Flow123d  last_with_con_2.0.0-663-gd0e2296
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 <boost/tokenizer.hpp>
32 #include "boost/lexical_cast.hpp"
33 
34 #include "mesh/mesh.h"
35 #include "mesh/ref_element.hh"
36 
37 // think about following dependencies
38 #include "mesh/boundaries.h"
39 #include "mesh/accessors.hh"
40 #include "mesh/partitioning.hh"
41 
42 #include "mesh/bih_tree.hh"
43 
44 
45 //TODO: sources, concentrations, initial condition and similarly boundary conditions should be
46 // instances of a Element valued field
47 // concentrations is in fact reimplemented in transport REMOVE it HERE
48 
49 // After removing non-geometrical things from mesh, this should be part of mash initializing.
50 #include "mesh/msh_gmshreader.h"
51 #include "mesh/region.hh"
52 
53 #define NDEF -1
54 
55 namespace IT = Input::Type;
56 
57 
59  return IT::Record("Mesh","Record with mesh related data." )
60  .allow_auto_conversion("mesh_file")
62  "Input file with mesh description.")
64  "List of additional region and region set definitions not contained in the mesh. "
65  "There are three region sets implicitly defined:\n\n"
66  "- ALL (all regions of the mesh)\n"
67  "- .BOUNDARY (all boundary regions)\n"
68  "- BULK (all bulk regions)")
69  .declare_key("partitioning", Partitioning::get_input_type(), IT::Default("\"any_neighboring\""), "Parameters of mesh partitioning algorithms.\n" )
70  .declare_key("print_regions", IT::Bool(), IT::Default("false"), "If true, print table of all used regions.")
71  .close();
72 }
73 
74 
75 
76 const unsigned int Mesh::undef_idx;
77 
79 : in_record_(in_record),
80  comm_(com),
81  row_4_el(nullptr),
82  el_4_loc(nullptr),
83  el_ds(nullptr)
84 {
85  // set in_record_, if input accessor is empty
86  if (in_record_.is_empty()) {
87  istringstream is("{mesh_file=\"\"}");
89  IT::Record &in_rec = const_cast<IT::Record &>(Mesh::get_input_type());
90  in_rec.finish();
91  reader.read_stream(is, in_rec, Input::FileFormat::format_JSON);
93  }
94 
96 }
97 
98 
99 
101 {
102 
103  n_insides = NDEF;
104  n_exsides = NDEF;
105  n_sides_ = NDEF;
106 
107  // number of element of particular dimension
108  n_lines = 0;
109  n_triangles = 0;
110  n_tetrahedras = 0;
111 
112  for (int d=0; d<3; d++) max_edge_sides_[d] = 0;
113 
114  // Initialize numbering of nodes on sides.
115  // This is temporary solution, until class Element is templated
116  // by dimension. Then we can replace Mesh::side_nodes by
117  // RefElement<dim>::side_nodes.
118 
119  // indices of side nodes in element node array
120  // Currently this is made ad libitum
121  // with some ordering here we can get sides with correct orientation.
122  // This speedup normal calculation.
123 
124  side_nodes.resize(3); // three side dimensions
125  for(int i=0; i < 3; i++) {
126  side_nodes[i].resize(i+2); // number of sides
127  for(int j=0; j < i+2; j++)
128  side_nodes[i][j].resize(i+1);
129  }
130 
131  for (unsigned int sid=0; sid<RefElement<1>::n_sides; sid++)
132  for (unsigned int nid=0; nid<RefElement<1>::n_nodes_per_side; nid++)
133  side_nodes[0][sid][nid] = RefElement<1>::side_nodes[sid][nid];
134 
135  for (unsigned int sid=0; sid<RefElement<2>::n_sides; sid++)
136  for (unsigned int nid=0; nid<RefElement<2>::n_nodes_per_side; nid++)
137  side_nodes[1][sid][nid] = RefElement<2>::side_nodes[sid][nid];
138 
139  for (unsigned int sid=0; sid<RefElement<3>::n_sides; sid++)
140  for (unsigned int nid=0; nid<RefElement<3>::n_nodes_per_side; nid++)
141  side_nodes[2][sid][nid] = RefElement<3>::side_nodes[sid][nid];
142 }
143 
144 
146  for(Edge &edg : this->edges)
147  if (edg.side_) delete[] edg.side_;
148 
149  FOR_ELEMENTS( this, ele ) {
150  if (ele->node) delete[] ele->node;
151  if (ele->edge_idx_) delete[] ele->edge_idx_;
152  if (ele->permutation_idx_) delete[] ele->permutation_idx_;
153  if (ele->boundary_idx_) delete[] ele->boundary_idx_;
154  }
155 
156  for(unsigned int idx=0; idx < this->bc_elements.size(); idx++) {
157  Element *ele=&(bc_elements[idx]);
158  if (ele->node) delete[] ele->node;
159  if (ele->edge_idx_) delete[] ele->edge_idx_;
160  if (ele->permutation_idx_) delete[] ele->permutation_idx_;
161  if (ele->boundary_idx_) delete[] ele->boundary_idx_;
162  }
163 
164  if (row_4_el != nullptr) delete[] row_4_el;
165  if (el_4_loc != nullptr) delete[] el_4_loc;
166  if (el_ds != nullptr) delete el_ds;
167 }
168 
169 
170 unsigned int Mesh::n_sides()
171 {
172  if (n_sides_ == NDEF) {
173  n_sides_=0;
174  FOR_ELEMENTS(this, ele) n_sides_ += ele->n_sides();
175  }
176  return n_sides_;
177 }
178 
179 unsigned int Mesh::n_corners() {
180  unsigned int li, count = 0;
181  FOR_ELEMENTS(this, ele) {
182  FOR_ELEMENT_NODES(ele, li) {
183  count++;
184  }
185  }
186  return count;
187 }
188 
190  return part_.get();
191 }
192 
193 
194 //=============================================================================
195 // COUNT ELEMENT TYPES
196 //=============================================================================
197 
199  FOR_ELEMENTS(this, elm)
200  switch (elm->dim()) {
201  case 1:
202  n_lines++;
203  break;
204  case 2:
205  n_triangles++;
206  break;
207  case 3:
208  n_tetrahedras++;
209  break;
210  }
211 }
212 
213 
214 void Mesh::read_gmsh_from_stream(istream &in) {
215 
216  START_TIMER("Reading mesh - from_stream");
217 
218  GmshMeshReader reader(in);
219  reader.read_physical_names(this);
220  reader.read_mesh(this);
221  setup_topology();
222  //close region_db_.
223  region_db_.close();
224 }
225 
226 
227 
229  START_TIMER("Reading mesh - init_from_input");
230 
231  try {
232  Input::Array region_list;
233  // read raw mesh, add regions from GMSH file
234  GmshMeshReader reader( in_record_.val<FilePath>("mesh_file") );
235  reader.read_physical_names(this);
236  // create regions from input
237  if (in_record_.opt_val("regions", region_list)) {
238  this->read_regions_from_input(region_list);
239  }
240  reader.read_mesh(this);
241  } INPUT_CATCH(FilePath::ExcFileOpen, FilePath::EI_Address_String, in_record_)
242  catch (ExceptionBase const &e) {
243  throw;
244  }
245  // possibly add implicit_boundary region.
246  setup_topology();
247  // finish mesh initialization
248  this->check_and_finish();
249 }
250 
251 
252 
253 
255  for (auto elem_to_region : map) {
256  ElementIter ele = this->element.find_id(elem_to_region.first);
257  ele->region_idx_ = region_db_.get_region( elem_to_region.second, ele->dim() );
259  }
260 }
261 
262 
263 
265  START_TIMER("MESH - setup topology");
266 
268 
269  // check mesh quality
270  FOR_ELEMENTS(this, ele)
271  if (ele->quality_measure_smooth() < 0.001) WarningOut().fmt("Bad quality (<0.001) of the element {}.\n", ele.id());
272 
277 
278  part_ = std::make_shared<Partitioning>(this, in_record_.val<Input::Record>("partitioning") );
279 
280  // create parallel distribution and numbering of elements
281  int *id_4_old = new int[element.size()];
282  int i = 0;
283  FOR_ELEMENTS(this, ele)
284  id_4_old[i++] = ele.index();
285  part_->id_maps(element.size(), id_4_old, el_ds, el_4_loc, row_4_el);
286 
287  delete[] id_4_old;
288 }
289 
290 
291 //
293 {
294 
295  n_insides = 0;
296  n_exsides = 0;
297  FOR_SIDES(this, sde ) {
298  if (sde->is_external()) n_exsides++;
299  else n_insides++;
300  }
301 }
302 
303 
304 
306  // for each node we make a list of elements that use this node
307  node_elements.resize(node_vector.size());
308 
309  FOR_ELEMENTS( this, e )
310  for (unsigned int n=0; n<e->n_nodes(); n++)
311  node_elements[node_vector.index(e->node[n])].push_back(e->index());
312 
313  for (vector<vector<unsigned int> >::iterator n=node_elements.begin(); n!=node_elements.end(); n++)
314  stable_sort(n->begin(), n->end());
315 }
316 
317 
318 void Mesh::intersect_element_lists(vector<unsigned int> const &nodes_list, vector<unsigned int> &intersection_element_list)
319 {
320  if (nodes_list.size() == 0) {
321  intersection_element_list.clear();
322  } else if (nodes_list.size() == 1) {
323  intersection_element_list = node_elements[ nodes_list[0] ];
324  } else {
325  vector<unsigned int>::const_iterator it1=nodes_list.begin();
327  intersection_element_list.resize( node_elements[*it1].size() ); // make enough space
328 
329  it1=set_intersection(
330  node_elements[*it1].begin(), node_elements[*it1].end(),
331  node_elements[*it2].begin(), node_elements[*it2].end(),
332  intersection_element_list.begin());
333  intersection_element_list.resize(it1-intersection_element_list.begin()); // resize to true size
334 
335  for(;it2<nodes_list.end();++it2) {
336  it1=set_intersection(
337  intersection_element_list.begin(), intersection_element_list.end(),
338  node_elements[*it2].begin(), node_elements[*it2].end(),
339  intersection_element_list.begin());
340  intersection_element_list.resize(it1-intersection_element_list.begin()); // resize to true size
341  }
342  }
343 }
344 
345 
347  unsigned int dim, unsigned int &element_idx) {
348  bool is_neighbour = false;
349 
350  vector<unsigned int>::iterator e_dest=element_list.begin();
351  for( vector<unsigned int>::iterator ele = element_list.begin(); ele!=element_list.end(); ++ele)
352  if (elements[*ele].dim() == dim) { // keep only indexes of elements of same dimension
353  *e_dest=*ele;
354  ++e_dest;
355  } else if (elements[*ele].dim() == dim-1) { // get only first element of lower dimension
356  if (is_neighbour) xprintf(UsrErr, "Too matching elements id: %d and id: %d in the same mesh.\n",
357  elements(*ele).id(), elements(element_idx).id() );
358 
359  is_neighbour = true;
360  element_idx = *ele;
361  }
362  element_list.resize( e_dest - element_list.begin());
363  return is_neighbour;
364 }
365 
367  // 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)
368  unsigned int ni=0;
369  while ( ni < si->n_nodes()
370  && find(side_nodes.begin(), side_nodes.end(), node_vector.index( si->node(ni) ) ) != side_nodes.end() ) ni++;
371  return ( ni == si->n_nodes() );
372 }
373 
374 /**
375  * TODO:
376  * - use std::is_any for setting is_neigbour
377  * - possibly make appropriate constructors for Edge and Neighbour
378  * - check side!=-1 when searching neigbouring element
379  * - process bc_elements first, there should be no Neigh, but check it
380  * set Edge and boundary there
381  */
382 
384 {
385  Neighbour neighbour;
386  Edge *edg;
387  unsigned int ngh_element_idx, last_edge_idx;
388 
390 
391  // pointers to created edges
392  //vector<Edge *> tmp_edges;
393  edges.resize(0); // be sure that edges are empty
394 
396  vector<unsigned int> intersection_list; // list of elements in intersection of node element lists
397 
398  for( ElementFullIter bc_ele = bc_elements.begin(); bc_ele != bc_elements.end(); ++bc_ele) {
399  // Find all elements that share this side.
400  side_nodes.resize(bc_ele->n_nodes());
401  for (unsigned n=0; n<bc_ele->n_nodes(); n++) side_nodes[n] = node_vector.index(bc_ele->node[n]);
402  intersect_element_lists(side_nodes, intersection_list);
403  bool is_neighbour = find_lower_dim_element(element, intersection_list, bc_ele->dim() +1, ngh_element_idx);
404  if (is_neighbour) {
405  xprintf(UsrErr, "Boundary element (id: %d) match a regular element (id: %d) of lower dimension.\n",
406  bc_ele.id(), element(ngh_element_idx).id());
407  } else {
408  if (intersection_list.size() == 0) {
409  // no matching dim+1 element found
410  WarningOut().fmt("Lonely boundary element, id: {}, region: {}, dimension {}.\n",
411  bc_ele.id(), bc_ele->region().id(), bc_ele->dim());
412  continue; // skip the boundary element
413  }
414  last_edge_idx=edges.size();
415  edges.resize(last_edge_idx+1);
416  edg = &( edges.back() );
417  edg->n_sides = 0;
418  edg->side_ = new struct SideIter[ intersection_list.size() ];
419 
420  // common boundary object
421  unsigned int bdr_idx=boundary_.size();
422  boundary_.resize(bdr_idx+1);
423  Boundary &bdr=boundary_.back();
424  bdr.bc_ele_idx_ = bc_ele.index();
425  bdr.edge_idx_ = last_edge_idx;
426  bdr.mesh_=this;
427 
428  // for 1d boundaries there can be more then one 1d elements connected to the boundary element
429  // we do not detect this case later in the main search over bulk elements
430  for( vector<unsigned int>::iterator isect = intersection_list.begin(); isect!=intersection_list.end(); ++isect) {
431  Element *elem = &(element[*isect]);
432  for (unsigned int ecs=0; ecs<elem->n_sides(); ecs++) {
433  SideIter si = elem->side(ecs);
434  if ( same_sides( si, side_nodes) ) {
435  if (elem->edge_idx_[ecs] != Mesh::undef_idx) {
436  OLD_ASSERT(elem->boundary_idx_!=nullptr, "Null boundary idx array.\n");
437  int last_bc_ele_idx=this->boundary_[elem->boundary_idx_[ecs]].bc_ele_idx_;
438  int new_bc_ele_idx=bc_ele.index();
439  THROW( ExcDuplicateBoundary()
440  << EI_ElemLast(this->bc_elements.get_id(last_bc_ele_idx))
441  << EI_RegLast(this->bc_elements[last_bc_ele_idx].region().label())
442  << EI_ElemNew(this->bc_elements.get_id(new_bc_ele_idx))
443  << EI_RegNew(this->bc_elements[new_bc_ele_idx].region().label())
444  );
445  }
446  elem->edge_idx_[ecs] = last_edge_idx;
447  edg->side_[ edg->n_sides++ ] = si;
448 
449  if (elem->boundary_idx_ == NULL) {
450  elem->boundary_idx_ = new unsigned int [ elem->n_sides() ];
451  std::fill( elem->boundary_idx_, elem->boundary_idx_ + elem->n_sides(), Mesh::undef_idx);
452  }
453  elem->boundary_idx_[ecs] = bdr_idx;
454  break; // next element in intersection list
455  }
456  }
457  }
458 
459  }
460 
461  }
462  // Now we go through all element sides and create edges and neighbours
463  FOR_ELEMENTS( this, e )
464  {
465  for (unsigned int s=0; s<e->n_sides(); s++)
466  {
467  // skip sides that were already found
468  if (e->edge_idx_[s] != Mesh::undef_idx) continue;
469 
470 
471  // Find all elements that share this side.
472  side_nodes.resize(e->side(s)->n_nodes());
473  for (unsigned n=0; n<e->side(s)->n_nodes(); n++) side_nodes[n] = node_vector.index(e->side(s)->node(n));
474  intersect_element_lists(side_nodes, intersection_list);
475 
476  bool is_neighbour = find_lower_dim_element(element, intersection_list, e->dim(), ngh_element_idx);
477 
478  if (is_neighbour) { // edge connects elements of different dimensions
479  neighbour.element_ = &(element[ngh_element_idx]);
480  } else { // edge connects only elements of the same dimension
481  // Allocate the array of sides.
482  last_edge_idx=edges.size();
483  edges.resize(last_edge_idx+1);
484  edg = &( edges.back() );
485  edg->n_sides = 0;
486  edg->side_ = new struct SideIter[ intersection_list.size() ];
487  if (intersection_list.size() > max_edge_sides_[e->dim()-1])
488  max_edge_sides_[e->dim()-1] = intersection_list.size();
489 
490  if (intersection_list.size() == 1) { // outer edge, create boundary object as well
491  edg->n_sides=1;
492  edg->side_[0] = e->side(s);
493  e->edge_idx_[s] = last_edge_idx;
494 
495  if (e->boundary_idx_ == NULL) {
496  e->boundary_idx_ = new unsigned int [ e->n_sides() ];
497  std::fill( e->boundary_idx_, e->boundary_idx_ + e->n_sides(), Mesh::undef_idx);
498  }
499 
500  unsigned int bdr_idx=boundary_.size();
501  boundary_.resize(bdr_idx+1);
502  Boundary &bdr=boundary_.back();
503  e->boundary_idx_[s] = bdr_idx;
504 
505  // fill boundary element
506  ElementFullIter bc_ele = bc_elements.add_item( -bdr_idx ); // use negative bcd index as ID,
507  bc_ele->init(e->dim()-1, this, region_db_.implicit_boundary_region() );
508  region_db_.mark_used_region( bc_ele->region_idx_.idx() );
509  for(unsigned int ni = 0; ni< side_nodes.size(); ni++) bc_ele->node[ni] = &( node_vector[side_nodes[ni]] );
510 
511  // fill Boundary object
512  bdr.edge_idx_ = last_edge_idx;
513  bdr.bc_ele_idx_ = bc_ele.index();
514  bdr.mesh_=this;
515 
516  continue; // next side of element e
517  }
518  }
519 
520  // go through the elements connected to the edge or neighbour
521  for( vector<unsigned int>::iterator isect = intersection_list.begin(); isect!=intersection_list.end(); ++isect) {
522  Element *elem = &(element[*isect]);
523  for (unsigned int ecs=0; ecs<elem->n_sides(); ecs++) {
524  if (elem->edge_idx_[ecs] != Mesh::undef_idx) continue;
525  SideIter si = elem->side(ecs);
526  if ( same_sides( si, side_nodes) ) {
527  if (is_neighbour) {
528  // create a new edge and neighbour for this side, and element to the edge
529  last_edge_idx=edges.size();
530  edges.resize(last_edge_idx+1);
531  edg = &( edges.back() );
532  edg->n_sides = 1;
533  edg->side_ = new struct SideIter[1];
534  edg->side_[0] = si;
535  elem->edge_idx_[ecs] = last_edge_idx;
536 
537  neighbour.edge_idx_ = last_edge_idx;
538 
539  vb_neighbours_.push_back(neighbour); // copy neighbour with this edge setting
540  } else {
541  // connect the side to the edge, and side to the edge
542  edg->side_[ edg->n_sides++ ] = si;
543  elem->edge_idx_[ecs] = last_edge_idx;
544  }
545  break; // next element from intersection list
546  }
547  } // search for side of other connected element
548  } // connected elements
549  OLD_ASSERT( is_neighbour || ( (unsigned int) edg->n_sides ) == intersection_list.size(), "Some connected sides were not found.\n");
550  } // for element sides
551  } // for elements
552 
553  MessageOut().fmt( "Created {} edges and {} neighbours.\n", edges.size(), vb_neighbours_.size() );
554 }
555 
556 
557 
559 {
560  for (EdgeVector::iterator edg=edges.begin(); edg!=edges.end(); edg++)
561  {
562  // side 0 is reference, so its permutation is 0
563  edg->side(0)->element()->permutation_idx_[edg->side(0)->el_idx()] = 0;
564 
565  if (edg->n_sides > 1)
566  {
567  map<const Node*,unsigned int> node_numbers;
568  unsigned int permutation[edg->side(0)->n_nodes()];
569 
570  for (unsigned int i=0; i<edg->side(0)->n_nodes(); i++)
571  node_numbers[edg->side(0)->node(i)] = i;
572 
573  for (int sid=1; sid<edg->n_sides; sid++)
574  {
575  for (unsigned int i=0; i<edg->side(0)->n_nodes(); i++)
576  permutation[node_numbers[edg->side(sid)->node(i)]] = i;
577 
578  switch (edg->side(0)->dim())
579  {
580  case 0:
581  edg->side(sid)->element()->permutation_idx_[edg->side(sid)->el_idx()] = RefElement<1>::permutation_index(permutation);
582  break;
583  case 1:
584  edg->side(sid)->element()->permutation_idx_[edg->side(sid)->el_idx()] = RefElement<2>::permutation_index(permutation);
585  break;
586  case 2:
587  edg->side(sid)->element()->permutation_idx_[edg->side(sid)->el_idx()] = RefElement<3>::permutation_index(permutation);
588  break;
589  }
590  }
591  }
592  }
593 
594  for (vector<Neighbour>::iterator nb=vb_neighbours_.begin(); nb!=vb_neighbours_.end(); nb++)
595  {
596  map<const Node*,unsigned int> node_numbers;
597  unsigned int permutation[nb->element()->n_nodes()];
598 
599  // element of lower dimension is reference, so
600  // we calculate permutation for the adjacent side
601  for (unsigned int i=0; i<nb->element()->n_nodes(); i++)
602  node_numbers[nb->element()->node[i]] = i;
603 
604  for (unsigned int i=0; i<nb->side()->n_nodes(); i++)
605  permutation[node_numbers[nb->side()->node(i)]] = i;
606 
607  switch (nb->side()->dim())
608  {
609  case 0:
610  nb->side()->element()->permutation_idx_[nb->side()->el_idx()] = RefElement<1>::permutation_index(permutation);
611  break;
612  case 1:
613  nb->side()->element()->permutation_idx_[nb->side()->el_idx()] = RefElement<2>::permutation_index(permutation);
614  break;
615  case 2:
616  nb->side()->element()->permutation_idx_[nb->side()->el_idx()] = RefElement<3>::permutation_index(permutation);
617  break;
618  }
619  }
620 }
621 
622 
623 
624 
625 
626 //=============================================================================
627 //
628 //=============================================================================
630 {
631 
632  //MessageOut() << "Element to neighbours of vb2 type... "/*orig verb 5*/;
633 
634  FOR_ELEMENTS(this,ele) ele->n_neighs_vb =0;
635 
636  // count vb neighs per element
637  FOR_NEIGHBOURS(this, ngh ) ngh->element_->n_neighs_vb++;
638 
639  // Allocation of the array per element
640  FOR_ELEMENTS(this, ele )
641  if( ele->n_neighs_vb > 0 ) {
642  ele->neigh_vb = new struct Neighbour* [ele->n_neighs_vb];
643  ele->n_neighs_vb=0;
644  }
645 
646  // fill
647  ElementIter ele;
648  FOR_NEIGHBOURS(this, ngh ) {
649  ele = ngh->element();
650  ele->neigh_vb[ ele->n_neighs_vb++ ] = &( *ngh );
651  }
652 
653  //MessageOut() << "... O.K.\n"/*orig verb 6*/;
654 }
655 
656 
657 
658 
662 
663 
665  /* Algorithm:
666  *
667  * 1) create BIH tree
668  * 2) for every 1D, find list of candidates
669  * 3) compute intersections for 1d, store it to master_elements
670  *
671  */
672  const BIHTree &bih_tree =get_bih_tree();
673  master_elements.resize(n_elements());
674 
675  for(unsigned int i_ele=0; i_ele<n_elements(); i_ele++) {
676  Element &ele = this->element[i_ele];
677 
678  if (ele.dim() == 1) {
679  vector<unsigned int> candidate_list;
680  bih_tree.find_bounding_box(ele.bounding_box(), candidate_list);
681 
682  //for(unsigned int i_elm=0; i_elm<n_elements(); i_elm++) {
683  for(unsigned int i_elm : candidate_list) {
684  ElementFullIter elm = this->element( i_elm );
685  if (elm->dim() == 2) {
686  IntersectionLocal *intersection;
687  GetIntersection( TAbscissa(ele), TTriangle(*elm), intersection);
688  if (intersection && intersection->get_type() == IntersectionLocal::line) {
689 
690  master_elements[i_ele].push_back( intersections.size() );
691  intersections.push_back( Intersection(this->element(i_ele), elm, intersection) );
692  }
693  }
694 
695  }
696  }
697  }
698 
699 }
700 
701 
702 
703 ElementAccessor<3> Mesh::element_accessor(unsigned int idx, bool boundary) {
704  return ElementAccessor<3>(this, idx, boundary);
705 }
706 
707 
708 
709 vector<int> const & Mesh::elements_id_maps( bool boundary_domain) const
710 {
711  if (bulk_elements_id_.size() ==0) {
713  int last_id;
714 
715  bulk_elements_id_.resize(n_elements());
716  map_it = bulk_elements_id_.begin();
717  last_id = -1;
718  for(unsigned int idx=0; idx < element.size(); idx++, ++map_it) {
719  int id = element.get_id(idx);
720  if (last_id >= id) xprintf(UsrErr, "Element IDs in non-increasing order, ID: %d\n", id);
721  last_id=*map_it = id;
722  }
723 
725  map_it = boundary_elements_id_.begin();
726  last_id = -1;
727  for(unsigned int idx=0; idx < bc_elements.size(); idx++, ++map_it) {
728  int id = bc_elements.get_id(idx);
729  // We set ID for boundary elements created by the mesh itself to "-1"
730  // this force gmsh reader to skip all remaining entries in boundary_elements_id_
731  // and thus report error for any remaining data lines
732  if (id < 0) last_id=*map_it=-1;
733  else {
734  if (last_id >= id) xprintf(UsrErr, "Element IDs in non-increasing order, ID: %d\n", id);
735  last_id=*map_it = id;
736  }
737  }
738  }
739 
740  if (boundary_domain) return boundary_elements_id_;
741  return bulk_elements_id_;
742 }
743 
745 {
747  it != region_list.end();
748  ++it) {
749  // constructor has side effect in the mesh - create new region or set and store them to Mesh::region_db_
750  (*it).factory< RegionSetBase, const Input::Record &, Mesh * >(*it, this);
751  }
752 }
753 
755 {
757  region_db_.el_to_reg_map_.clear();
758  region_db_.close();
760 
761  if ( in_record_.val<bool>("print_regions") ) {
763  }
764 }
765 
766 
768  if (! this->bih_tree_)
769  bih_tree_ = std::make_shared<BIHTree>(this);
770  return *bih_tree_;
771 }
772 
773 //-----------------------------------------------------------------------------
774 // vim: set cindent:
int n_triangles
Definition: mesh.h:253
Distribution * el_ds
Parallel distribution of elements.
Definition: mesh.h:392
Class for the mesh partitioning. This should provide:
Definition: partitioning.hh:39
Mesh(Input::Record in_record, MPI_Comm com=MPI_COMM_WORLD)
Definition: mesh.cc:78
Iterator< ValueType > begin() const
void read_stream(istream &in, const Type::TypeBase &root_type, FileFormat format)
This method actually reads the given stream in.
#define FOR_ELEMENT_NODES(i, j)
Definition: elements.h:188
vector< vector< unsigned int > > node_elements
Definition: mesh.h:275
Accessor to input data conforming to declared Array.
Definition: accessors.hh:567
unsigned int * boundary_idx_
Definition: elements.h:83
int get_id(const T *it) const
Definition: sys_vector.hh:458
void count_side_types()
Definition: mesh.cc:292
int MPI_Comm
Definition: mpi.h:141
void make_intersec_elements()
Definition: mesh.cc:664
void check_and_finish()
Definition: mesh.cc:754
MapElementIDToRegionID el_to_reg_map_
Definition: region.hh:570
Reader for (slightly) modified input files.
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:50
#define FOR_ELEMENTS(_mesh_, __i)
Definition: mesh.h:405
Class for declaration of the input of type Bool.
Definition: type_base.hh:465
#define MessageOut()
Macro defining &#39;message&#39; record of log.
Definition: logger.hh:231
unsigned int * permutation_idx_
Definition: elements.h:92
int n_lines
Definition: mesh.h:252
static const Input::Type::Record & get_input_type()
Definition: partitioning.cc:45
static const unsigned int undef_idx
Definition: mesh.h:108
#define INPUT_CATCH(ExceptionType, AddressEITag, input_accessor)
Definition: accessors.hh:57
void create_node_element_lists()
Definition: mesh.cc:305
static unsigned int permutation_index(unsigned int p[n_nodes_per_side])
Definition: ref_element.cc:346
static const Input::Type::Record & get_input_type()
Definition: mesh.cc:58
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:99
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:355
Definition: mesh.h:95
int index() const
Definition: sys_vector.hh:78
Input::Record in_record_
Definition: mesh.h:374
Node ** node
Definition: elements.h:79
FullIter add_item(int id)
Definition: sys_vector.hh:359
int n_sides
Definition: edges.h:36
vector< vector< vector< unsigned int > > > side_nodes
Definition: mesh.h:265
FullIter find_id(const int id)
Definition: sys_vector.hh:434
Definition: edges.h:26
int n_tetrahedras
Definition: mesh.h:254
ElementAccessor< 3 > element_accessor(unsigned int idx, bool boundary=false)
Definition: mesh.cc:703
vector< Boundary > boundary_
Definition: mesh.h:220
Partitioning * get_part()
Definition: mesh.cc:189
unsigned int n_sides()
Definition: mesh.cc:170
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:348
int * el_4_loc
Index set assigning to local element index its global index.
Definition: mesh.h:390
int * row_4_el
Index set assigning to global element index the local index used in parallel vectors.
Definition: mesh.h:388
int n_exsides
Definition: mesh.h:249
IteratorBase end() const
Region get_region(unsigned int id, unsigned int dim)
Definition: region.cc:150
unsigned int dim() const
unsigned int size() const
Returns size of the container. This is independent of the allocated space.
Definition: sys_vector.hh:391
static Default optional()
The factory function to make an empty default value which is optional.
Definition: type_record.hh:113
#define OLD_ASSERT(...)
Definition: global_defs.h:131
bool opt_val(const string &key, Ret &value) const
unsigned int edge_idx_
Definition: neighbours.h:135
unsigned int * edge_idx_
Definition: elements.h:82
void read_regions_from_input(Input::Array region_list)
Definition: mesh.cc:744
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
unsigned int n_elements() const
Definition: mesh.h:133
bool same_sides(const SideIter &si, vector< unsigned int > &side_nodes)
Definition: mesh.cc:366
static FileName input()
The factory function for declaring type FileName for input files.
Definition: type_base.hh:637
std::vector< T >::iterator iterator
Definition: sys_vector.hh:215
vector< int > const & elements_id_maps(bool boundary_domain) const
Definition: mesh.cc:709
void setup_topology()
Definition: mesh.cc:264
Neighbour ** neigh_vb
Definition: elements.h:166
T get_root_interface() const
Returns the root accessor.
friend class RegionSetBase
Definition: mesh.h:382
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
#define NDEF
Definition: mesh.cc:53
const Ret val(const string &key) const
unsigned int n_sides() const
#define xprintf(...)
Definition: system.hh:87
#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:36
unsigned int n_corners()
Definition: mesh.cc:179
unsigned int index(const T *pointer) const
Definition: sys_vector.hh:373
#define FOR_NEIGHBOURS(_mesh_, it)
Definition: mesh.h:456
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:484
ElementVector bc_elements
Definition: mesh.h:224
void close()
Definition: region.cc:249
SideIter side(const unsigned int loc_index)
FullIter begin()
Definition: sys_vector.hh:383
bool is_empty() const
Definition: accessors.hh:366
Region implicit_boundary_region()
Definition: region.cc:75
std::shared_ptr< BIHTree > bih_tree_
Definition: mesh.h:370
unsigned int edge_idx_
Definition: boundaries.h:75
void read_mesh(Mesh *mesh)
Accessor to the polymorphic input data of a type given by an AbstracRecord object.
Definition: accessors.hh:459
vector< int > boundary_elements_id_
Definition: mesh.h:350
unsigned int bc_ele_idx_
Definition: boundaries.h:76
void find_bounding_box(const BoundingBox &boundingBox, std::vector< unsigned int > &result_list) const
Definition: bih_tree.cc:193
void count_element_types()
Definition: mesh.cc:198
const BIHTree & get_bih_tree()
Definition: mesh.cc:767
void check_regions()
Definition: region.cc:461
Dedicated class for storing path to input and output files.
Definition: file_path.hh:48
Definition: system.hh:59
void read_gmsh_from_stream(istream &in)
Definition: mesh.cc:214
void mark_used_region(unsigned int idx)
Definition: region.cc:235
void print_region_table(ostream &stream) const
Definition: region.cc:409
Support classes for parallel programing.
vector< vector< unsigned int > > master_elements
Definition: mesh.h:241
vector< Neighbour > vb_neighbours_
Definition: mesh.h:246
FinishStatus finish(FinishStatus finish_type=FinishStatus::regular_) override
Finish declaration of the Record type.
Definition: type_record.cc:242
int n_sides_
Definition: mesh.h:250
RegionIdx region_idx_
Definition: elements.h:174
void read_physical_names(Mesh *mesh)
RegionDB region_db_
Definition: mesh.h:361
std::vector< Edge > edges
Vector of MH edges, this should not be part of the geometrical mesh.
Definition: mesh.h:227
vector< Intersection > intersections
Definition: mesh.h:235
MPI_Comm comm_
Definition: mesh.h:379
unsigned int n_neighs_vb
Definition: elements.h:164
int n_insides
Definition: mesh.h:248
void intersect_element_lists(vector< unsigned int > const &nodes_list, vector< unsigned int > &intersection_element_list)
Definition: mesh.cc:318
void modify_element_ids(const RegionDB::MapElementIDToRegionID &map)
Definition: mesh.cc:254
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:234
void make_edge_permutations()
Definition: mesh.cc:558
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc...
static Input::Type::Abstract & get_input_type()
Definition: region_set.cc:25
Mesh * mesh_
Definition: boundaries.h:77
unsigned int n_nodes() const
Definition: mesh.h:129
Record type proxy class.
Definition: type_record.hh:171
ElementIter element_
Definition: neighbours.h:136
Base of exceptions used in Flow123d.
Definition: exceptions.hh:67
SideIter * side_
Definition: edges.h:37
void element_to_neigh_vb()
Definition: mesh.cc:629
const Node * node(unsigned int i) const
Definition: side_impl.hh:46
void reinit(Input::Record in_record)
Definition: mesh.cc:100
FullIter end()
Returns FullFullIterer of the fictions past the end element.
Definition: sys_vector.hh:387
#define FOR_SIDES(_mesh_, it)
Definition: mesh.h:449
void init_from_input()
Definition: mesh.cc:228
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
Definition: exceptions.hh:45
void GetIntersection(const TBisector &, const TBisector &, TPosition &, double &, double &)
bool find_lower_dim_element(ElementVector &elements, vector< unsigned int > &element_list, unsigned int dim, unsigned int &element_idx)
Definition: mesh.cc:346
vector< int > bulk_elements_id_
Definition: mesh.h:350
IntersectionType get_type() const
std::shared_ptr< Partitioning > part_
Definition: mesh.h:365
void make_neighbours_and_edges()
Definition: mesh.cc:383
NodeVector node_vector
Vector of nodes of the mesh.
Definition: mesh.h:214
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:216
unsigned int n_nodes() const
Definition: side_impl.hh:33
~Mesh()
Destructor.
Definition: mesh.cc:145
unsigned int idx() const
Returns a global index of the region.
Definition: region.hh:81
BoundingBox bounding_box()
Definition: elements.h:102