Flow123d  JS_before_hm-62-ge56f9d5
mesh.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file mesh.cc
15  * @ingroup mesh
16  * @brief Mesh construction
17  */
18 
19 #include <unistd.h>
20 #include <set>
21 
22 
23 #include "system/system.hh"
24 #include "system/exceptions.hh"
26 #include "input/input_type.hh"
27 #include "input/accessors.hh"
28 #include "system/sys_profiler.hh"
29 #include "la/distribution.hh"
30 
31 #include "mesh/long_idx.hh"
32 #include "mesh/mesh.h"
33 #include "mesh/bc_mesh.hh"
34 #include "mesh/ref_element.hh"
35 #include "mesh/region_set.hh"
36 #include "mesh/range_wrapper.hh"
37 
38 // think about following dependencies
39 #include "mesh/accessors.hh"
40 #include "mesh/node_accessor.hh"
41 #include "mesh/partitioning.hh"
42 #include "mesh/neighbours.h"
43 
44 
45 #include "mesh/bih_tree.hh"
46 #include "mesh/duplicate_nodes.h"
47 
49 
50 
51 //TODO: sources, concentrations, initial condition and similarly boundary conditions should be
52 // instances of a Element valued field
53 // concentrations is in fact reimplemented in transport REMOVE it HERE
54 
55 // After removing non-geometrical things from mesh, this should be part of mash initializing.
56 #include "mesh/region.hh"
57 
58 #define NDEF -1
59 
60 namespace IT = Input::Type;
61 
63  return Input::Type::Selection("Types of search algorithm for finding intersection candidates.")
64  .add_value(Mesh::BIHsearch, "BIHsearch",
65  "Use BIH for finding initial candidates, then continue by prolongation.")
66  .add_value(Mesh::BIHonly, "BIHonly",
67  "Use BIH for finding all candidates.")
68  .add_value(Mesh::BBsearch, "BBsearch",
69  "Use bounding boxes for finding initial candidates, then continue by prolongation.")
70  .close();
71 }
72 
74  return IT::Record("Mesh","Record with mesh related data." )
75  .allow_auto_conversion("mesh_file")
77  "Input file with mesh description.")
79  "List of additional region and region set definitions not contained in the mesh. "
80  "There are three region sets implicitly defined:\n\n"
81  "- ALL (all regions of the mesh)\n"
82  "- .BOUNDARY (all boundary regions)\n"
83  "- BULK (all bulk regions)")
84  .declare_key("partitioning", Partitioning::get_input_type(), IT::Default("\"any_neighboring\""), "Parameters of mesh partitioning algorithms.\n" )
85  .declare_key("print_regions", IT::Bool(), IT::Default("true"), "If true, print table of all used regions.")
86  .declare_key("intersection_search", Mesh::get_input_intersection_variant(),
87  IT::Default("\"BIHsearch\""), "Search algorithm for element intersections.")
88  .declare_key("global_snap_radius", IT::Double(0.0), IT::Default("1E-3"),
89  "Maximal snapping distance from the mesh in various search operations. In particular, it is used "
90  "to find the closest mesh element of an observe point; and in FieldFormula to find closest surface "
91  "element in plan view (Z projection).")
93  "Output file with neighboring data from mesh.")
94  .close();
95 }
96 
97 const unsigned int Mesh::undef_idx;
98 
100 : tree(nullptr),
101  bulk_size_(0),
102  row_4_el(nullptr),
103  el_4_loc(nullptr),
104  el_ds(nullptr),
105  node_4_loc_(nullptr),
106  node_ds_(nullptr),
107  bc_mesh_(nullptr)
108 
109 {}
110 
111 
112 
114 : tree(nullptr),
115  in_record_(in_record),
116  comm_(com),
117  row_4_el(nullptr),
118  el_4_loc(nullptr),
119  el_ds(nullptr),
120  node_4_loc_(nullptr),
121  node_ds_(nullptr),
122  bc_mesh_(nullptr)
123 {
124  // set in_record_, if input accessor is empty
125  if (in_record_.is_empty()) {
126  istringstream is("{mesh_file=\"\"}");
127  Input::ReaderToStorage reader;
128  IT::Record &in_rec = const_cast<IT::Record &>(Mesh::get_input_type());
129  in_rec.finish();
130  reader.read_stream(is, in_rec, Input::FileFormat::format_JSON);
132  }
133 
134  int rank;
136  if (rank == 0) {
137  // optionally open raw output file
138  FilePath raw_output_file_path;
139  if (in_record_.opt_val("raw_ngh_output", raw_output_file_path)) {
140  MessageOut() << "Opening raw ngh output: " << raw_output_file_path << "\n";
141  try {
142  raw_output_file_path.open_stream(raw_ngh_output_file);
143  } INPUT_CATCH(FilePath::ExcFileOpen, FilePath::EI_Address_String, (in_record_))
144  }
145 
146  }
147  init();
148 }
149 
151 {
152  return in_record_.val<Mesh::IntersectionSearch>("intersection_search");
153 }
154 
155 
157 {
158 
159  n_insides = NDEF;
160  n_exsides = NDEF;
161  n_sides_ = NDEF;
162 
163  // number of element of particular dimension
164  n_lines = 0;
165  n_triangles = 0;
166  n_tetrahedras = 0;
167 
168  for (int d=0; d<3; d++) max_edge_sides_[d] = 0;
169 
170  // Initialize numbering of nodes on sides.
171  // This is temporary solution, until class Element is templated
172  // by dimension. Then we can replace Mesh::side_nodes by
173  // RefElement<dim>::side_nodes.
174 
175  // indices of side nodes in element node array
176  // Currently this is made ad libitum
177  // with some ordering here we can get sides with correct orientation.
178  // This speedup normal calculation.
179 
180  side_nodes.resize(3); // three side dimensions
181  for(int i=0; i < 3; i++) {
182  side_nodes[i].resize(i+2); // number of sides
183  for(int j=0; j < i+2; j++)
184  side_nodes[i][j].resize(i+1);
185  }
186 
187  for (unsigned int sid=0; sid<RefElement<1>::n_sides; sid++)
188  for (unsigned int nid=0; nid<RefElement<1>::n_nodes_per_side; nid++)
189  side_nodes[0][sid][nid] = RefElement<1>::interact(Interaction<0,0>(sid))[nid];
190 
191  for (unsigned int sid=0; sid<RefElement<2>::n_sides; sid++)
192  for (unsigned int nid=0; nid<RefElement<2>::n_nodes_per_side; nid++)
193  side_nodes[1][sid][nid] = RefElement<2>::interact(Interaction<0,1>(sid))[nid];
194 
195  for (unsigned int sid=0; sid<RefElement<3>::n_sides; sid++)
196  for (unsigned int nid=0; nid<RefElement<3>::n_nodes_per_side; nid++)
197  side_nodes[2][sid][nid] = RefElement<3>::interact(Interaction<0,2>(sid))[nid];
198 }
199 
200 
202  for(EdgeData &edg : this->edges)
203  if (edg.side_) delete[] edg.side_;
204 
205  for (unsigned int idx=0; idx < bulk_size_; idx++) {
206  Element *ele=&(element_vec_[idx]);
207  if (ele->boundary_idx_) delete[] ele->boundary_idx_;
208  if (ele->neigh_vb) delete[] ele->neigh_vb;
209  }
210 
211  for(unsigned int idx=bulk_size_; idx < element_vec_.size(); idx++) {
212  Element *ele=&(element_vec_[idx]);
213  if (ele->boundary_idx_) delete[] ele->boundary_idx_;
214  }
215 
216  if (row_4_el != nullptr) delete[] row_4_el;
217  if (el_4_loc != nullptr) delete[] el_4_loc;
218  if (el_ds != nullptr) delete el_ds;
219  if (node_4_loc_ != nullptr) delete[] node_4_loc_;
220  if (node_ds_ != nullptr) delete node_ds_;
221  if (bc_mesh_ != nullptr) delete bc_mesh_;
222  if (tree != nullptr) delete tree;
223 }
224 
225 
226 unsigned int Mesh::n_sides() const
227 {
228  if (n_sides_ == NDEF) {
229  n_sides_=0;
230  for (auto ele : this->elements_range()) n_sides_ += ele->n_sides();
231  }
232  return n_sides_;
233 }
234 
235 unsigned int Mesh::n_vb_neighbours() const {
236  return vb_neighbours_.size();
237  }
238 
239 
240 unsigned int Mesh::n_corners() {
241  unsigned int li, count = 0;
242  for (auto ele : this->elements_range()) {
243  for (li=0; li<ele->n_nodes(); li++) {
244  count++;
245  }
246  }
247  return count;
248 }
249 
250 Edge Mesh::edge(uint edge_idx) const
251 {
252  ASSERT_LT_DBG(edge_idx, edges.size());
253  return Edge(this, edge_idx);
254 }
255 
257 {
258  ASSERT_LT_DBG(bc_idx, boundary_.size());
259  return Boundary(&boundary_[bc_idx]);
260 }
261 
263  return part_.get();
264 }
265 
267  return (LongIdx*)this->get_part()->get_loc_part();
268 }
269 
270 
271 //=============================================================================
272 // COUNT ELEMENT TYPES
273 //=============================================================================
274 
276  for (auto elm : this->elements_range())
277  switch (elm->dim()) {
278  case 1:
279  n_lines++;
280  break;
281  case 2:
282  n_triangles++;
283  break;
284  case 3:
285  n_tetrahedras++;
286  break;
287  }
288 }
289 
290 
291 
293  for (auto elem_to_region : map) {
294  Element &ele = element_vec_[ elem_index(elem_to_region.first) ];
295  ele.region_idx_ = region_db_.get_region( elem_to_region.second, ele.dim() );
297  }
298 }
299 
300 
301 
303  START_TIMER("MESH - setup topology");
304 
306 
307  // check mesh quality
308  for (auto ele : this->elements_range()) {
309  double quality = ele.quality_measure_smooth(ele.side(0));
310  if ( quality< 0.001) WarningOut().fmt("Bad quality (<0.001) of the element {}.\n", ele.idx());
311  }
312 
317 
318  tree = new DuplicateNodes(this);
319 
320  part_ = std::make_shared<Partitioning>(this, in_record_.val<Input::Record>("partitioning") );
321 
322  // create parallel distribution and numbering of elements
323  LongIdx *id_4_old = new LongIdx[n_elements()];
324  int i = 0;
325  for (auto ele : this->elements_range())
326  id_4_old[i++] = ele.idx();
327  part_->id_maps(n_elements(), id_4_old, el_ds, el_4_loc, row_4_el);
328 
329  delete[] id_4_old;
330 
331  this->distribute_nodes();
332 
334 }
335 
336 
337 //
339 {
340 
341  n_insides = 0;
342  n_exsides = 0;
343  for (auto ele : this->elements_range())
344  for(SideIter sde = ele.side(0); sde->side_idx() < ele->n_sides(); ++sde) {
345  if (sde->is_external()) n_exsides++;
346  else n_insides++;
347  }
348 }
349 
350 
351 
353  // for each node we make a list of elements that use this node
354  node_elements_.resize( this->n_nodes() );
355 
356  for (auto ele : this->elements_range())
357  for (unsigned int n=0; n<ele->n_nodes(); n++)
358  node_elements_[ele.node_accessor(n).idx()].push_back(ele.idx());
359 
360  for (vector<vector<unsigned int> >::iterator n=node_elements_.begin(); n!=node_elements_.end(); n++)
361  stable_sort(n->begin(), n->end());
362 }
363 
364 
365 void Mesh::intersect_element_lists(vector<unsigned int> const &nodes_list, vector<unsigned int> &intersection_element_list)
366 {
367  if (node_elements_.size() == 0) {
369  }
370 
371  if (nodes_list.size() == 0) {
372  intersection_element_list.clear();
373  } else if (nodes_list.size() == 1) {
374  intersection_element_list = node_elements_[ nodes_list[0] ];
375  } else {
376  vector<unsigned int>::const_iterator it1=nodes_list.begin();
378  intersection_element_list.resize( node_elements_[*it1].size() ); // make enough space
379 
380  it1=set_intersection(
381  node_elements_[*it1].begin(), node_elements_[*it1].end(),
382  node_elements_[*it2].begin(), node_elements_[*it2].end(),
383  intersection_element_list.begin());
384  intersection_element_list.resize(it1-intersection_element_list.begin()); // resize to true size
385 
386  for(;it2<nodes_list.end();++it2) {
387  it1=set_intersection(
388  intersection_element_list.begin(), intersection_element_list.end(),
389  node_elements_[*it2].begin(), node_elements_[*it2].end(),
390  intersection_element_list.begin());
391  intersection_element_list.resize(it1-intersection_element_list.begin()); // resize to true size
392  }
393  }
394 }
395 
396 
397 bool Mesh::find_lower_dim_element( vector<unsigned int> &element_list, unsigned int dim, unsigned int &element_idx) {
398  bool is_neighbour = false;
399 
400  vector<unsigned int>::iterator e_dest=element_list.begin();
401  for( vector<unsigned int>::iterator ele = element_list.begin(); ele!=element_list.end(); ++ele) {
402  if (element_vec_[*ele].dim() == dim) { // keep only indexes of elements of same dimension
403  *e_dest=*ele;
404  ++e_dest;
405  } else if (element_vec_[*ele].dim() == dim-1) { // get only first element of lower dimension
406  if (is_neighbour) xprintf(UsrErr, "Too matching elements id: %d and id: %d in the same mesh.\n",
407  this->elem_index(*ele), this->elem_index(element_idx) );
408 
409  is_neighbour = true;
410  element_idx = *ele;
411  }
412  }
413  element_list.resize( e_dest - element_list.begin());
414  return is_neighbour;
415 }
416 
418  // check if nodes lists match (this is slow and will be faster only when we convert whole mesh into hierarchical design like in deal.ii)
419  unsigned int ni=0;
420  while ( ni < si->n_nodes()
421  && find(side_nodes.begin(), side_nodes.end(), si->node(ni).idx() ) != side_nodes.end() ) ni++;
422  return ( ni == si->n_nodes() );
423 }
424 
425 /**
426  * TODO:
427  * - use std::is_any for setting is_neigbour
428  * - possibly make appropriate constructors for Edge and Neighbour
429  * - check side!=-1 when searching neigbouring element
430  * - process boundary elements first, there should be no Neigh, but check it
431  * set Edge and boundary there
432  */
433 
435 {
436  ASSERT(bc_element_tmp_.size()==0)
437  .error("Temporary structure of boundary element data is not empty. Did you call create_boundary_elements?");
438 
439  Neighbour neighbour;
440  EdgeData *edg = nullptr;
441  unsigned int ngh_element_idx;
442  unsigned int last_edge_idx = Mesh::undef_idx;
443 
444  neighbour.mesh_ = this;
445 
447 
448  // pointers to created edges
449  //vector<Edge *> tmp_edges;
450  edges.resize(0); // be sure that edges are empty
451 
453  vector<unsigned int> intersection_list; // list of elements in intersection of node element lists
454 
455  for( unsigned int i=bulk_size_; i<element_vec_.size(); ++i) {
456 
457  ElementAccessor<3> bc_ele = this->element_accessor(i);
458  ASSERT(bc_ele.region().is_boundary());
459  // Find all elements that share this side.
460  side_nodes.resize(bc_ele->n_nodes());
461  for (unsigned n=0; n<bc_ele->n_nodes(); n++) side_nodes[n] = bc_ele->node_idx(n);
462  intersect_element_lists(side_nodes, intersection_list);
463  bool is_neighbour = find_lower_dim_element(intersection_list, bc_ele->dim() +1, ngh_element_idx);
464  if (is_neighbour) {
465  xprintf(UsrErr, "Boundary element (id: %d) match a regular element (id: %d) of lower dimension.\n",
466  bc_ele.idx(), this->elem_index(ngh_element_idx));
467  } else {
468  if (intersection_list.size() == 0) {
469  // no matching dim+1 element found
470  WarningOut().fmt("Lonely boundary element, id: {}, region: {}, dimension {}.\n",
471  bc_ele.idx(), bc_ele.region().id(), bc_ele->dim());
472  continue; // skip the boundary element
473  }
474  last_edge_idx=edges.size();
475  edges.resize(last_edge_idx+1);
476  edg = &( edges.back() );
477  edg->n_sides = 0;
478  edg->side_ = new struct SideIter[ intersection_list.size() ];
479 
480  // common boundary object
481  unsigned int bdr_idx=boundary_.size();
482  boundary_.resize(bdr_idx+1);
483  BoundaryData &bdr=boundary_.back();
484  bdr.bc_ele_idx_ = i;
485  bdr.edge_idx_ = last_edge_idx;
486  bdr.mesh_=this;
487 
488  // for 1d boundaries there can be more then one 1d elements connected to the boundary element
489  // we do not detect this case later in the main search over bulk elements
490  for( vector<unsigned int>::iterator isect = intersection_list.begin(); isect!=intersection_list.end(); ++isect) {
491  ElementAccessor<3> elem = this->element_accessor(*isect);
492  for (unsigned int ecs=0; ecs<elem->n_sides(); ecs++) {
493  SideIter si = elem.side(ecs);
494  if ( same_sides( si, side_nodes) ) {
495  if (elem->edge_idx(ecs) != Mesh::undef_idx) {
496  OLD_ASSERT(elem->boundary_idx_!=nullptr, "Null boundary idx array.\n");
497  int last_bc_ele_idx=this->boundary_[elem->boundary_idx_[ecs]].bc_ele_idx_;
498  int new_bc_ele_idx=i;
499  THROW( ExcDuplicateBoundary()
500  << EI_ElemLast(this->find_elem_id(last_bc_ele_idx))
501  << EI_RegLast(this->element_accessor(last_bc_ele_idx).region().label())
502  << EI_ElemNew(this->find_elem_id(new_bc_ele_idx))
503  << EI_RegNew(this->element_accessor(new_bc_ele_idx).region().label())
504  );
505  }
506  element_vec_[*isect].edge_idx_[ecs] = last_edge_idx;
507  edg->side_[ edg->n_sides++ ] = si;
508 
509  if (elem->boundary_idx_ == NULL) {
510  Element *el = &(element_vec_[*isect]);
511  el->boundary_idx_ = new unsigned int [ el->n_sides() ];
512  std::fill( el->boundary_idx_, el->boundary_idx_ + el->n_sides(), Mesh::undef_idx);
513  }
514  elem->boundary_idx_[ecs] = bdr_idx;
515  break; // next element in intersection list
516  }
517  }
518  }
519 
520  }
521 
522  }
523  // Now we go through all element sides and create edges and neighbours
524  for (auto e : this->elements_range()) {
525  for (unsigned int s=0; s<e->n_sides(); s++)
526  {
527  // skip sides that were already found
528  if (e->edge_idx(s) != Mesh::undef_idx) continue;
529 
530 
531  // Find all elements that share this side.
532  side_nodes.resize(e.side(s)->n_nodes());
533  for (unsigned n=0; n<e.side(s)->n_nodes(); n++) side_nodes[n] = e.side(s)->node(n).idx();
534  intersect_element_lists(side_nodes, intersection_list);
535 
536  bool is_neighbour = find_lower_dim_element(intersection_list, e->dim(), ngh_element_idx);
537 
538  if (is_neighbour) { // edge connects elements of different dimensions
539  // Initialize for the neighbour case.
540  neighbour.elem_idx_ = ngh_element_idx;
541  } else { // edge connects only elements of the same dimension
542  // Initialize for the edge case.
543  last_edge_idx=edges.size();
544  edges.resize(last_edge_idx+1);
545  edg = &( edges.back() );
546  edg->n_sides = 0;
547  edg->side_ = new struct SideIter[ intersection_list.size() ];
548  if (intersection_list.size() > max_edge_sides_[e->dim()-1])
549  max_edge_sides_[e->dim()-1] = intersection_list.size();
550 
551  if (intersection_list.size() == 1) {
552  // outer edge, create boundary object as well
553  Element &elm = element_vec_[e.idx()];
554  edg->n_sides=1;
555  edg->side_[0] = e.side(s);
556  element_vec_[e.idx()].edge_idx_[s] = last_edge_idx;
557 
558  if (e->boundary_idx_ == NULL) {
559  elm.boundary_idx_ = new unsigned int [ e->n_sides() ];
560  std::fill( elm.boundary_idx_, elm.boundary_idx_ + e->n_sides(), Mesh::undef_idx);
561  }
562 
563  unsigned int bdr_idx=boundary_.size()+1; // need for VTK mesh that has no boundary elements
564  // and bulk elements are indexed from 0
565  boundary_.resize(bdr_idx+1);
566  BoundaryData &bdr=boundary_.back();
567  elm.boundary_idx_[s] = bdr_idx;
568 
569  // fill boundary element
570  Element * bc_ele = add_element_to_vector(-bdr_idx, true);
571  bc_ele->init(e->dim()-1, region_db_.implicit_boundary_region() );
573  for(unsigned int ni = 0; ni< side_nodes.size(); ni++) bc_ele->nodes_[ni] = side_nodes[ni];
574 
575  // fill Boundary object
576  bdr.edge_idx_ = last_edge_idx;
577  bdr.bc_ele_idx_ = elem_index(-bdr_idx);
578  bdr.mesh_=this;
579 
580  continue; // next side of element e
581  }
582  }
583 
584  // go through the elements connected to the edge or neighbour
585  // setup neigbour or edge
586  for( vector<unsigned int>::iterator isect = intersection_list.begin(); isect!=intersection_list.end(); ++isect) {
587  ElementAccessor<3> elem = this->element_accessor(*isect);
588  for (unsigned int ecs=0; ecs<elem->n_sides(); ecs++) {
589  if (elem->edge_idx(ecs) != Mesh::undef_idx) continue; // ??? This should not happen.
590  SideIter si = elem.side(ecs);
591  if ( same_sides( si, side_nodes) ) {
592  if (is_neighbour) {
593  // create a new edge and neighbour for this side, and element to the edge
594  last_edge_idx=edges.size();
595  edges.resize(last_edge_idx+1);
596  edg = &( edges.back() );
597  edg->n_sides = 1;
598  edg->side_ = new struct SideIter[1];
599  edg->side_[0] = si;
600  element_vec_[elem.idx()].edge_idx_[ecs] = last_edge_idx;
601 
602  neighbour.edge_idx_ = last_edge_idx;
603 
604  vb_neighbours_.push_back(neighbour); // copy neighbour with this edge setting
605  } else {
606  // connect the side to the edge, and side to the edge
607  ASSERT_PTR_DBG(edg);
608  edg->side_[ edg->n_sides++ ] = si;
609  ASSERT_DBG(last_edge_idx != Mesh::undef_idx);
610  element_vec_[elem.idx()].edge_idx_[ecs] = last_edge_idx;
611  }
612  break; // next element from intersection list
613  }
614  } // search for side of other connected element
615  } // connected elements
616 
617  if (! is_neighbour)
618  ASSERT_EQ( (unsigned int) edg->n_sides, intersection_list.size())(e.index())(s).error("Missing edge sides.");
619  } // for element sides
620  } // for elements
621 
622  MessageOut().fmt( "Created {} edges and {} neighbours.\n", edges.size(), vb_neighbours_.size() );
623 }
624 
625 
626 
628 {
629  for (auto edg : edge_range())
630  {
631  // side 0 is reference, so its permutation is 0
632  edg.side(0)->element()->permutation_idx_[edg.side(0)->side_idx()] = 0;
633 
634  if (edg.n_sides() > 1)
635  {
636  map<unsigned int,unsigned int> node_numbers;
637  unsigned int permutation[edg.side(0)->n_nodes()];
638 
639  for (unsigned int i=0; i<edg.side(0)->n_nodes(); i++)
640  node_numbers[edg.side(0)->node(i).idx()] = i;
641  //node_numbers[edg.side(0)->node(i).node()] = i;
642 
643  for (uint sid=1; sid<edg.n_sides(); sid++)
644  {
645  for (unsigned int i=0; i<edg.side(0)->n_nodes(); i++)
646  permutation[node_numbers[edg.side(sid)->node(i).idx()]] = i;
647  //permutation[node_numbers[edg.side(sid)->node(i).node()]] = i;
648 
649  switch (edg.side(0)->dim())
650  {
651  case 0:
652  edg.side(sid)->element()->permutation_idx_[edg.side(sid)->side_idx()] = RefElement<1>::permutation_index(permutation);
653  break;
654  case 1:
655  edg.side(sid)->element()->permutation_idx_[edg.side(sid)->side_idx()] = RefElement<2>::permutation_index(permutation);
656  break;
657  case 2:
658  edg.side(sid)->element()->permutation_idx_[edg.side(sid)->side_idx()] = RefElement<3>::permutation_index(permutation);
659  break;
660  }
661  }
662  }
663  }
664 
665  for (vector<Neighbour>::iterator nb=vb_neighbours_.begin(); nb!=vb_neighbours_.end(); nb++)
666  {
667  map<const Node*,unsigned int> node_numbers;
668  unsigned int permutation[nb->element()->n_nodes()];
669 
670  // element of lower dimension is reference, so
671  // we calculate permutation for the adjacent side
672  for (unsigned int i=0; i<nb->element()->n_nodes(); i++)
673  node_numbers[nb->element().node(i)] = i;
674 
675  for (unsigned int i=0; i<nb->side()->n_nodes(); i++)
676  permutation[node_numbers[nb->side()->node(i).node()]] = i;
677 
678  switch (nb->side()->dim())
679  {
680  case 0:
681  nb->side()->element()->permutation_idx_[nb->side()->side_idx()] = RefElement<1>::permutation_index(permutation);
682  break;
683  case 1:
684  nb->side()->element()->permutation_idx_[nb->side()->side_idx()] = RefElement<2>::permutation_index(permutation);
685  break;
686  case 2:
687  nb->side()->element()->permutation_idx_[nb->side()->side_idx()] = RefElement<3>::permutation_index(permutation);
688  break;
689  }
690  }
691 }
692 
693 
694 
695 
696 
697 //=============================================================================
698 //
699 //=============================================================================
701 {
702 
703  //MessageOut() << "Element to neighbours of vb2 type... "/*orig verb 5*/;
704 
705  for (vector<Element>::iterator ele = element_vec_.begin(); ele!= element_vec_.begin()+bulk_size_; ++ele)
706  ele->n_neighs_vb_ =0;
707 
708  // count vb neighs per element
709  for (auto & ngh : this->vb_neighbours_) ngh.element()->n_neighs_vb_++;
710 
711  // Allocation of the array per element
712  for (vector<Element>::iterator ele = element_vec_.begin(); ele!= element_vec_.begin()+bulk_size_; ++ele)
713  if( ele->n_neighs_vb() > 0 ) {
714  ele->neigh_vb = new struct Neighbour* [ele->n_neighs_vb()];
715  ele->n_neighs_vb_=0;
716  }
717 
718  // fill
719  ElementAccessor<3> ele;
720  for (auto & ngh : this->vb_neighbours_) {
721  ele = ngh.element();
722  ele->neigh_vb[ ele->n_neighs_vb_++ ] = &ngh;
723  }
724 
725  //MessageOut() << "... O.K.\n"/*orig verb 6*/;
726 }
727 
728 
729 
730 
731 
732 
734  /* Algorithm:
735  *
736  * 1) create BIH tree
737  * 2) for every 1D, find list of candidates
738  * 3) compute intersections for 1d, store it to master_elements
739  *
740  */
741  if (! intersections) {
742  intersections = std::make_shared<MixedMeshIntersections>(this);
743  intersections->compute_intersections();
744  }
745  return *intersections;
746 }
747 
748 
749 
751  return ElementAccessor<3>(this, idx);
752 }
753 
754 
755 
756 NodeAccessor<3> Mesh::node_accessor(unsigned int idx) const {
757  return NodeAccessor<3>(this, idx);
758 }
759 
760 
761 
762 void Mesh::elements_id_maps( vector<LongIdx> & bulk_elements_id, vector<LongIdx> & boundary_elements_id) const
763 {
764  if (bulk_elements_id.size() ==0) {
766  LongIdx last_id;
767 
768  bulk_elements_id.resize(n_elements());
769  map_it = bulk_elements_id.begin();
770  last_id = -1;
771  for(unsigned int idx=0; idx < n_elements(); idx++, ++map_it) {
772  LongIdx id = this->find_elem_id(idx);
773  if (last_id >= id) xprintf(UsrErr, "Element IDs in non-increasing order, ID: %d\n", id);
774  last_id=*map_it = id;
775  }
776 
777  boundary_elements_id.resize(n_elements(true));
778  map_it = boundary_elements_id.begin();
779  last_id = -1;
780  for(unsigned int idx=bulk_size_; idx<element_vec_.size(); idx++, ++map_it) {
781  LongIdx id = this->find_elem_id(idx);
782  // We set ID for boundary elements created by the mesh itself to "-1"
783  // this force gmsh reader to skip all remaining entries in boundary_elements_id
784  // and thus report error for any remaining data lines
785  if (id < 0) last_id=*map_it=-1;
786  else {
787  if (last_id >= id) xprintf(UsrErr, "Element IDs in non-increasing order, ID: %d\n", id);
788  last_id=*map_it = id;
789  }
790  }
791  }
792 }
793 
794 
795 bool compare_points(const arma::vec3 &p1, const arma::vec3 &p2) {
796  static const double point_tolerance = 1E-10;
797  return fabs(p1[0]-p2[0]) < point_tolerance
798  && fabs(p1[1]-p2[1]) < point_tolerance
799  && fabs(p1[2]-p2[2]) < point_tolerance;
800 }
801 
802 
803 bool Mesh::check_compatible_mesh( Mesh & mesh, vector<LongIdx> & bulk_elements_id, vector<LongIdx> & boundary_elements_id )
804 {
805  std::vector<unsigned int> node_ids; // allow mapping ids of nodes from source mesh to target mesh
806  std::vector<unsigned int> node_list;
807  std::vector<unsigned int> candidate_list; // returned by intersect_element_lists
808  std::vector<unsigned int> result_list; // list of elements with same dimension as vtk element
809  unsigned int i; // counter over vectors
810 
811  {
812  // iterates over node vector of \p this object
813  // to each node must be found just only one node in target \p mesh
814  // store orders (mapping between source and target meshes) into node_ids vector
815  std::vector<unsigned int> searched_elements; // for BIH tree
816  unsigned int i_node, i_elm_node;
817  const BIHTree &bih_tree=mesh.get_bih_tree();
818 
819  // create nodes of mesh
820  node_ids.resize( this->n_nodes() );
821  i=0;
822  for (auto nod : this->node_range()) {
823  arma::vec3 point = nod->point();
824  uint found_i_node = Mesh::undef_idx;
825  bih_tree.find_point(point, searched_elements);
826 
827  for (std::vector<unsigned int>::iterator it = searched_elements.begin(); it!=searched_elements.end(); it++) {
828  ElementAccessor<3> ele = mesh.element_accessor( *it );
829  for (i_node=0; i_node<ele->n_nodes(); i_node++)
830  {
831  if ( compare_points(ele.node(i_node)->point(), point) ) {
832  i_elm_node = ele.node_accessor(i_node).idx();
833  if (found_i_node == Mesh::undef_idx) found_i_node = i_elm_node;
834  else if (found_i_node != i_elm_node) {
835  // duplicate nodes in target mesh
836  this->elements_id_maps(bulk_elements_id, boundary_elements_id);
837  return false;
838  }
839  }
840  }
841  }
842  if (found_i_node == Mesh::undef_idx) {
843  // no node found in target mesh
844  this->elements_id_maps(bulk_elements_id, boundary_elements_id);
845  return false;
846  }
847  node_ids[i] = found_i_node;
848  searched_elements.clear();
849  i++;
850  }
851  }
852 
853  {
854  // iterates over bulk elements of \p this object
855  // elements in both meshes must be in ratio 1:1
856  // store orders (mapping between both mesh files) into bulk_elements_id vector
857  bulk_elements_id.clear();
858  bulk_elements_id.resize(this->n_elements());
859  // iterate trough bulk part of element vector, to each element in source mesh must exist only one element in target mesh
860  // fill bulk_elements_id vector
861  i=0;
862  for (auto elm : this->elements_range()) {
863  for (unsigned int j=0; j<elm->n_nodes(); j++) { // iterate trough all nodes of any element
864  node_list.push_back( node_ids[ elm->node_idx(j) ] );
865  }
866  mesh.intersect_element_lists(node_list, candidate_list);
867  for (auto i_elm : candidate_list) {
868  if ( mesh.element_accessor(i_elm)->dim() == elm.dim() ) result_list.push_back( elm.index() );
869  }
870  if (result_list.size() != 1) {
871  // intersect_element_lists must produce one element
872  this->elements_id_maps(bulk_elements_id, boundary_elements_id);
873  return false;
874  }
875  bulk_elements_id[i] = (LongIdx)result_list[0];
876  node_list.clear();
877  result_list.clear();
878  i++;
879  }
880  }
881 
882  {
883  // iterates over boundary elements of \p this object
884  // elements in both meshes must be in ratio 1:1
885  // store orders (mapping between both mesh files) into boundary_elements_id vector
886  auto bc_mesh = this->get_bc_mesh();
887  boundary_elements_id.clear();
888  boundary_elements_id.resize(bc_mesh->n_elements());
889  // iterate trough boundary part of element vector, to each element in source mesh must exist only one element in target mesh
890  // fill boundary_elements_id vector
891  i=0;
892  for (auto elm : bc_mesh->elements_range()) {
893  for (unsigned int j=0; j<elm->n_nodes(); j++) { // iterate trough all nodes of any element
894  node_list.push_back( node_ids[ elm->node_idx(j) ] );
895  }
896  mesh.get_bc_mesh()->intersect_element_lists(node_list, candidate_list);
897  for (auto i_elm : candidate_list) {
898  if ( mesh.get_bc_mesh()->element_accessor(i_elm)->dim() == elm.dim() ) result_list.push_back( elm.index() );
899  }
900  if (result_list.size() != 1) {
901  // intersect_element_lists must produce one element
902  this->elements_id_maps(bulk_elements_id, boundary_elements_id);
903  return false;
904  }
905  boundary_elements_id[i] = (LongIdx)result_list[0];
906  node_list.clear();
907  result_list.clear();
908  i++;
909  }
910  }
911 
912  return true;
913 }
914 
916 {
918  it != region_list.end();
919  ++it) {
920  // constructor has side effect in the mesh - create new region or set and store them to Mesh::region_db_
921  (*it).factory< RegionSetBase, const Input::Record &, Mesh * >(*it, this);
922  }
923 }
924 
926 {
928  region_db_.el_to_reg_map_.clear();
929  region_db_.close();
931 
932  if ( in_record_.val<bool>("print_regions") ) {
933  stringstream ss;
935  MessageOut() << ss.str();
936  }
937 }
938 
939 
941  START_TIMER("Mesh::compute_element_boxes");
943 
944  // make element boxes
945  unsigned int i=0;
946  boxes.resize(this->n_elements());
947  for (auto element : this->elements_range()) {
948  boxes[i] = element.bounding_box();
949  i++;
950  }
951 
952  return boxes;
953 }
954 
956  if (! this->bih_tree_) {
957  bih_tree_ = std::make_shared<BIHTree>();
958  bih_tree_->add_boxes( this->get_element_boxes() );
959  bih_tree_->construct();
960  }
961  return *bih_tree_;
962 }
963 
964 double Mesh::global_snap_radius() const {
965  return in_record_.val<double>("global_snap_radius");
966 }
967 
968 void Mesh::add_physical_name(unsigned int dim, unsigned int id, std::string name) {
969  region_db_.add_region(id, name, dim, "$PhysicalNames");
970 }
971 
972 
973 void Mesh::add_node(unsigned int node_id, arma::vec3 coords) {
974  node_vec_.push_back( Node() );
975  Node &node = node_vec_[node_vec_.size()-1];
976  node.point() = coords;
977  node_ids_.add_item(node_id);
978 }
979 
980 
981 void Mesh::add_element(unsigned int elm_id, unsigned int dim, unsigned int region_id, unsigned int partition_id,
982  std::vector<unsigned int> node_ids) {
983  RegionIdx region_idx = region_db_.get_region( region_id, dim );
984  if ( !region_idx.is_valid() ) {
985  region_idx = region_db_.add_region( region_id, region_db_.create_label_from_id(region_id), dim, "$Element" );
986  }
987  region_db_.mark_used_region(region_idx.idx());
988 
989  if (region_idx.is_boundary()) {
990  bc_element_tmp_.push_back( ElementTmpData(elm_id, dim, region_idx, partition_id, node_ids) );
991  } else {
992  if(dim == 0 ) {
993  WarningOut().fmt("Bulk elements of zero size(dim=0) are not supported. Element ID: {}.\n", elm_id);
994  }
995  else {
996  Element *ele = add_element_to_vector(elm_id);
997  this->init_element(ele, elm_id, dim, region_idx, partition_id, node_ids);
998  }
999  }
1000 }
1001 
1002 
1003 void Mesh::init_element(Element *ele, unsigned int elm_id, unsigned int dim, RegionIdx region_idx, unsigned int partition_id,
1004  std::vector<unsigned int> node_ids) {
1005  ele->init(dim, region_idx);
1006  ele->pid_ = partition_id;
1007 
1008  for (unsigned int ni=0; ni<ele->n_nodes(); ni++) {
1009  ele->nodes_[ni] = this->node_index(node_ids[ni]);
1010  }
1011 
1012  // check that tetrahedron element is numbered correctly and is not degenerated
1013  if(ele->dim() == 3)
1014  {
1015  double jac = this->element_accessor( this->elem_index(elm_id) ).tetrahedron_jacobian();
1016  if( ! (jac > 0) )
1017  WarningOut().fmt("Tetrahedron element with id {} has wrong numbering or is degenerated (Jacobian = {}).",elm_id,jac);
1018  }
1019 }
1020 
1021 
1023  if (node_elements_.size() == 0) {
1024  this->create_node_element_lists();
1025  }
1026  return node_elements_;
1027 }
1028 
1029 
1030 void Mesh::init_element_vector(unsigned int size) {
1031  element_vec_.clear();
1032  element_vec_.resize(size);
1033  element_ids_.reinit(size);
1034  bc_element_tmp_.clear();
1035  bc_element_tmp_.reserve(size);
1036  bulk_size_ = 0;
1038 }
1039 
1040 
1041 void Mesh::init_node_vector(unsigned int size) {
1042  node_vec_.clear();
1043  node_vec_.reserve(size);
1044  node_ids_.reinit(0);
1045 }
1046 
1047 
1049  Element * elem=nullptr;
1050  if (boundary) {
1051  ASSERT_DBG(id<0)(id).error("Add boundary element from mesh file trough temporary structure!");
1052  element_vec_.push_back( Element() );
1053  elem = &element_vec_[element_vec_.size()-1];
1054  element_ids_.add_item(id);
1055  } else {
1056  elem = &element_vec_[bulk_size_];
1058  bulk_size_++;
1059  }
1060 
1061  return elem;
1062 }
1063 
1065  auto bgn_it = make_iter<ElementAccessor<3>>( ElementAccessor<3>(this, 0) );
1066  auto end_it = make_iter<ElementAccessor<3>>( ElementAccessor<3>(this, bulk_size_) );
1067  return Range<ElementAccessor<3>>(bgn_it, end_it);
1068 }
1069 
1071  auto bgn_it = make_iter<NodeAccessor<3>>( NodeAccessor<3>(this, 0) );
1072  auto end_it = make_iter<NodeAccessor<3>>( NodeAccessor<3>(this, node_vec_.size()) );
1073  return Range<NodeAccessor<3>>(bgn_it, end_it);
1074 }
1075 
1077  auto bgn_it = make_iter<Edge>( Edge(this, 0) );
1078  auto end_it = make_iter<Edge>( Edge(this, edges.size()) );
1079  return Range<Edge>(bgn_it, end_it);
1080 }
1081 
1082 inline void Mesh::check_element_size(unsigned int elem_idx) const
1083 {
1084  ASSERT(elem_idx < element_vec_.size())(elem_idx)(element_vec_.size()).error("Index of element is out of bound of element vector!");
1085 }
1086 
1087 /*
1088  * Output of internal flow data.
1089  */
1091 {
1092  START_TIMER("Mesh::output_internal_ngh_data");
1093 
1094  if (! raw_ngh_output_file.is_open()) return;
1095 
1096  // header
1097  raw_ngh_output_file << "// fields:\n//ele_id n_sides ns_side_neighbors[n] neighbors[n*ns] n_vb_neighbors vb_neighbors[n_vb]\n";
1099 
1100  int cit = 0;
1101 
1102  // map from higher dim elements to its lower dim neighbors, using gmsh IDs: ele->id()
1103  unsigned int undefined_ele_id = -1;
1105  for (auto ele : this->elements_range()) {
1106  if(ele->n_neighs_vb() > 0){
1107  for (unsigned int i = 0; i < ele->n_neighs_vb(); i++){
1108  ElementAccessor<3> higher_ele = ele->neigh_vb[i]->side()->element();
1109 
1110  auto search = neigh_vb_map.find(higher_ele.idx());
1111  if(search != neigh_vb_map.end()){
1112  // if found, add id to correct local side idx
1113  search->second[ele->neigh_vb[i]->side()->side_idx()] = ele.idx();
1114  }
1115  else{
1116  // if not found, create new vector, each side can have one vb neighbour
1117  std::vector<unsigned int> higher_ele_side_ngh(higher_ele->n_sides(), undefined_ele_id);
1118  higher_ele_side_ngh[ele->neigh_vb[i]->side()->side_idx()] = ele.idx();
1119  neigh_vb_map[higher_ele.idx()] = higher_ele_side_ngh;
1120  }
1121  }
1122  }
1123  }
1124 
1125  for (auto ele : this->elements_range()) {
1126  raw_ngh_output_file << ele.idx() << " ";
1127  raw_ngh_output_file << ele->n_sides() << " ";
1128 
1129  auto search_neigh = neigh_vb_map.end();
1130  for (unsigned int i = 0; i < ele->n_sides(); i++) {
1131  unsigned int n_side_neighs = ele.side(i)->edge().n_sides()-1; //n_sides - the current one
1132  // check vb neighbors (lower dimension)
1133  if(n_side_neighs == 0){
1134  //update search
1135  if(search_neigh == neigh_vb_map.end())
1136  search_neigh = neigh_vb_map.find(ele.idx());
1137 
1138  if(search_neigh != neigh_vb_map.end())
1139  if(search_neigh->second[i] != undefined_ele_id)
1140  n_side_neighs = 1;
1141  }
1142  raw_ngh_output_file << n_side_neighs << " ";
1143  }
1144 
1145  for (unsigned int i = 0; i < ele->n_sides(); i++) {
1146  Edge edge = ele.side(i)->edge();
1147  if(edge.n_sides() > 1){
1148  for (uint j = 0; j < edge.n_sides(); j++) {
1149  if(edge.side(j) != ele.side(i))
1150  raw_ngh_output_file << edge.side(j)->element().idx() << " ";
1151  }
1152  }
1153  //check vb neighbour
1154  else if(search_neigh != neigh_vb_map.end()
1155  && search_neigh->second[i] != undefined_ele_id){
1156  raw_ngh_output_file << search_neigh->second[i] << " ";
1157  }
1158  }
1159 
1160  // list higher dim neighbours
1161  raw_ngh_output_file << ele->n_neighs_vb() << " ";
1162  for (unsigned int i = 0; i < ele->n_neighs_vb(); i++)
1163  raw_ngh_output_file << ele->neigh_vb[i]->side()->element().idx() << " ";
1164 
1165  raw_ngh_output_file << endl;
1166  cit ++;
1167  }
1168  raw_ngh_output_file << "$EndFlowField\n" << endl;
1169 }
1170 
1171 
1173  unsigned int i, pos;
1175  for (i=0, pos=bulk_size_; i<bc_element_tmp_.size(); ++i, ++pos) {
1176  Element *ele = &element_vec_[pos];
1177  element_ids_.set_item(bc_element_tmp_[i].elm_id, pos);
1178  this->init_element(ele, bc_element_tmp_[i].elm_id, bc_element_tmp_[i].dim, bc_element_tmp_[i].region_idx,
1179  bc_element_tmp_[i].partition_id, bc_element_tmp_[i].node_ids);
1180 
1181  }
1182 
1183  element_vec_.resize(pos); // remove empty element (count is equal with zero-dimensional bulk elements)
1184  bc_element_tmp_.clear();
1185  bc_element_tmp_.reserve(0);
1186 }
1187 
1188 
1189 void Mesh::permute_tetrahedron(unsigned int elm_idx, std::vector<unsigned int> permutation_vec)
1190 {
1191  ASSERT_LT_DBG(elm_idx, element_vec_.size());
1192  ASSERT_EQ_DBG(permutation_vec.size(), 4);
1193 
1194  std::array<unsigned int, 4> tmp_nodes;
1195  Element &elem = element_vec_[elm_idx];
1196  ASSERT_EQ_DBG(elem.dim(), 3);
1197 
1198  for(unsigned int i=0; i<elem.n_nodes(); i++)
1199  {
1200  tmp_nodes[i] = elem.nodes_[permutation_vec[i]];
1201  }
1202  elem.nodes_ = tmp_nodes;
1203 }
1204 
1205 
1206 void Mesh::permute_triangle(unsigned int elm_idx, std::vector<unsigned int> permutation_vec)
1207 {
1208  ASSERT_LT_DBG(elm_idx, element_vec_.size());
1209  ASSERT_EQ_DBG(permutation_vec.size(), 3);
1210 
1211  std::array<unsigned int, 4> tmp_nodes;
1212  Element &elem = element_vec_[elm_idx];
1213  ASSERT_EQ_DBG(elem.dim(), 2);
1214 
1215  for(unsigned int i=0; i<elem.n_nodes(); i++)
1216  {
1217  tmp_nodes[i] = elem.nodes_[permutation_vec[i]];
1218  }
1219  elem.nodes_ = tmp_nodes;
1220 }
1221 
1222 
1224  if (bc_mesh_ == nullptr) bc_mesh_ = new BCMesh(this);
1225  return bc_mesh_;
1226 }
1227 
1228 
1230  ASSERT_PTR(el_4_loc).error("Array 'el_4_loc' is not initialized. Did you call Partitioning::id_maps?\n");
1231 
1232  unsigned int i_proc, i_node, i_ghost_node, elm_node;
1233  unsigned int my_proc = el_ds->myp();
1234  unsigned int n_proc = el_ds->np();
1235 
1236  // distribute nodes between processes, every node is assigned to minimal process of elements that own node
1237  // fill min_node_proc vector with same values on all processes
1238  std::vector<unsigned int> node_proc( this->n_nodes(), n_proc );
1239  std::vector<bool> local_node_flag( this->n_nodes(), false );
1240 
1241  for ( auto elm : this->elements_range() ) {
1242  i_proc = elm.proc();
1243  for (elm_node=0; elm_node<elm->n_nodes(); elm_node++) {
1244  i_node = elm->node_idx(elm_node);
1245  if (i_proc == my_proc) local_node_flag[i_node] = true;
1246  if (i_proc < node_proc[i_node]) node_proc[i_node] = i_proc;
1247  }
1248  }
1249 
1250  unsigned int n_own_nodes=0, n_local_nodes=0; // number of own and ghost nodes
1251  for(uint i_proc : node_proc) if (i_proc == my_proc) n_own_nodes++;
1252  for(uint loc_flag : local_node_flag) if (loc_flag) n_local_nodes++;
1253 
1254  //DebugOut() << print_var(n_own_nodes) << print_var(n_local_nodes) << this->n_nodes();
1255  // create and fill node_4_loc_ (mapping local to global indexes)
1256  node_4_loc_ = new LongIdx [ n_local_nodes ];
1257  i_node=0;
1258  i_ghost_node = n_own_nodes;
1259  for (unsigned int i=0; i<this->n_nodes(); ++i) {
1260  if (local_node_flag[i]) {
1261  if (node_proc[i]==my_proc)
1262  node_4_loc_[i_node++] = i;
1263  else
1264  node_4_loc_[i_ghost_node++] = i;
1265  }
1266  }
1267 
1268  // Construct node distribution object, set number of local nodes (own+ghost)
1269  node_ds_ = new Distribution(n_own_nodes, PETSC_COMM_WORLD);
1270  node_ds_->get_lsizes_array(); // need to initialize lsizes data member
1272 
1273 }
1274 
1275 //-----------------------------------------------------------------------------
1276 // vim: set cindent:
int n_triangles
Definition: mesh.h:275
Distribution * el_ds
Parallel distribution of elements.
Definition: mesh.h:560
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: long_idx.hh:22
Class for the mesh partitioning. This should provide:
Definition: partitioning.hh:52
Iterator< ValueType > begin() const
void read_stream(istream &in, const Type::TypeBase &root_type, FileFormat format)
This method actually reads the given stream in.
unsigned int n_sides() const
Returns number of sides aligned with the edge.
Definition: accessors.hh:298
void output_internal_ngh_data()
Output of neighboring data into raw output.
Definition: mesh.cc:1090
const Element * element() const
Definition: accessors.hh:159
unsigned int n_nodes() const
Definition: elements.h:128
std::array< unsigned int, 4 > nodes_
indices to element&#39;s nodes
Definition: elements.h:113
vector< Element > element_vec_
Definition: mesh.h:512
BidirectionalMap< int > node_ids_
Maps node ids to indexes into vector node_vec_.
Definition: mesh.h:532
vector< vector< unsigned int > > const & node_elements()
Definition: mesh.cc:1022
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:332
Accessor to input data conforming to declared Array.
Definition: accessors.hh:567
unsigned int n_local_nodes_
Hold number of local nodes (own + ghost), value is equal with size of node_4_loc array.
Definition: mesh.h:566
unsigned int * boundary_idx_
Definition: elements.h:82
void permute_tetrahedron(unsigned int elm_idx, std::vector< unsigned int > permutation_vec)
Permute nodes of 3D elements of given elm_idx.
Definition: mesh.cc:1189
bool is_boundary() const
Returns true if it is a Boundary region and false if it is a Bulk region.
Definition: region.hh:74
unsigned int uint
unsigned int bc_ele_idx_
Definition: mesh_data.hh:53
friend class Element
Definition: mesh.h:541
void count_side_types()
Definition: mesh.cc:338
NodeAccessor< 3 > node_accessor(unsigned int ni) const
Definition: accessors.hh:202
int MPI_Comm
Definition: mpi.h:141
BCMesh * bc_mesh_
Boundary mesh, object is created only if it&#39;s necessary.
Definition: mesh.h:568
Definition: nodes.hh:31
unsigned int side_idx() const
Returns local index of the side on the element.
Definition: accessors.hh:413
void check_and_finish()
Definition: mesh.cc:925
MapElementIDToRegionID el_to_reg_map_
Definition: region.hh:571
Reader for (slightly) modified input files.
void init_element_vector(unsigned int size)
Initialize element_vec_, set size and reset counters of boundary and bulk elements.
Definition: mesh.cc:1030
virtual const LongIdx * get_local_part()
Definition: mesh.cc:266
MixedMeshIntersections & mixed_intersections()
Definition: mesh.cc:733
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
BidirectionalMap< int > element_ids_
Maps element ids to indexes into vector element_vec_.
Definition: mesh.h:524
Class for declaration of the input of type Bool.
Definition: type_base.hh:459
Edge edge() const
Returns pointer to the edge connected to the side.
void add_node(unsigned int node_id, arma::vec3 coords)
Add new node of given id and coordinates to mesh.
Definition: mesh.cc:973
LongIdx * node_4_loc_
Index set assigning to local node index its global index.
Definition: mesh.h:562
#define MessageOut()
Macro defining &#39;message&#39; record of log.
Definition: logger.hh:243
DuplicateNodes * tree
Definition: mesh.h:263
unsigned int elem_idx_
Index of element in Mesh::element_vec_.
Definition: neighbours.h:137
int n_lines
Definition: mesh.h:274
static const Input::Type::Record & get_input_type()
Definition: partitioning.cc:49
ElementAccessor< 3 > element_accessor(unsigned int idx) const override
Overwrite Mesh::element_accessor()
Definition: bc_mesh.cc:89
static const unsigned int undef_idx
Definition: mesh.h:102
std::vector< BoundingBox > get_element_boxes()
Compute bounding boxes of elements contained in mesh.
Definition: mesh.cc:940
void init()
Definition: mesh.cc:156
ofstream raw_ngh_output_file
Definition: mesh.h:570
std::string format(CStringRef format_str, ArgList args)
Definition: format.h:3141
NodeAccessor< 3 > node_accessor(unsigned int idx) const
Create and return NodeAccessor to node of given idx.
Definition: mesh.cc:756
#define INPUT_CATCH(ExceptionType, AddressEITag, input_accessor)
Definition: accessors.hh:64
void create_node_element_lists()
Definition: mesh.cc:352
static unsigned int permutation_index(unsigned int p[n_nodes_per_side])
Definition: ref_element.cc:552
BCMesh * get_bc_mesh()
Create boundary mesh if doesn&#39;t exist and return it.
Definition: mesh.cc:1223
static const Input::Type::Record & get_input_type()
Definition: mesh.cc:73
int elem_index(int elem_id) const
For element of given elem_id returns index in element_vec_ or (-1) if element doesn&#39;t exist...
Definition: mesh.h:353
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
std::shared_ptr< MixedMeshIntersections > intersections
Definition: mesh.h:255
double global_snap_radius() const
Maximal distance of observe point from Mesh relative to its size.
Definition: mesh.cc:964
unsigned int max_edge_sides_[3]
Maximal number of sides per one edge in the actual mesh (set in make_neighbours_and_edges()).
Definition: mesh.h:476
Definition: mesh.h:76
unsigned int n_local_nodes() const
Definition: mesh.h:177
std::vector< EdgeData > edges
Vector of MH edges, this should not be part of the geometrical mesh.
Definition: mesh.h:535
Input::Record in_record_
Definition: mesh.h:500
std::string create_label_from_id(unsigned int id) const
Definition: region.cc:338
vector< vector< vector< unsigned int > > > side_nodes
Definition: mesh.h:287
SideIter side(const unsigned int loc_index)
int n_tetrahedras
Definition: mesh.h:276
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:347
virtual unsigned int n_nodes() const
Definition: mesh.h:124
vector< BoundaryData > boundary_
Definition: mesh.h:247
Boundary boundary(uint edge_idx) const
Definition: mesh.cc:256
void add_physical_name(unsigned int dim, unsigned int id, std::string name)
Add new node of given id and coordinates to mesh.
Definition: mesh.cc:968
int node_index(int node_id) const
For node of given node_id returns index in element_vec_ or (-1) if node doesn&#39;t exist.
Definition: mesh.h:365
ElementAccessor< 3 > element() const
Returns iterator to the element of the side.
virtual Partitioning * get_part()
Definition: mesh.cc:262
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
Class for declaration of inputs sequences.
Definition: type_base.hh:346
virtual ElementAccessor< 3 > element_accessor(unsigned int idx) const
Create and return ElementAccessor to element of given idx.
Definition: mesh.cc:750
int n_exsides
Definition: mesh.h:271
void init(unsigned int dim, RegionIdx reg)
Definition: elements.cc:52
IteratorBase end() const
const unsigned int * get_lsizes_array()
get local sizes array
Region get_region(unsigned int id, unsigned int dim)
Definition: region.cc:151
unsigned int dim() const
Definition: elements.h:123
static Default optional()
The factory function to make an empty default value which is optional.
Definition: type_record.hh:124
#define OLD_ASSERT(...)
Definition: global_defs.h:131
unsigned int n_vb_neighbours() const
Definition: mesh.cc:235
ElementAccessor< 3 > element()
Definition: neighbours.h:161
void open_stream(Stream &stream) const
Definition: file_path.cc:211
virtual bool check_compatible_mesh(Mesh &mesh, vector< LongIdx > &bulk_elements_id, vector< LongIdx > &boundary_elements_id)
Definition: mesh.cc:803
bool opt_val(const string &key, Ret &value) const
unsigned int edge_idx_
Index of Edge in Mesh.
Definition: neighbours.h:138
void read_regions_from_input(Input::Array region_list)
Definition: mesh.cc:915
unsigned int boundary_loaded_size_
Count of boundary elements loaded from mesh file.
Definition: mesh.h:521
virtual Record & allow_auto_conversion(const string &from_key)
Allows shorter input of the Record providing only value of the from_key given as the parameter...
Definition: type_record.cc:133
bool same_sides(const SideIter &si, vector< unsigned int > &side_nodes)
Definition: mesh.cc:417
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:541
int find_elem_id(unsigned int pos) const
Return element id (in GMSH file) of element of given position in element vector.
Definition: mesh.h:359
unsigned int add_item(T val)
Add new item at the end position of map.
Mesh * mesh_
Definition: mesh_data.hh:54
unsigned int node_idx(unsigned int ni) const
Return index (in Mesh::node_vec) of ni-th node.
Definition: elements.h:74
void check_element_size(unsigned int elem_idx) const
Check if given index is in element_vec_.
Definition: mesh.cc:1082
void setup_topology()
Definition: mesh.cc:302
IntersectionSearch get_intersection_search()
Getter for input type selection for intersection search algorithm.
Definition: mesh.cc:150
Neighbour ** neigh_vb
Definition: elements.h:86
LongIdx * el_4_loc
Index set assigning to local element index its global index.
Definition: mesh.h:558
Range helper class.
friend class BCMesh
Definition: mesh.h:544
T get_root_interface() const
Returns the root accessor.
unsigned int edge_idx(unsigned int edg_idx) const
Return edge_idx of given index.
Definition: elements.h:138
friend class RegionSetBase
Definition: mesh.h:540
static const Input::Type::Selection & get_input_intersection_variant()
The definition of input record for selection of variant of file format.
Definition: mesh.cc:62
static FileName input()
The factory function for declaring type FileName for input files.
Definition: type_base.cc:526
SideIter side()
Definition: neighbours.h:145
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
unsigned int edge_idx_
Definition: mesh_data.hh:49
#define NDEF
Definition: mesh.cc:58
const Ret val(const string &key) const
unsigned int n_sides() const
Definition: elements.h:134
#define xprintf(...)
Definition: system.hh:93
#define START_TIMER(tag)
Starts a timer with specified tag.
Class for O(log N) lookup for intersections with a set of bounding boxes.
Definition: bih_tree.hh:38
Selection & add_value(const int value, const std::string &key, const std::string &description="", TypeBase::attribute_map attributes=TypeBase::attribute_map())
Adds one new value with name given by key to the Selection.
unsigned int n_corners()
Definition: mesh.cc:240
void set_item(T val, unsigned int pos)
LongIdx * row_4_el
Index set assigning to global element index the local index used in parallel vectors.
Definition: mesh.h:556
BoundingBox bounding_box() const
Definition: accessors.hh:210
Record & declare_key(const string &key, std::shared_ptr< TypeBase > type, const Default &default_value, const string &description, TypeBase::attribute_map key_attributes=TypeBase::attribute_map())
Declares a new key of the Record.
Definition: type_record.cc:503
void close()
Definition: region.cc:250
bool is_empty() const
Definition: accessors.hh:366
vector< Node > node_vec_
Definition: mesh.h:529
Region implicit_boundary_region()
Definition: region.cc:76
virtual Range< ElementAccessor< 3 > > elements_range() const
Returns range of bulk elements.
Definition: mesh.cc:1064
unsigned int np() const
get num of processors
unsigned int n_sides() const
Definition: mesh.cc:226
std::shared_ptr< BIHTree > bih_tree_
Definition: mesh.h:494
Element * add_element_to_vector(int id, bool boundary=false)
Adds element to mesh data structures (element_vec_, element_ids_), returns pointer to this element...
Definition: mesh.cc:1048
vector< vector< unsigned int > > node_elements_
For each node the vector contains a list of elements that use this node.
Definition: mesh.h:350
void reinit(unsigned int init_size=0)
Reset data of map, allow reserve size.
NodeAccessor< 3 > node(unsigned int i) const
Returns node for given local index i on the side.
friend class Boundary
Definition: mesh.h:543
Accessor to the polymorphic input data of a type given by an AbstracRecord object.
Definition: accessors.hh:459
Region region() const
Definition: accessors.hh:164
void elements_id_maps(vector< LongIdx > &bulk_elements_id, vector< LongIdx > &boundary_elements_id) const
Definition: mesh.cc:762
void count_element_types()
Definition: mesh.cc:275
bool compare_points(const arma::vec3 &p1, const arma::vec3 &p2)
Definition: mesh.cc:795
const BIHTree & get_bih_tree()
Getter for BIH. Creates and compute BIH at first call.
Definition: mesh.cc:955
void create_boundary_elements()
Create boundary elements from data of temporary structure, this method MUST be call after read mesh f...
Definition: mesh.cc:1172
#define MPI_Comm_rank
Definition: mpi.h:236
void check_regions()
Definition: region.cc:462
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
unsigned int myp() const
get my processor
Definition: system.hh:65
void mark_used_region(unsigned int idx)
Definition: region.cc:236
Support classes for parallel programing.
Region add_region(unsigned int id, const std::string &label, unsigned int dim, const std::string &address="implicit")
Definition: region.cc:86
vector< Neighbour > vb_neighbours_
Definition: mesh.h:268
Distribution * node_ds_
Parallel distribution of nodes. Depends on elements distribution.
Definition: mesh.h:564
FinishStatus finish(FinishStatus finish_type=FinishStatus::regular_) override
Finish declaration of the Record type.
Definition: type_record.cc:243
virtual unsigned int n_elements(bool boundary=false) const
Returns count of boundary or bulk elements.
Definition: mesh.h:344
void distribute_nodes()
Fill array node_4_loc_ and create object node_ds_ according to element distribution.
Definition: mesh.cc:1229
void find_point(const Space< 3 >::Point &point, std::vector< unsigned int > &result_list, bool full_list=false) const
Definition: bih_tree.cc:287
int n_sides_
Definition: mesh.h:272
RegionIdx region_idx_
Definition: elements.h:109
Class represents boundary part of mesh.
Definition: bc_mesh.hh:35
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
Definition: asserts.hh:336
const Selection & close() const
Close the Selection, no more values can be added.
RegionDB region_db_
Definition: mesh.h:485
MPI_Comm comm_
Definition: mesh.h:505
#define ASSERT_DBG(expr)
bool find_lower_dim_element(vector< unsigned int > &element_list, unsigned int dim, unsigned int &element_idx)
Definition: mesh.cc:397
int pid_
Id # of mesh partition.
Definition: elements.h:91
int n_insides
Definition: mesh.h:270
void print_region_table(std::ostream &stream) const
Definition: region.cc:410
void intersect_element_lists(vector< unsigned int > const &nodes_list, vector< unsigned int > &intersection_element_list)
Definition: mesh.cc:365
void modify_element_ids(const RegionDB::MapElementIDToRegionID &map)
Definition: mesh.cc:292
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:246
#define MPI_COMM_WORLD
Definition: mpi.h:123
void make_edge_permutations()
Definition: mesh.cc:627
double tetrahedron_jacobian() const
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc...
static Input::Type::Abstract & get_input_type()
Definition: region_set.cc:25
Range< NodeAccessor< 3 > > node_range() const
Returns range of nodes.
Definition: mesh.cc:1070
SideIter * side_
Definition: mesh_data.hh:32
Edge edge(uint edge_idx) const
Definition: mesh.cc:250
const LongIdx * get_loc_part() const
Definition: partitioning.cc:85
void add_element(unsigned int elm_id, unsigned int dim, unsigned int region_id, unsigned int partition_id, std::vector< unsigned int > node_ids)
Add new element of given id to mesh.
Definition: mesh.cc:981
Record type proxy class.
Definition: type_record.hh:182
bool is_valid() const
Returns false if the region has undefined/invalid value.
Definition: region.hh:78
IntersectionSearch
Types of search algorithm for finding intersection candidates.
Definition: mesh.h:91
void element_to_neigh_vb()
Definition: mesh.cc:700
void init_element(Element *ele, unsigned int elm_id, unsigned int dim, RegionIdx region_idx, unsigned int partition_id, std::vector< unsigned int > node_ids)
Initialize element.
Definition: mesh.cc:1003
friend class Edge
Definition: mesh.h:538
Mesh()
Definition: mesh.cc:99
arma::vec3 & point()
Definition: nodes.hh:67
static FileName output()
The factory function for declaring type FileName for input files.
Definition: type_base.cc:533
unsigned int n_neighs_vb_
of neighbours, V-B type (comp.)
Definition: elements.h:93
void permute_triangle(unsigned int elm_idx, std::vector< unsigned int > permutation_vec)
Permute nodes of 2D elements of given elm_idx.
Definition: mesh.cc:1206
#define ASSERT_PTR_DBG(ptr)
Definition of assert macro checking non-null pointer (PTR) only for debug mode.
Definition: asserts.hh:340
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
Definition: accessors.hh:180
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
Definition: exceptions.hh:53
SideIter side(const unsigned int i) const
Gets side iterator of the i -th side.
Range< Edge > edge_range() const
Returns range of edges.
Definition: mesh.cc:1076
Template for classes storing finite set of named values.
Implementation of range helper class.
unsigned int bulk_size_
Count of bulk elements.
Definition: mesh.h:518
std::shared_ptr< Partitioning > part_
Definition: mesh.h:489
void make_neighbours_and_edges()
Definition: mesh.cc:434
vector< ElementTmpData > bc_element_tmp_
Hold data of boundary elements during reading mesh (allow to preserve correct order during reading of...
Definition: mesh.h:515
Mesh * mesh_
Pointer to Mesh to which belonged.
Definition: neighbours.h:136
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual)
Definition: asserts.hh:328
void init_node_vector(unsigned int size)
Initialize node_vec_, set size.
Definition: mesh.cc:1041
unsigned int id() const
Returns id of the region (using RegionDB)
Definition: region.cc:38
Main class for computation of intersection of meshes of combined dimensions.
unsigned int n_nodes() const
Returns number of nodes of the side.
Definition: accessors.hh:405
virtual ~Mesh()
Destructor.
Definition: mesh.cc:201
unsigned int idx() const
Returns a global index of the region.
Definition: region.hh:82
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:300
const Node * node(unsigned int ni) const
Definition: accessors.hh:198
unsigned int n_sides
Definition: mesh_data.hh:29