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