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