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