Flow123d  release_3.0.0-916-g95df358
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 (bc_mesh_ != nullptr) delete bc_mesh_;
222  if (tree != nullptr) delete tree;
223 }
224 
225 
226 unsigned int Mesh::n_sides() const
227 {
228  if (n_sides_ == NDEF) {
229  n_sides_=0;
230  for (auto ele : this->elements_range()) n_sides_ += ele->n_sides();
231  }
232  return n_sides_;
233 }
234 
235 unsigned int Mesh::n_vb_neighbours() const {
236  return vb_neighbours_.size();
237  }
238 
239 
240 unsigned int Mesh::n_corners() {
241  unsigned int li, count = 0;
242  for (auto ele : this->elements_range()) {
243  for (li=0; li<ele->n_nodes(); li++) {
244  count++;
245  }
246  }
247  return count;
248 }
249 
251  return part_.get();
252 }
253 
255  return (LongIdx*)this->get_part()->get_loc_part();
256 }
257 
258 
259 //=============================================================================
260 // COUNT ELEMENT TYPES
261 //=============================================================================
262 
264  for (auto elm : this->elements_range())
265  switch (elm->dim()) {
266  case 1:
267  n_lines++;
268  break;
269  case 2:
270  n_triangles++;
271  break;
272  case 3:
273  n_tetrahedras++;
274  break;
275  }
276 }
277 
278 
279 
281  for (auto elem_to_region : map) {
282  Element &ele = element_vec_[ elem_index(elem_to_region.first) ];
283  ele.region_idx_ = region_db_.get_region( elem_to_region.second, ele.dim() );
285  }
286 }
287 
288 
289 
291  START_TIMER("MESH - setup topology");
292 
294 
295  // check mesh quality
296  for (auto ele : this->elements_range())
297  if (ele.quality_measure_smooth(ele.side(0)) < 0.001) WarningOut().fmt("Bad quality (<0.001) of the element {}.\n", ele.idx());
298 
303 
304  tree = new DuplicateNodes(this);
305 
306  part_ = std::make_shared<Partitioning>(this, in_record_.val<Input::Record>("partitioning") );
307 
308  // create parallel distribution and numbering of elements
309  LongIdx *id_4_old = new LongIdx[n_elements()];
310  int i = 0;
311  for (auto ele : this->elements_range())
312  id_4_old[i++] = ele.idx();
313  part_->id_maps(n_elements(), id_4_old, el_ds, el_4_loc, row_4_el);
314 
315  delete[] id_4_old;
316 
317  this->distribute_nodes();
318 
320 }
321 
322 
323 //
325 {
326 
327  n_insides = 0;
328  n_exsides = 0;
329  for (auto ele : this->elements_range())
330  for(SideIter sde = ele.side(0); sde->side_idx() < ele->n_sides(); ++sde) {
331  if (sde->is_external()) n_exsides++;
332  else n_insides++;
333  }
334 }
335 
336 
337 
339  // for each node we make a list of elements that use this node
340  node_elements_.resize( this->n_nodes() );
341 
342  for (auto ele : this->elements_range())
343  for (unsigned int n=0; n<ele->n_nodes(); n++)
344  node_elements_[ele.node_accessor(n).idx()].push_back(ele.idx());
345 
346  for (vector<vector<unsigned int> >::iterator n=node_elements_.begin(); n!=node_elements_.end(); n++)
347  stable_sort(n->begin(), n->end());
348 }
349 
350 
351 void Mesh::intersect_element_lists(vector<unsigned int> const &nodes_list, vector<unsigned int> &intersection_element_list)
352 {
353  if (node_elements_.size() == 0) {
355  }
356 
357  if (nodes_list.size() == 0) {
358  intersection_element_list.clear();
359  } else if (nodes_list.size() == 1) {
360  intersection_element_list = node_elements_[ nodes_list[0] ];
361  } else {
362  vector<unsigned int>::const_iterator it1=nodes_list.begin();
364  intersection_element_list.resize( node_elements_[*it1].size() ); // make enough space
365 
366  it1=set_intersection(
367  node_elements_[*it1].begin(), node_elements_[*it1].end(),
368  node_elements_[*it2].begin(), node_elements_[*it2].end(),
369  intersection_element_list.begin());
370  intersection_element_list.resize(it1-intersection_element_list.begin()); // resize to true size
371 
372  for(;it2<nodes_list.end();++it2) {
373  it1=set_intersection(
374  intersection_element_list.begin(), intersection_element_list.end(),
375  node_elements_[*it2].begin(), node_elements_[*it2].end(),
376  intersection_element_list.begin());
377  intersection_element_list.resize(it1-intersection_element_list.begin()); // resize to true size
378  }
379  }
380 }
381 
382 
383 bool Mesh::find_lower_dim_element( vector<unsigned int> &element_list, unsigned int dim, unsigned int &element_idx) {
384  bool is_neighbour = false;
385 
386  vector<unsigned int>::iterator e_dest=element_list.begin();
387  for( vector<unsigned int>::iterator ele = element_list.begin(); ele!=element_list.end(); ++ele) {
388  if (element_vec_[*ele].dim() == dim) { // keep only indexes of elements of same dimension
389  *e_dest=*ele;
390  ++e_dest;
391  } else if (element_vec_[*ele].dim() == dim-1) { // get only first element of lower dimension
392  if (is_neighbour) xprintf(UsrErr, "Too matching elements id: %d and id: %d in the same mesh.\n",
393  this->elem_index(*ele), this->elem_index(element_idx) );
394 
395  is_neighbour = true;
396  element_idx = *ele;
397  }
398  }
399  element_list.resize( e_dest - element_list.begin());
400  return is_neighbour;
401 }
402 
404  // 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)
405  unsigned int ni=0;
406  while ( ni < si->n_nodes()
407  && find(side_nodes.begin(), side_nodes.end(), si->node(ni).idx() ) != side_nodes.end() ) ni++;
408  return ( ni == si->n_nodes() );
409 }
410 
411 /**
412  * TODO:
413  * - use std::is_any for setting is_neigbour
414  * - possibly make appropriate constructors for Edge and Neighbour
415  * - check side!=-1 when searching neigbouring element
416  * - process boundary elements first, there should be no Neigh, but check it
417  * set Edge and boundary there
418  */
419 
421 {
422  ASSERT(bc_element_tmp_.size()==0)
423  .error("Temporary structure of boundary element data is not empty. Did you call create_boundary_elements?");
424 
425  Neighbour neighbour;
426  Edge *edg;
427  unsigned int ngh_element_idx, last_edge_idx;
428 
429  neighbour.mesh_ = this;
430 
432 
433  // pointers to created edges
434  //vector<Edge *> tmp_edges;
435  edges.resize(0); // be sure that edges are empty
436 
438  vector<unsigned int> intersection_list; // list of elements in intersection of node element lists
439 
440  for( unsigned int i=bulk_size_; i<element_vec_.size(); ++i) {
441  ElementAccessor<3> bc_ele = this->element_accessor(i);
442  // Find all elements that share this side.
443  side_nodes.resize(bc_ele->n_nodes());
444  for (unsigned n=0; n<bc_ele->n_nodes(); n++) side_nodes[n] = bc_ele->node_idx(n);
445  intersect_element_lists(side_nodes, intersection_list);
446  bool is_neighbour = find_lower_dim_element(intersection_list, bc_ele->dim() +1, ngh_element_idx);
447  if (is_neighbour) {
448  xprintf(UsrErr, "Boundary element (id: %d) match a regular element (id: %d) of lower dimension.\n",
449  bc_ele.idx(), this->elem_index(ngh_element_idx));
450  } else {
451  if (intersection_list.size() == 0) {
452  // no matching dim+1 element found
453  WarningOut().fmt("Lonely boundary element, id: {}, region: {}, dimension {}.\n",
454  bc_ele.idx(), bc_ele.region().id(), bc_ele->dim());
455  continue; // skip the boundary element
456  }
457  last_edge_idx=edges.size();
458  edges.resize(last_edge_idx+1);
459  edg = &( edges.back() );
460  edg->n_sides = 0;
461  edg->side_ = new struct SideIter[ intersection_list.size() ];
462 
463  // common boundary object
464  unsigned int bdr_idx=boundary_.size();
465  boundary_.resize(bdr_idx+1);
466  Boundary &bdr=boundary_.back();
467  bdr.bc_ele_idx_ = i;
468  bdr.edge_idx_ = last_edge_idx;
469  bdr.mesh_=this;
470 
471  // for 1d boundaries there can be more then one 1d elements connected to the boundary element
472  // we do not detect this case later in the main search over bulk elements
473  for( vector<unsigned int>::iterator isect = intersection_list.begin(); isect!=intersection_list.end(); ++isect) {
474  ElementAccessor<3> elem = this->element_accessor(*isect);
475  for (unsigned int ecs=0; ecs<elem->n_sides(); ecs++) {
476  SideIter si = elem.side(ecs);
477  if ( same_sides( si, side_nodes) ) {
478  if (elem->edge_idx(ecs) != Mesh::undef_idx) {
479  OLD_ASSERT(elem->boundary_idx_!=nullptr, "Null boundary idx array.\n");
480  int last_bc_ele_idx=this->boundary_[elem->boundary_idx_[ecs]].bc_ele_idx_;
481  int new_bc_ele_idx=i;
482  THROW( ExcDuplicateBoundary()
483  << EI_ElemLast(this->elem_index(last_bc_ele_idx))
484  << EI_RegLast(this->element_accessor(last_bc_ele_idx).region().label())
485  << EI_ElemNew(this->elem_index(new_bc_ele_idx))
486  << EI_RegNew(this->element_accessor(new_bc_ele_idx).region().label())
487  );
488  }
489  element_vec_[*isect].edge_idx_[ecs] = last_edge_idx;
490  edg->side_[ edg->n_sides++ ] = si;
491 
492  if (elem->boundary_idx_ == NULL) {
493  Element *el = &(element_vec_[*isect]);
494  el->boundary_idx_ = new unsigned int [ el->n_sides() ];
495  std::fill( el->boundary_idx_, el->boundary_idx_ + el->n_sides(), Mesh::undef_idx);
496  }
497  elem->boundary_idx_[ecs] = bdr_idx;
498  break; // next element in intersection list
499  }
500  }
501  }
502 
503  }
504 
505  }
506  // Now we go through all element sides and create edges and neighbours
507  for (auto e : this->elements_range()) {
508  for (unsigned int s=0; s<e->n_sides(); s++)
509  {
510  // skip sides that were already found
511  if (e->edge_idx(s) != Mesh::undef_idx) continue;
512 
513 
514  // Find all elements that share this side.
515  side_nodes.resize(e.side(s)->n_nodes());
516  for (unsigned n=0; n<e.side(s)->n_nodes(); n++) side_nodes[n] = e.side(s)->node(n).idx();
517  intersect_element_lists(side_nodes, intersection_list);
518 
519  bool is_neighbour = find_lower_dim_element(intersection_list, e->dim(), ngh_element_idx);
520 
521  if (is_neighbour) { // edge connects elements of different dimensions
522  neighbour.elem_idx_ = ngh_element_idx;
523  } else { // edge connects only elements of the same dimension
524  // Allocate the array of sides.
525  last_edge_idx=edges.size();
526  edges.resize(last_edge_idx+1);
527  edg = &( edges.back() );
528  edg->n_sides = 0;
529  edg->side_ = new struct SideIter[ intersection_list.size() ];
530  if (intersection_list.size() > max_edge_sides_[e->dim()-1])
531  max_edge_sides_[e->dim()-1] = intersection_list.size();
532 
533  if (intersection_list.size() == 1) { // outer edge, create boundary object as well
534  Element &elm = element_vec_[e.idx()];
535  edg->n_sides=1;
536  edg->side_[0] = e.side(s);
537  element_vec_[e.idx()].edge_idx_[s] = last_edge_idx;
538 
539  if (e->boundary_idx_ == NULL) {
540  elm.boundary_idx_ = new unsigned int [ e->n_sides() ];
541  std::fill( elm.boundary_idx_, elm.boundary_idx_ + e->n_sides(), Mesh::undef_idx);
542  }
543 
544  unsigned int bdr_idx=boundary_.size()+1; // need for VTK mesh that has no boundary elements
545  // and bulk elements are indexed from 0
546  boundary_.resize(bdr_idx+1);
547  Boundary &bdr=boundary_.back();
548  elm.boundary_idx_[s] = bdr_idx;
549 
550  // fill boundary element
551  Element * bc_ele = add_element_to_vector(-bdr_idx, true);
552  bc_ele->init(e->dim()-1, region_db_.implicit_boundary_region() );
554  for(unsigned int ni = 0; ni< side_nodes.size(); ni++) bc_ele->nodes_[ni] = side_nodes[ni];
555 
556  // fill Boundary object
557  bdr.edge_idx_ = last_edge_idx;
558  bdr.bc_ele_idx_ = elem_index(-bdr_idx);
559  bdr.mesh_=this;
560 
561  continue; // next side of element e
562  }
563  }
564 
565  // go through the elements connected to the edge or neighbour
566  for( vector<unsigned int>::iterator isect = intersection_list.begin(); isect!=intersection_list.end(); ++isect) {
567  ElementAccessor<3> elem = this->element_accessor(*isect);
568  for (unsigned int ecs=0; ecs<elem->n_sides(); ecs++) {
569  if (elem->edge_idx(ecs) != Mesh::undef_idx) continue;
570  SideIter si = elem.side(ecs);
571  if ( same_sides( si, side_nodes) ) {
572  if (is_neighbour) {
573  // create a new edge and neighbour for this side, and element to the edge
574  last_edge_idx=edges.size();
575  edges.resize(last_edge_idx+1);
576  edg = &( edges.back() );
577  edg->n_sides = 1;
578  edg->side_ = new struct SideIter[1];
579  edg->side_[0] = si;
580  element_vec_[elem.idx()].edge_idx_[ecs] = last_edge_idx;
581 
582  neighbour.edge_idx_ = last_edge_idx;
583 
584  vb_neighbours_.push_back(neighbour); // copy neighbour with this edge setting
585  } else {
586  // connect the side to the edge, and side to the edge
587  edg->side_[ edg->n_sides++ ] = si;
588  element_vec_[elem.idx()].edge_idx_[ecs] = last_edge_idx;
589  }
590  break; // next element from intersection list
591  }
592  } // search for side of other connected element
593  } // connected elements
594  OLD_ASSERT( is_neighbour || ( (unsigned int) edg->n_sides ) == intersection_list.size(), "Some connected sides were not found.\n");
595  } // for element sides
596  } // for elements
597 
598  MessageOut().fmt( "Created {} edges and {} neighbours.\n", edges.size(), vb_neighbours_.size() );
599 }
600 
601 
602 
604 {
605  for (std::vector<Edge>::iterator edg=edges.begin(); edg!=edges.end(); edg++)
606  {
607  // side 0 is reference, so its permutation is 0
608  edg->side(0)->element()->permutation_idx_[edg->side(0)->side_idx()] = 0;
609 
610  if (edg->n_sides > 1)
611  {
612  map<unsigned int,unsigned int> node_numbers;
613  unsigned int permutation[edg->side(0)->n_nodes()];
614 
615  for (unsigned int i=0; i<edg->side(0)->n_nodes(); i++)
616  node_numbers[edg->side(0)->node(i).idx()] = i;
617  //node_numbers[edg->side(0)->node(i).node()] = i;
618 
619  for (int sid=1; sid<edg->n_sides; sid++)
620  {
621  for (unsigned int i=0; i<edg->side(0)->n_nodes(); i++)
622  permutation[node_numbers[edg->side(sid)->node(i).idx()]] = i;
623  //permutation[node_numbers[edg->side(sid)->node(i).node()]] = i;
624 
625  switch (edg->side(0)->dim())
626  {
627  case 0:
628  edg->side(sid)->element()->permutation_idx_[edg->side(sid)->side_idx()] = RefElement<1>::permutation_index(permutation);
629  break;
630  case 1:
631  edg->side(sid)->element()->permutation_idx_[edg->side(sid)->side_idx()] = RefElement<2>::permutation_index(permutation);
632  break;
633  case 2:
634  edg->side(sid)->element()->permutation_idx_[edg->side(sid)->side_idx()] = RefElement<3>::permutation_index(permutation);
635  break;
636  }
637  }
638  }
639  }
640 
641  for (vector<Neighbour>::iterator nb=vb_neighbours_.begin(); nb!=vb_neighbours_.end(); nb++)
642  {
643  map<const Node*,unsigned int> node_numbers;
644  unsigned int permutation[nb->element()->n_nodes()];
645 
646  // element of lower dimension is reference, so
647  // we calculate permutation for the adjacent side
648  for (unsigned int i=0; i<nb->element()->n_nodes(); i++)
649  node_numbers[nb->element().node(i)] = i;
650 
651  for (unsigned int i=0; i<nb->side()->n_nodes(); i++)
652  permutation[node_numbers[nb->side()->node(i).node()]] = i;
653 
654  switch (nb->side()->dim())
655  {
656  case 0:
657  nb->side()->element()->permutation_idx_[nb->side()->side_idx()] = RefElement<1>::permutation_index(permutation);
658  break;
659  case 1:
660  nb->side()->element()->permutation_idx_[nb->side()->side_idx()] = RefElement<2>::permutation_index(permutation);
661  break;
662  case 2:
663  nb->side()->element()->permutation_idx_[nb->side()->side_idx()] = RefElement<3>::permutation_index(permutation);
664  break;
665  }
666  }
667 }
668 
669 
670 
671 
672 
673 //=============================================================================
674 //
675 //=============================================================================
677 {
678 
679  //MessageOut() << "Element to neighbours of vb2 type... "/*orig verb 5*/;
680 
681  for (vector<Element>::iterator ele = element_vec_.begin(); ele!= element_vec_.begin()+bulk_size_; ++ele)
682  ele->n_neighs_vb_ =0;
683 
684  // count vb neighs per element
685  for (auto & ngh : this->vb_neighbours_) ngh.element()->n_neighs_vb_++;
686 
687  // Allocation of the array per element
688  for (vector<Element>::iterator ele = element_vec_.begin(); ele!= element_vec_.begin()+bulk_size_; ++ele)
689  if( ele->n_neighs_vb() > 0 ) {
690  ele->neigh_vb = new struct Neighbour* [ele->n_neighs_vb()];
691  ele->n_neighs_vb_=0;
692  }
693 
694  // fill
695  ElementAccessor<3> ele;
696  for (auto & ngh : this->vb_neighbours_) {
697  ele = ngh.element();
698  ele->neigh_vb[ ele->n_neighs_vb_++ ] = &ngh;
699  }
700 
701  //MessageOut() << "... O.K.\n"/*orig verb 6*/;
702 }
703 
704 
705 
706 
707 
708 
710  /* Algorithm:
711  *
712  * 1) create BIH tree
713  * 2) for every 1D, find list of candidates
714  * 3) compute intersections for 1d, store it to master_elements
715  *
716  */
717  if (! intersections) {
718  intersections = std::make_shared<MixedMeshIntersections>(this);
719  intersections->compute_intersections();
720  }
721  return *intersections;
722 }
723 
724 
725 
727  return ElementAccessor<3>(this, idx);
728 }
729 
730 
731 
732 NodeAccessor<3> Mesh::node_accessor(unsigned int idx) const {
733  return NodeAccessor<3>(this, idx);
734 }
735 
736 
737 
738 void Mesh::elements_id_maps( vector<LongIdx> & bulk_elements_id, vector<LongIdx> & boundary_elements_id) const
739 {
740  if (bulk_elements_id.size() ==0) {
742  LongIdx last_id;
743 
744  bulk_elements_id.resize(n_elements());
745  map_it = bulk_elements_id.begin();
746  last_id = -1;
747  for(unsigned int idx=0; idx < n_elements(); idx++, ++map_it) {
748  LongIdx id = this->find_elem_id(idx);
749  if (last_id >= id) xprintf(UsrErr, "Element IDs in non-increasing order, ID: %d\n", id);
750  last_id=*map_it = id;
751  }
752 
753  boundary_elements_id.resize(n_elements(true));
754  map_it = boundary_elements_id.begin();
755  last_id = -1;
756  for(unsigned int idx=bulk_size_; idx<element_vec_.size(); idx++, ++map_it) {
757  LongIdx id = this->find_elem_id(idx);
758  // We set ID for boundary elements created by the mesh itself to "-1"
759  // this force gmsh reader to skip all remaining entries in boundary_elements_id
760  // and thus report error for any remaining data lines
761  if (id < 0) last_id=*map_it=-1;
762  else {
763  if (last_id >= id) xprintf(UsrErr, "Element IDs in non-increasing order, ID: %d\n", id);
764  last_id=*map_it = id;
765  }
766  }
767  }
768 }
769 
770 
771 bool compare_points(const arma::vec3 &p1, const arma::vec3 &p2) {
772  static const double point_tolerance = 1E-10;
773  return fabs(p1[0]-p2[0]) < point_tolerance
774  && fabs(p1[1]-p2[1]) < point_tolerance
775  && fabs(p1[2]-p2[2]) < point_tolerance;
776 }
777 
778 
779 bool Mesh::check_compatible_mesh( Mesh & mesh, vector<LongIdx> & bulk_elements_id, vector<LongIdx> & boundary_elements_id )
780 {
781  std::vector<unsigned int> node_ids; // allow mapping ids of nodes from source mesh to target mesh
782  std::vector<unsigned int> node_list;
783  std::vector<unsigned int> candidate_list; // returned by intersect_element_lists
784  std::vector<unsigned int> result_list; // list of elements with same dimension as vtk element
785  unsigned int i; // counter over vectors
786 
787  {
788  // iterates over node vector of \p this object
789  // to each node must be found just only one node in target \p mesh
790  // store orders (mapping between source and target meshes) into node_ids vector
791  std::vector<unsigned int> searched_elements; // for BIH tree
792  unsigned int i_node, i_elm_node;
793  const BIHTree &bih_tree=mesh.get_bih_tree();
794 
795  // create nodes of mesh
796  node_ids.resize( this->n_nodes() );
797  i=0;
798  for (auto nod : this->node_range()) {
799  arma::vec3 point = nod->point();
800  int found_i_node = -1;
801  bih_tree.find_point(point, searched_elements);
802 
803  for (std::vector<unsigned int>::iterator it = searched_elements.begin(); it!=searched_elements.end(); it++) {
804  ElementAccessor<3> ele = mesh.element_accessor( *it );
805  for (i_node=0; i_node<ele->n_nodes(); i_node++)
806  {
807  if ( compare_points(ele.node(i_node)->point(), point) ) {
808  i_elm_node = ele.node_accessor(i_node).idx();
809  if (found_i_node == -1) found_i_node = i_elm_node;
810  else if (found_i_node != i_elm_node) {
811  // duplicate nodes in target mesh
812  this->elements_id_maps(bulk_elements_id, boundary_elements_id);
813  return false;
814  }
815  }
816  }
817  }
818  if (found_i_node == -1) {
819  // no node found in target mesh
820  this->elements_id_maps(bulk_elements_id, boundary_elements_id);
821  return false;
822  }
823  node_ids[i] = (unsigned int)found_i_node;
824  searched_elements.clear();
825  i++;
826  }
827  }
828 
829  {
830  // iterates over bulk elements of \p this object
831  // elements in both meshes must be in ratio 1:1
832  // store orders (mapping between both mesh files) into bulk_elements_id vector
833  bulk_elements_id.clear();
834  bulk_elements_id.resize(this->n_elements());
835  // iterate trough bulk part of element vector, to each element in source mesh must exist only one element in target mesh
836  // fill bulk_elements_id vector
837  i=0;
838  for (auto elm : this->elements_range()) {
839  for (unsigned int j=0; j<elm->n_nodes(); j++) { // iterate trough all nodes of any element
840  node_list.push_back( node_ids[ elm->node_idx(j) ] );
841  }
842  mesh.intersect_element_lists(node_list, candidate_list);
843  for (auto i_elm : candidate_list) {
844  if ( mesh.element_accessor(i_elm)->dim() == elm.dim() ) result_list.push_back( elm.index() );
845  }
846  if (result_list.size() != 1) {
847  // intersect_element_lists must produce one element
848  this->elements_id_maps(bulk_elements_id, boundary_elements_id);
849  return false;
850  }
851  bulk_elements_id[i] = (LongIdx)result_list[0];
852  node_list.clear();
853  result_list.clear();
854  i++;
855  }
856  }
857 
858  {
859  // iterates over boundary elements of \p this object
860  // elements in both meshes must be in ratio 1:1
861  // store orders (mapping between both mesh files) into boundary_elements_id vector
862  auto bc_mesh = this->get_bc_mesh();
863  boundary_elements_id.clear();
864  boundary_elements_id.resize(bc_mesh->n_elements());
865  // iterate trough boundary part of element vector, to each element in source mesh must exist only one element in target mesh
866  // fill boundary_elements_id vector
867  i=0;
868  for (auto elm : bc_mesh->elements_range()) {
869  for (unsigned int j=0; j<elm->n_nodes(); j++) { // iterate trough all nodes of any element
870  node_list.push_back( node_ids[ elm->node_idx(j) ] );
871  }
872  mesh.get_bc_mesh()->intersect_element_lists(node_list, candidate_list);
873  for (auto i_elm : candidate_list) {
874  if ( mesh.get_bc_mesh()->element_accessor(i_elm)->dim() == elm.dim() ) result_list.push_back( elm.index() );
875  }
876  if (result_list.size() != 1) {
877  // intersect_element_lists must produce one element
878  this->elements_id_maps(bulk_elements_id, boundary_elements_id);
879  return false;
880  }
881  boundary_elements_id[i] = (LongIdx)result_list[0];
882  node_list.clear();
883  result_list.clear();
884  i++;
885  }
886  }
887 
888  return true;
889 }
890 
892 {
894  it != region_list.end();
895  ++it) {
896  // constructor has side effect in the mesh - create new region or set and store them to Mesh::region_db_
897  (*it).factory< RegionSetBase, const Input::Record &, Mesh * >(*it, this);
898  }
899 }
900 
902 {
904  region_db_.el_to_reg_map_.clear();
905  region_db_.close();
907 
908  if ( in_record_.val<bool>("print_regions") ) {
909  stringstream ss;
911  MessageOut() << ss.str();
912  }
913 }
914 
915 
917  START_TIMER("Mesh::compute_element_boxes");
919 
920  // make element boxes
921  unsigned int i=0;
922  boxes.resize(this->n_elements());
923  for (auto element : this->elements_range()) {
924  boxes[i] = element.bounding_box();
925  i++;
926  }
927 
928  return boxes;
929 }
930 
932  if (! this->bih_tree_) {
933  bih_tree_ = std::make_shared<BIHTree>();
934  bih_tree_->add_boxes( this->get_element_boxes() );
935  bih_tree_->construct();
936  }
937  return *bih_tree_;
938 }
939 
940 double Mesh::global_snap_radius() const {
941  return in_record_.val<double>("global_snap_radius");
942 }
943 
944 void Mesh::add_physical_name(unsigned int dim, unsigned int id, std::string name) {
945  region_db_.add_region(id, name, dim, "$PhysicalNames");
946 }
947 
948 
949 void Mesh::add_node(unsigned int node_id, arma::vec3 coords) {
950  node_vec_.push_back( Node() );
951  Node &node = node_vec_[node_vec_.size()-1];
952  node.point() = coords;
953  node_ids_.add_item(node_id);
954 }
955 
956 
957 void Mesh::add_element(unsigned int elm_id, unsigned int dim, unsigned int region_id, unsigned int partition_id,
958  std::vector<unsigned int> node_ids) {
959  RegionIdx region_idx = region_db_.get_region( region_id, dim );
960  if ( !region_idx.is_valid() ) {
961  region_idx = region_db_.add_region( region_id, region_db_.create_label_from_id(region_id), dim, "$Element" );
962  }
963  region_db_.mark_used_region(region_idx.idx());
964 
965  if (region_idx.is_boundary()) {
966  bc_element_tmp_.push_back( ElementTmpData(elm_id, dim, region_idx, partition_id, node_ids) );
967  } else {
968  if(dim == 0 ) {
969  WarningOut().fmt("Bulk elements of zero size(dim=0) are not supported. Element ID: {}.\n", elm_id);
970  }
971  else {
972  Element *ele = add_element_to_vector(elm_id);
973  this->init_element(ele, elm_id, dim, region_idx, partition_id, node_ids);
974  }
975  }
976 }
977 
978 
979 void Mesh::init_element(Element *ele, unsigned int elm_id, unsigned int dim, RegionIdx region_idx, unsigned int partition_id,
980  std::vector<unsigned int> node_ids) {
981  ele->init(dim, region_idx);
982  ele->pid_ = partition_id;
983 
984  for (unsigned int ni=0; ni<ele->n_nodes(); ni++) {
985  ele->nodes_[ni] = this->node_index(node_ids[ni]);
986  }
987 
988  // check that tetrahedron element is numbered correctly and is not degenerated
989  if(ele->dim() == 3)
990  {
991  double jac = this->element_accessor( this->elem_index(elm_id) ).tetrahedron_jacobian();
992  if( ! (jac > 0) )
993  WarningOut().fmt("Tetrahedron element with id {} has wrong numbering or is degenerated (Jacobian = {}).",elm_id,jac);
994  }
995 }
996 
997 
999  if (node_elements_.size() == 0) {
1000  this->create_node_element_lists();
1001  }
1002  return node_elements_;
1003 }
1004 
1005 
1006 void Mesh::init_element_vector(unsigned int size) {
1007  element_vec_.clear();
1008  element_vec_.resize(size);
1009  element_ids_.reinit(size);
1010  bc_element_tmp_.clear();
1011  bc_element_tmp_.reserve(size);
1012  bulk_size_ = 0;
1014 }
1015 
1016 
1017 void Mesh::init_node_vector(unsigned int size) {
1018  node_vec_.clear();
1019  node_vec_.reserve(size);
1020  node_ids_.reinit(0);
1021 }
1022 
1023 
1024 Element * Mesh::add_element_to_vector(int id, bool boundary) {
1025  Element * elem=nullptr;
1026  if (boundary) {
1027  ASSERT_DBG(id<0)(id).error("Add boundary element from mesh file trough temporary structure!");
1028  element_vec_.push_back( Element() );
1029  elem = &element_vec_[element_vec_.size()-1];
1030  element_ids_.add_item(id);
1031  } else {
1032  elem = &element_vec_[bulk_size_];
1034  bulk_size_++;
1035  }
1036 
1037  return elem;
1038 }
1039 
1041  auto bgn_it = make_iter<ElementAccessor<3>>( ElementAccessor<3>(this, 0) );
1042  auto end_it = make_iter<ElementAccessor<3>>( ElementAccessor<3>(this, bulk_size_) );
1043  return Range<ElementAccessor<3>>(bgn_it, end_it);
1044 }
1045 
1047  auto bgn_it = make_iter<NodeAccessor<3>>( NodeAccessor<3>(this, 0) );
1048  auto end_it = make_iter<NodeAccessor<3>>( NodeAccessor<3>(this, node_vec_.size()) );
1049  return Range<NodeAccessor<3>>(bgn_it, end_it);
1050 }
1051 
1052 inline void Mesh::check_element_size(unsigned int elem_idx) const
1053 {
1054  ASSERT(elem_idx < element_vec_.size())(elem_idx)(element_vec_.size()).error("Index of element is out of bound of element vector!");
1055 }
1056 
1057 /*
1058  * Output of internal flow data.
1059  */
1061 {
1062  START_TIMER("Mesh::output_internal_ngh_data");
1063 
1064  if (! raw_ngh_output_file.is_open()) return;
1065 
1066  // header
1067  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";
1069 
1070  int cit = 0;
1071 
1072  // map from higher dim elements to its lower dim neighbors, using gmsh IDs: ele->id()
1073  unsigned int undefined_ele_id = -1;
1075  for (auto ele : this->elements_range()) {
1076  if(ele->n_neighs_vb() > 0){
1077  for (unsigned int i = 0; i < ele->n_neighs_vb(); i++){
1078  ElementAccessor<3> higher_ele = ele->neigh_vb[i]->side()->element();
1079 
1080  auto search = neigh_vb_map.find(higher_ele.idx());
1081  if(search != neigh_vb_map.end()){
1082  // if found, add id to correct local side idx
1083  search->second[ele->neigh_vb[i]->side()->side_idx()] = ele.idx();
1084  }
1085  else{
1086  // if not found, create new vector, each side can have one vb neighbour
1087  std::vector<unsigned int> higher_ele_side_ngh(higher_ele->n_sides(), undefined_ele_id);
1088  higher_ele_side_ngh[ele->neigh_vb[i]->side()->side_idx()] = ele.idx();
1089  neigh_vb_map[higher_ele.idx()] = higher_ele_side_ngh;
1090  }
1091  }
1092  }
1093  }
1094 
1095  for (auto ele : this->elements_range()) {
1096  raw_ngh_output_file << ele.idx() << " ";
1097  raw_ngh_output_file << ele->n_sides() << " ";
1098 
1099  auto search_neigh = neigh_vb_map.end();
1100  for (unsigned int i = 0; i < ele->n_sides(); i++) {
1101  unsigned int n_side_neighs = ele.side(i)->edge()->n_sides-1; //n_sides - the current one
1102  // check vb neighbors (lower dimension)
1103  if(n_side_neighs == 0){
1104  //update search
1105  if(search_neigh == neigh_vb_map.end())
1106  search_neigh = neigh_vb_map.find(ele.idx());
1107 
1108  if(search_neigh != neigh_vb_map.end())
1109  if(search_neigh->second[i] != undefined_ele_id)
1110  n_side_neighs = 1;
1111  }
1112  raw_ngh_output_file << n_side_neighs << " ";
1113  }
1114 
1115  for (unsigned int i = 0; i < ele->n_sides(); i++) {
1116  const Edge* edge = ele.side(i)->edge();
1117  if(ele.side(i)->edge()->n_sides > 1){
1118  for (int j = 0; j < edge->n_sides; j++) {
1119  if(edge->side(j) != ele.side(i))
1120  raw_ngh_output_file << edge->side(j)->element().idx() << " ";
1121  }
1122  }
1123  //check vb neighbour
1124  else if(search_neigh != neigh_vb_map.end()
1125  && search_neigh->second[i] != undefined_ele_id){
1126  raw_ngh_output_file << search_neigh->second[i] << " ";
1127  }
1128  }
1129 
1130  // list higher dim neighbours
1131  raw_ngh_output_file << ele->n_neighs_vb() << " ";
1132  for (unsigned int i = 0; i < ele->n_neighs_vb(); i++)
1133  raw_ngh_output_file << ele->neigh_vb[i]->side()->element().idx() << " ";
1134 
1135  raw_ngh_output_file << endl;
1136  cit ++;
1137  }
1138  raw_ngh_output_file << "$EndFlowField\n" << endl;
1139 }
1140 
1141 
1143  unsigned int i, pos;
1145  for (i=0, pos=bulk_size_; i<bc_element_tmp_.size(); ++i, ++pos) {
1146  Element *ele = &element_vec_[pos];
1147  element_ids_.set_item(bc_element_tmp_[i].elm_id, pos);
1148  this->init_element(ele, bc_element_tmp_[i].elm_id, bc_element_tmp_[i].dim, bc_element_tmp_[i].region_idx,
1149  bc_element_tmp_[i].partition_id, bc_element_tmp_[i].node_ids);
1150 
1151  }
1152 
1153  element_vec_.resize(pos); // remove empty element (count is equal with zero-dimensional bulk elements)
1154  bc_element_tmp_.clear();
1155  bc_element_tmp_.reserve(0);
1156 }
1157 
1158 
1159 void Mesh::permute_tetrahedron(unsigned int elm_idx, std::vector<unsigned int> permutation_vec)
1160 {
1161  ASSERT_LT_DBG(elm_idx, element_vec_.size());
1162  ASSERT_EQ_DBG(permutation_vec.size(), 4);
1163 
1164  std::array<unsigned int, 4> tmp_nodes;
1165  Element &elem = element_vec_[elm_idx];
1166  ASSERT_EQ_DBG(elem.dim(), 3);
1167 
1168  for(unsigned int i=0; i<elem.n_nodes(); i++)
1169  {
1170  tmp_nodes[i] = elem.nodes_[permutation_vec[i]];
1171  }
1172  elem.nodes_ = tmp_nodes;
1173 }
1174 
1175 
1176 void Mesh::permute_triangle(unsigned int elm_idx, std::vector<unsigned int> permutation_vec)
1177 {
1178  ASSERT_LT_DBG(elm_idx, element_vec_.size());
1179  ASSERT_EQ_DBG(permutation_vec.size(), 3);
1180 
1181  std::array<unsigned int, 4> tmp_nodes;
1182  Element &elem = element_vec_[elm_idx];
1183  ASSERT_EQ_DBG(elem.dim(), 2);
1184 
1185  for(unsigned int i=0; i<elem.n_nodes(); i++)
1186  {
1187  tmp_nodes[i] = elem.nodes_[permutation_vec[i]];
1188  }
1189  elem.nodes_ = tmp_nodes;
1190 }
1191 
1192 
1194  if (bc_mesh_ == nullptr) bc_mesh_ = new BCMesh(this);
1195  return bc_mesh_;
1196 }
1197 
1198 
1200  ASSERT_PTR(el_4_loc).error("Array 'el_4_loc' is not initialized. Did you call Partitioning::id_maps?\n");
1201 
1202  unsigned int i_proc, i_node, i_ghost_node, elm_node;
1203  unsigned int my_proc = el_ds->myp();
1204 
1205  // distribute nodes between processes, every node is assigned to minimal process of elements that own node
1206  // fill min_node_proc vector with same values on all processes
1207  std::vector<unsigned int> min_node_proc( this->n_nodes(), Mesh::undef_idx );
1208  std::vector<bool> ghost_node_flag( this->n_nodes(), false );
1209  unsigned int n_own_nodes=0, n_ghost_nodes=0; // number of own and ghost nodes
1210  for ( elm : this->elements_range() ) {
1211  i_proc = elm.proc();
1212  for (elm_node=0; elm_node<elm->n_nodes(); elm_node++) {
1213  i_node = elm->node_idx(elm_node);
1214  if ( (min_node_proc[i_node]==Mesh::undef_idx) || (min_node_proc[i_node]>i_proc) ) {
1215  if (i_proc==my_proc) n_own_nodes++;
1216  else if (min_node_proc[i_node]==my_proc) { n_own_nodes--; n_ghost_nodes++; ghost_node_flag[i_node] = true; }
1217  min_node_proc[i_node] = i_proc;
1218  } else if ( !ghost_node_flag[i_node] && (i_proc==my_proc) && (min_node_proc[i_node]!=my_proc) ) {
1219  n_ghost_nodes++;
1220  ghost_node_flag[i_node] = true;
1221  }
1222  }
1223  }
1224 
1225  // create and fill node_4_loc_ (mapping local to global indexes)
1226  node_4_loc_ = new LongIdx [ n_own_nodes+n_ghost_nodes ];
1227  i_node=0;
1228  i_ghost_node = n_own_nodes;
1229  for (unsigned int i=0; i<this->n_nodes(); ++i) {
1230  if (min_node_proc[i]==my_proc) node_4_loc_[i_node++] = i;
1231  if (ghost_node_flag[i]) node_4_loc_[i_ghost_node++] = i;
1232  }
1233 
1234  // Construct node distribution object, set number of local nodes (own+ghost)
1235  node_ds_ = new Distribution(n_own_nodes, PETSC_COMM_WORLD);
1236  node_ds_->get_lsizes_array(); // need to initialize lsizes data member
1237  n_local_nodes_ = n_own_nodes+n_ghost_nodes;
1238 }
1239 
1240 //-----------------------------------------------------------------------------
1241 // 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:1060
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:998
#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:1159
bool is_boundary() const
Returns true if it is a Boundary region and false if it is a Bulk region.
Definition: region.hh:74
friend class Element
Definition: mesh.h:534
void count_side_types()
Definition: mesh.cc:324
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:901
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:1006
virtual const LongIdx * get_local_part()
Definition: mesh.cc:254
MixedMeshIntersections & mixed_intersections()
Definition: mesh.cc:709
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:949
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:916
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:732
#define INPUT_CATCH(ExceptionType, AddressEITag, input_accessor)
Definition: accessors.hh:64
void create_node_element_lists()
Definition: mesh.cc:338
static unsigned int permutation_index(unsigned int p[n_nodes_per_side])
Definition: ref_element.cc:533
BCMesh * get_bc_mesh()
Create boundary mesh if doesn&#39;t exist and return it.
Definition: mesh.cc:1193
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:940
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
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:944
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:250
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:303
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:726
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:235
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:779
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:891
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:132
bool same_sides(const SideIter &si, vector< unsigned int > &side_nodes)
Definition: mesh.cc:403
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:1052
void setup_topology()
Definition: mesh.cc:290
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:240
void set_item(T val, unsigned int pos)
LongIdx * row_4_el
Index set assigning to global element index the local index used in parallel vectors.
Definition: mesh.h: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:501
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:1040
unsigned int n_sides() const
Definition: mesh.cc:226
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:1024
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:738
unsigned int bc_ele_idx_
Definition: boundaries.h:79
void count_element_types()
Definition: mesh.cc:263
bool compare_points(const arma::vec3 &p1, const arma::vec3 &p2)
Definition: mesh.cc:771
const BIHTree & get_bih_tree()
Getter for BIH. Creates and compute BIH at first call.
Definition: mesh.cc:931
void create_boundary_elements()
Create boundary elements from data of temporary structure, this method MUST be call after read mesh f...
Definition: mesh.cc:1142
#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:242
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:1199
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:383
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:351
void modify_element_ids(const RegionDB::MapElementIDToRegionID &map)
Definition: mesh.cc:280
#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:603
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:1046
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:957
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:676
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:979
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:1176
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:420
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:1017
unsigned int id() const
Returns id of the region (using RegionDB)
Definition: region.cc:39
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