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