Flow123d  jenkins-Flow123d-linux-release-multijob-282
mesh.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @ingroup mesh
27  * @brief Mesh construction
28  *
29  */
30 
31 #include <unistd.h>
32 #include <set>
33 
34 
35 #include "system/system.hh"
36 #include "system/xio.h"
37 #include "input/json_to_storage.hh"
38 #include "input/input_type.hh"
39 #include "system/sys_profiler.hh"
40 
41 #include <boost/tokenizer.hpp>
42 #include "boost/lexical_cast.hpp"
43 #include <boost/make_shared.hpp>
44 
45 #include "mesh/mesh.h"
46 #include "mesh/ref_element.hh"
47 
48 // think about following dependencies
49 #include "mesh/boundaries.h"
50 #include "mesh/accessors.hh"
51 #include "mesh/partitioning.hh"
52 
53 #include "mesh/bih_tree.hh"
54 
55 
56 //TODO: sources, concentrations, initial condition and similarly boundary conditions should be
57 // instances of a Element valued field
58 // concentrations is in fact reimplemented in transport REMOVE it HERE
59 
60 // After removing non-geometrical things from mesh, this should be part of mash initializing.
61 #include "mesh/msh_gmshreader.h"
62 #include "mesh/region.hh"
63 
64 #define NDEF -1
65 
66 namespace IT = Input::Type;
67 
68 
70  = IT::Record("Mesh","Record with mesh related data." )
72  "Input file with mesh description.")
74  "List of additional region definitions not contained in the mesh.")
76  "List of region set definitions. There are three region sets implicitly defined:\n"
77  "ALL (all regions of the mesh), BOUNDARY (all boundary regions), and BULK (all bulk regions)")
78  .declare_key("partitioning", Partitioning::input_type, IT::Default("any_neighboring"), "Parameters of mesh partitioning algorithms.\n" )
79  .close();
80 
81 
82 
83 const unsigned int Mesh::undef_idx;
84 
85 Mesh::Mesh(const std::string &input_str, MPI_Comm comm)
86 :comm_(comm)
87 {
88 
89  Input::JSONToStorage reader( input_str, Mesh::input_type );
91 
93 }
94 
95 
96 
98 : in_record_(in_record),
99  comm_(com)
100 {
102 }
103 
104 
105 
107 {
108 
109  n_insides = NDEF;
110  n_exsides = NDEF;
111  n_sides_ = NDEF;
112 
113  // number of element of particular dimension
114  n_lines = 0;
115  n_triangles = 0;
116  n_tetrahedras = 0;
117 
118  for (int d=0; d<3; d++) max_edge_sides_[d] = 0;
119 
120  // Initialize numbering of nodes on sides.
121  // This is temporary solution, until class Element is templated
122  // by dimension. Then we can replace Mesh::side_nodes by
123  // RefElement<dim>::side_nodes.
124 
125  // indices of side nodes in element node array
126  // Currently this is made ad libitum
127  // with some ordering here we can get sides with correct orientation.
128  // This speedup normal calculation.
129 
130  side_nodes.resize(3); // three side dimensions
131  for(int i=0; i < 3; i++) {
132  side_nodes[i].resize(i+2); // number of sides
133  for(int j=0; j < i+2; j++)
134  side_nodes[i][j].resize(i+1);
135  }
136 
137  for (unsigned int sid=0; sid<RefElement<1>::n_sides; sid++)
138  for (unsigned int nid=0; nid<RefElement<1>::n_nodes_per_side; nid++)
139  side_nodes[0][sid][nid] = RefElement<1>::side_nodes[sid][nid];
140 
141  for (unsigned int sid=0; sid<RefElement<2>::n_sides; sid++)
142  for (unsigned int nid=0; nid<RefElement<2>::n_nodes_per_side; nid++)
143  side_nodes[1][sid][nid] = RefElement<2>::side_nodes[sid][nid];
144 
145  for (unsigned int sid=0; sid<RefElement<3>::n_sides; sid++)
146  for (unsigned int nid=0; nid<RefElement<3>::n_nodes_per_side; nid++)
147  side_nodes[2][sid][nid] = RefElement<3>::side_nodes[sid][nid];
148 }
149 
150 
152  for(Edge &edg : this->edges)
153  if (edg.side_) delete[] edg.side_;
154 
155  FOR_ELEMENTS( this, ele ) {
156  if (ele->node) delete[] ele->node;
157  if (ele->edge_idx_) delete[] ele->edge_idx_;
158  if (ele->permutation_idx_) delete[] ele->permutation_idx_;
159  if (ele->boundary_idx_) delete[] ele->boundary_idx_;
160  }
161 
162  for(unsigned int idx=0; idx < this->bc_elements.size(); idx++) {
163  Element *ele=&(bc_elements[idx]);
164  if (ele->node) delete[] ele->node;
165  if (ele->edge_idx_) delete[] ele->edge_idx_;
166  if (ele->permutation_idx_) delete[] ele->permutation_idx_;
167  if (ele->boundary_idx_) delete[] ele->boundary_idx_;
168  }
169 }
170 
171 
172 unsigned int Mesh::n_sides()
173 {
174  if (n_sides_ == NDEF) {
175  n_sides_=0;
176  FOR_ELEMENTS(this, ele) n_sides_ += ele->n_sides();
177  }
178  return n_sides_;
179 }
180 
181 unsigned int Mesh::n_corners() {
182  unsigned int li, count = 0;
183  FOR_ELEMENTS(this, ele) {
184  FOR_ELEMENT_NODES(ele, li) {
185  count++;
186  }
187  }
188  return count;
189 }
190 
192  return part_.get();
193 }
194 
195 
196 //=============================================================================
197 // COUNT ELEMENT TYPES
198 //=============================================================================
199 
201  FOR_ELEMENTS(this, elm)
202  switch (elm->dim()) {
203  case 1:
204  n_lines++;
205  break;
206  case 2:
207  n_triangles++;
208  break;
209  case 3:
210  n_tetrahedras++;
211  break;
212  }
213 }
214 
215 
216 void Mesh::read_gmsh_from_stream(istream &in) {
217 
218  START_TIMER("Reading mesh - from_stream");
219 
220  GmshMeshReader reader(in);
221  reader.read_mesh(this);
222  setup_topology();
223 }
224 
225 
226 
228  START_TIMER("Reading mesh - init_from_input");
229 
230  Input::Array region_list;
231  RegionDB::MapElementIDToRegionID el_to_reg_map;
232 
233  // create regions from our input
234  if (in_record_.opt_val("regions", region_list)) {
235  region_db_.read_regions_from_input(region_list, el_to_reg_map);
236  }
237  // read raw mesh, add regions from GMSH file
238  GmshMeshReader reader( in_record_.val<FilePath>("mesh_file") );
239  reader.read_mesh(this, &el_to_reg_map);
240  // possibly add implicit_boundary region, close region_db_.
241  setup_topology();
242  // create sets
243  Input::Array set_list;
244  if (in_record_.opt_val("sets", set_list)) {
246  }
247 }
248 
249 
250 
251 
253  START_TIMER("MESH - setup topology");
254 
256 
257  // check mesh quality
258  FOR_ELEMENTS(this, ele)
259  if (ele->quality_measure_smooth() < 0.001) xprintf(Warn, "Bad quality (<0.001) of the element %u.\n", ele.id());
260 
265 
266  region_db_.close();
267  part_ = boost::make_shared<Partitioning>(this, in_record_.val<Input::Record>("partitioning") );
268 }
269 
270 
271 //
273 {
274 
275  n_insides = 0;
276  n_exsides = 0;
277  FOR_SIDES(this, sde ) {
278  if (sde->is_external()) n_exsides++;
279  else n_insides++;
280  }
281 }
282 
283 
284 
286  // for each node we make a list of elements that use this node
287  node_elements.resize(node_vector.size());
288 
289  FOR_ELEMENTS( this, e )
290  for (unsigned int n=0; n<e->n_nodes(); n++)
291  node_elements[node_vector.index(e->node[n])].push_back(e->index());
292 
293  for (vector<vector<unsigned int> >::iterator n=node_elements.begin(); n!=node_elements.end(); n++)
294  stable_sort(n->begin(), n->end());
295 }
296 
297 
298 void Mesh::intersect_element_lists(vector<unsigned int> const &nodes_list, vector<unsigned int> &intersection_element_list)
299 {
300  if (nodes_list.size() == 0) {
301  intersection_element_list.clear();
302  } else if (nodes_list.size() == 1) {
303  intersection_element_list = node_elements[ nodes_list[0] ];
304  } else {
305  vector<unsigned int>::const_iterator it1=nodes_list.begin();
307  intersection_element_list.resize( node_elements[*it1].size() ); // make enough space
308 
309  it1=set_intersection(
310  node_elements[*it1].begin(), node_elements[*it1].end(),
311  node_elements[*it2].begin(), node_elements[*it2].end(),
312  intersection_element_list.begin());
313  intersection_element_list.resize(it1-intersection_element_list.begin()); // resize to true size
314 
315  for(;it2<nodes_list.end();++it2) {
316  it1=set_intersection(
317  intersection_element_list.begin(), intersection_element_list.end(),
318  node_elements[*it2].begin(), node_elements[*it2].end(),
319  intersection_element_list.begin());
320  intersection_element_list.resize(it1-intersection_element_list.begin()); // resize to true size
321  }
322  }
323 }
324 
325 
327  unsigned int dim, unsigned int &element_idx) {
328  bool is_neighbour = false;
329 
330  vector<unsigned int>::iterator e_dest=element_list.begin();
331  for( vector<unsigned int>::iterator ele = element_list.begin(); ele!=element_list.end(); ++ele)
332  if (elements[*ele].dim() == dim) { // keep only indexes of elements of same dimension
333  *e_dest=*ele;
334  ++e_dest;
335  } else if (elements[*ele].dim() == dim-1) { // get only first element of lower dimension
336  if (is_neighbour) xprintf(UsrErr, "Too matching elements id: %d and id: %d in the same mesh.\n",
337  elements(*ele).id(), elements(element_idx).id() );
338 
339  is_neighbour = true;
340  element_idx = *ele;
341  }
342  element_list.resize( e_dest - element_list.begin());
343  return is_neighbour;
344 }
345 
346 bool Mesh::same_sides(const SideIter &si, vector<unsigned int> &side_nodes) {
347  // 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)
348  unsigned int ni=0;
349  while ( ni < si->n_nodes()
350  && find(side_nodes.begin(), side_nodes.end(), node_vector.index( si->node(ni) ) ) != side_nodes.end() ) ni++;
351  return ( ni == si->n_nodes() );
352 }
353 
354 /**
355  * TODO:
356  * - use std::is_any for setting is_neigbour
357  * - possibly make appropriate constructors for Edge and Neighbour
358  * - check side!=-1 when searching neigbouring element
359  * - process bc_elements first, there should be no Neigh, but check it
360  * set Edge and boundary there
361  */
362 
364 {
365  Neighbour neighbour;
366  Edge *edg;
367  unsigned int ngh_element_idx, last_edge_idx;
368 
370 
371  // pointers to created edges
372  //vector<Edge *> tmp_edges;
373  edges.resize(0); // be sure that edges are empty
374 
376  vector<unsigned int> intersection_list; // list of elements in intersection of node element lists
377 
378  for( ElementFullIter bc_ele = bc_elements.begin(); bc_ele != bc_elements.end(); ++bc_ele) {
379  // Find all elements that share this side.
380  side_nodes.resize(bc_ele->n_nodes());
381  for (unsigned n=0; n<bc_ele->n_nodes(); n++) side_nodes[n] = node_vector.index(bc_ele->node[n]);
382  intersect_element_lists(side_nodes, intersection_list);
383  bool is_neighbour = find_lower_dim_element(element, intersection_list, bc_ele->dim() +1, ngh_element_idx);
384  if (is_neighbour) {
385  xprintf(UsrErr, "Boundary element (id: %d) match a regular element (id: %d) of lower dimension.\n",
386  bc_ele.id(), element(ngh_element_idx).id());
387  } else {
388  if (intersection_list.size() == 0) {
389  // no matching dim+1 element found
390  xprintf(Warn, "Lonely boundary element, id: %d, region: %d, dimension %d.\n", bc_ele.id(), bc_ele->region().id(), bc_ele->dim());
391  continue; // skip the boundary element
392  }
393  last_edge_idx=edges.size();
394  edges.resize(last_edge_idx+1);
395  edg = &( edges.back() );
396  edg->n_sides = 0;
397  edg->side_ = new struct SideIter[ intersection_list.size() ];
398 
399  // common boundary object
400  unsigned int bdr_idx=boundary_.size();
401  boundary_.resize(bdr_idx+1);
402  Boundary &bdr=boundary_.back();
403  bdr.bc_ele_idx_ = bc_ele.index();
404  bdr.edge_idx_ = last_edge_idx;
405  bdr.mesh_=this;
406 
407  // for 1d boundaries there can be more then one 1d elements connected to the boundary element
408  // we do not detect this case later in the main search over bulk elements
409  for( vector<unsigned int>::iterator isect = intersection_list.begin(); isect!=intersection_list.end(); ++isect) {
410  Element *elem = &(element[*isect]);
411  for (unsigned int ecs=0; ecs<elem->n_sides(); ecs++) {
412  SideIter si = elem->side(ecs);
413  if ( same_sides( si, side_nodes) ) {
414  edg->side_[ edg->n_sides++ ] = si;
415  elem->edge_idx_[ecs] = last_edge_idx;
416 
417  if (elem->boundary_idx_ == NULL) {
418  elem->boundary_idx_ = new unsigned int [ elem->n_sides() ];
419  std::fill( elem->boundary_idx_, elem->boundary_idx_ + elem->n_sides(), Mesh::undef_idx);
420  }
421  elem->boundary_idx_[ecs] = bdr_idx;
422  break; // next element in intersection list
423  }
424  }
425  }
426 
427  }
428 
429  }
430  // Now we go through all element sides and create edges and neighbours
431  FOR_ELEMENTS( this, e )
432  {
433  for (unsigned int s=0; s<e->n_sides(); s++)
434  {
435  // skip sides that were already found
436  if (e->edge_idx_[s] != Mesh::undef_idx) continue;
437 
438 
439  // Find all elements that share this side.
440  side_nodes.resize(e->side(s)->n_nodes());
441  for (unsigned n=0; n<e->side(s)->n_nodes(); n++) side_nodes[n] = node_vector.index(e->side(s)->node(n));
442  intersect_element_lists(side_nodes, intersection_list);
443 
444  bool is_neighbour = find_lower_dim_element(element, intersection_list, e->dim(), ngh_element_idx);
445 
446  if (is_neighbour) { // edge connects elements of different dimensions
447  neighbour.element_ = &(element[ngh_element_idx]);
448  } else { // edge connects only elements of the same dimension
449  // Allocate the array of sides.
450  last_edge_idx=edges.size();
451  edges.resize(last_edge_idx+1);
452  edg = &( edges.back() );
453  edg->n_sides = 0;
454  edg->side_ = new struct SideIter[ intersection_list.size() ];
455  if (intersection_list.size() > max_edge_sides_[e->dim()-1])
456  max_edge_sides_[e->dim()-1] = intersection_list.size();
457 
458  if (intersection_list.size() == 1) { // outer edge, create boundary object as well
459  edg->n_sides=1;
460  edg->side_[0] = e->side(s);
461  e->edge_idx_[s] = last_edge_idx;
462 
463  if (e->boundary_idx_ == NULL) {
464  e->boundary_idx_ = new unsigned int [ e->n_sides() ];
465  std::fill( e->boundary_idx_, e->boundary_idx_ + e->n_sides(), Mesh::undef_idx);
466  }
467 
468  unsigned int bdr_idx=boundary_.size();
469  boundary_.resize(bdr_idx+1);
470  Boundary &bdr=boundary_.back();
471  e->boundary_idx_[s] = bdr_idx;
472 
473  // fill boundary element
474  ElementFullIter bc_ele = bc_elements.add_item( -bdr_idx ); // use negative bcd index as ID,
475  bc_ele->init(e->dim()-1, this, region_db_.implicit_boundary_region() );
476  for(unsigned int ni = 0; ni< side_nodes.size(); ni++) bc_ele->node[ni] = &( node_vector[side_nodes[ni]] );
477 
478  // fill Boundary object
479  bdr.edge_idx_ = last_edge_idx;
480  bdr.bc_ele_idx_ = bc_ele.index();
481  bdr.mesh_=this;
482 
483  continue; // next side of element e
484  }
485  }
486 
487  // go through the elements connected to the edge or neighbour
488  for( vector<unsigned int>::iterator isect = intersection_list.begin(); isect!=intersection_list.end(); ++isect) {
489  Element *elem = &(element[*isect]);
490  for (unsigned int ecs=0; ecs<elem->n_sides(); ecs++) {
491  if (elem->edge_idx_[ecs] != Mesh::undef_idx) continue;
492  SideIter si = elem->side(ecs);
493  if ( same_sides( si, side_nodes) ) {
494  if (is_neighbour) {
495  // create a new edge and neighbour for this side, and element to the edge
496  last_edge_idx=edges.size();
497  edges.resize(last_edge_idx+1);
498  edg = &( edges.back() );
499  edg->n_sides = 1;
500  edg->side_ = new struct SideIter[1];
501  edg->side_[0] = si;
502  elem->edge_idx_[ecs] = last_edge_idx;
503 
504  neighbour.edge_idx_ = last_edge_idx;
505 
506  vb_neighbours_.push_back(neighbour); // copy neighbour with this edge setting
507  } else {
508  // connect the side to the edge, and side to the edge
509  edg->side_[ edg->n_sides++ ] = si;
510  elem->edge_idx_[ecs] = last_edge_idx;
511  }
512  break; // next element from intersection list
513  }
514  } // search for side of other connected element
515  } // connected elements
516  ASSERT( is_neighbour || ( (unsigned int) edg->n_sides ) == intersection_list.size(), "Some connected sides were not found.\n");
517  } // for element sides
518  } // for elements
519 
520  xprintf( Msg, "Created %d edges and %d neighbours.\n", edges.size(), vb_neighbours_.size() );
521 }
522 
523 
524 
526 {
527  for (EdgeVector::iterator edg=edges.begin(); edg!=edges.end(); edg++)
528  {
529  // side 0 is reference, so its permutation is 0
530  edg->side(0)->element()->permutation_idx_[edg->side(0)->el_idx()] = 0;
531 
532  if (edg->n_sides > 1)
533  {
534  map<const Node*,unsigned int> node_numbers;
535  unsigned int permutation[edg->side(0)->n_nodes()];
536 
537  for (unsigned int i=0; i<edg->side(0)->n_nodes(); i++)
538  node_numbers[edg->side(0)->node(i)] = i;
539 
540  for (int sid=1; sid<edg->n_sides; sid++)
541  {
542  for (unsigned int i=0; i<edg->side(0)->n_nodes(); i++)
543  permutation[node_numbers[edg->side(sid)->node(i)]] = i;
544 
545  switch (edg->side(0)->dim())
546  {
547  case 0:
548  edg->side(sid)->element()->permutation_idx_[edg->side(sid)->el_idx()] = RefElement<1>::permutation_index(permutation);
549  break;
550  case 1:
551  edg->side(sid)->element()->permutation_idx_[edg->side(sid)->el_idx()] = RefElement<2>::permutation_index(permutation);
552  break;
553  case 2:
554  edg->side(sid)->element()->permutation_idx_[edg->side(sid)->el_idx()] = RefElement<3>::permutation_index(permutation);
555  break;
556  }
557  }
558  }
559  }
560 
561  for (vector<Neighbour>::iterator nb=vb_neighbours_.begin(); nb!=vb_neighbours_.end(); nb++)
562  {
563  map<const Node*,unsigned int> node_numbers;
564  unsigned int permutation[nb->element()->n_nodes()];
565 
566  // element of lower dimension is reference, so
567  // we calculate permutation for the adjacent side
568  for (unsigned int i=0; i<nb->element()->n_nodes(); i++)
569  node_numbers[nb->element()->node[i]] = i;
570 
571  for (unsigned int i=0; i<nb->side()->n_nodes(); i++)
572  permutation[node_numbers[nb->side()->node(i)]] = i;
573 
574  switch (nb->side()->dim())
575  {
576  case 0:
577  nb->side()->element()->permutation_idx_[nb->side()->el_idx()] = RefElement<1>::permutation_index(permutation);
578  break;
579  case 1:
580  nb->side()->element()->permutation_idx_[nb->side()->el_idx()] = RefElement<2>::permutation_index(permutation);
581  break;
582  case 2:
583  nb->side()->element()->permutation_idx_[nb->side()->el_idx()] = RefElement<3>::permutation_index(permutation);
584  break;
585  }
586  }
587 }
588 
589 
590 
591 
592 
593 //=============================================================================
594 //
595 //=============================================================================
597 {
598 
599  xprintf( MsgVerb, " Element to neighbours of vb2 type... ")/*orig verb 5*/;
600 
601  FOR_ELEMENTS(this,ele) ele->n_neighs_vb =0;
602 
603  // count vb neighs per element
604  FOR_NEIGHBOURS(this, ngh ) ngh->element_->n_neighs_vb++;
605 
606  // Allocation of the array per element
607  FOR_ELEMENTS(this, ele )
608  if( ele->n_neighs_vb > 0 ) {
609  ele->neigh_vb = new struct Neighbour* [ele->n_neighs_vb];
610  ele->n_neighs_vb=0;
611  }
612 
613  // fill
614  ElementIter ele;
615  FOR_NEIGHBOURS(this, ngh ) {
616  ele = ngh->element();
617  ele->neigh_vb[ ele->n_neighs_vb++ ] = &( *ngh );
618  }
619 
620  xprintf( MsgVerb, "O.K.\n")/*orig verb 6*/;
621 }
622 
623 
624 
625 
629 
630 
632  /* Algorithm:
633  *
634  * 1) create BIH tree
635  * 2) for every 1D, find list of candidates
636  * 3) compute intersections for 1d, store it to master_elements
637  *
638  */
639  BIHTree bih_tree( this );
640  master_elements.resize(n_elements());
641 
642  for(unsigned int i_ele=0; i_ele<n_elements(); i_ele++) {
643  Element &ele = this->element[i_ele];
644 
645  if (ele.dim() == 1) {
646  vector<unsigned int> candidate_list;
647  bih_tree.find_bounding_box(ele.bounding_box(), candidate_list);
648 
649  //for(unsigned int i_elm=0; i_elm<n_elements(); i_elm++) {
650  for(unsigned int i_elm : candidate_list) {
651  ElementFullIter elm = this->element( i_elm );
652  if (elm->dim() == 2) {
653  IntersectionLocal *intersection;
654  GetIntersection( TAbscissa(ele), TTriangle(*elm), intersection);
655  if (intersection && intersection->get_type() == IntersectionLocal::line) {
656 
657  master_elements[i_ele].push_back( intersections.size() );
658  intersections.push_back( Intersection(this->element(i_ele), elm, intersection) );
659  }
660  }
661 
662  }
663  }
664  }
665 
666 }
667 
668 
669 
670 ElementAccessor<3> Mesh::element_accessor(unsigned int idx, bool boundary) {
671  return ElementAccessor<3>(this, idx, boundary);
672 }
673 
674 
675 
676 vector<int> const & Mesh::elements_id_maps( bool boundary_domain) const
677 {
678  if (bulk_elements_id_.size() ==0) {
680  int last_id;
681 
682  bulk_elements_id_.resize(n_elements());
683  map_it = bulk_elements_id_.begin();
684  last_id = -1;
685  for(unsigned int idx=0; idx < element.size(); idx++, ++map_it) {
686  int id = element.get_id(idx);
687  if (last_id >= id) xprintf(UsrErr, "Element IDs in non-increasing order, ID: %d\n", id);
688  last_id=*map_it = id;
689  }
690 
692  map_it = boundary_elements_id_.begin();
693  last_id = -1;
694  for(unsigned int idx=0; idx < bc_elements.size(); idx++, ++map_it) {
695  int id = bc_elements.get_id(idx);
696  // We set ID for boundary elements created by the mesh itself to "-1"
697  // this force gmsh reader to skip all remaining entries in boundary_elements_id_
698  // and thus report error for any remaining data lines
699  if (id < 0) last_id=*map_it=-1;
700  else {
701  if (last_id >= id) xprintf(UsrErr, "Element IDs in non-increasing order, ID: %d\n", id);
702  last_id=*map_it = id;
703  }
704  }
705  }
706 
707  if (boundary_domain) return boundary_elements_id_;
708  return bulk_elements_id_;
709 }
710 
711 //-----------------------------------------------------------------------------
712 // vim: set cindent:
int n_triangles
Definition: mesh.h:242
Class for the mesh partitioning. This should provide:
Definition: partitioning.hh:29
Mesh(const std::string &input_str="{mesh_file=\"\"}", MPI_Comm com=MPI_COMM_WORLD)
Definition: mesh.cc:85
void read_regions_from_input(Input::Array region_list, MapElementIDToRegionID &map)
Definition: region.cc:455
#define FOR_ELEMENT_NODES(i, j)
Definition: elements.h:150
vector< vector< unsigned int > > node_elements
Definition: mesh.h:327
Accessor to input data conforming to declared Array.
Definition: accessors.hh:558
unsigned int * boundary_idx_
Definition: elements.h:94
int get_id(const T *it) const
Definition: sys_vector.hh:468
Definition: system.hh:72
void read_sets_from_input(Input::Array arr)
Definition: region.cc:370
void count_side_types()
Definition: mesh.cc:272
int MPI_Comm
Definition: mpi.h:141
void find_bounding_box(const BoundingBox &boundingBox, std::vector< unsigned int > &result_list)
Definition: bih_tree.cc:203
void make_intersec_elements()
Definition: mesh.cc:631
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:39
#define FOR_ELEMENTS(_mesh_, __i)
Definition: mesh.h:364
unsigned int * permutation_idx_
Definition: elements.h:103
int n_lines
Definition: mesh.h:241
static const unsigned int undef_idx
Definition: mesh.h:111
void create_node_element_lists()
Definition: mesh.cc:285
static unsigned int permutation_index(unsigned int p[n_nodes_per_side])
Definition: ref_element.cc:146
static Default obligatory()
Definition: type_record.hh:89
???
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:330
boost::shared_ptr< Partitioning > part_
Definition: mesh.h:340
int index() const
Definition: sys_vector.hh:88
Input::Record in_record_
Definition: mesh.h:344
Node ** node
Definition: elements.h:90
FullIter add_item(int id)
Definition: sys_vector.hh:369
int n_sides
Definition: edges.h:48
vector< vector< vector< unsigned int > > > side_nodes
Definition: mesh.h:254
Definition: edges.h:38
int n_tetrahedras
Definition: mesh.h:243
ElementAccessor< 3 > element_accessor(unsigned int idx, bool boundary=false)
Definition: mesh.cc:670
vector< Boundary > boundary_
Definition: mesh.h:209
I/O functions with filename storing, able to track current line in opened file. All standard stdio fu...
Partitioning * get_part()
Definition: mesh.cc:191
unsigned int n_sides()
Definition: mesh.cc:172
Class for declaration of inputs sequences.
Definition: type_base.hh:239
int n_exsides
Definition: mesh.h:238
unsigned int dim() const
unsigned int size() const
Returns size of the container. This is independent of the allocated space.
Definition: sys_vector.hh:401
static Default optional()
Definition: type_record.hh:102
bool opt_val(const string &key, Ret &value) const
unsigned int edge_idx_
Definition: neighbours.h:147
unsigned int * edge_idx_
Definition: elements.h:93
unsigned int n_elements() const
Definition: mesh.h:141
bool same_sides(const SideIter &si, vector< unsigned int > &side_nodes)
Definition: mesh.cc:346
#define ASSERT(...)
Definition: global_defs.h:121
Definition: system.hh:72
static FileName input()
Definition: type_base.hh:464
static Input::Type::Record input_type
Definition: partitioning.hh:35
std::vector< T >::iterator iterator
Definition: sys_vector.hh:225
vector< int > const & elements_id_maps(bool boundary_domain) const
Definition: mesh.cc:676
void setup_topology()
Definition: mesh.cc:252
Neighbour ** neigh_vb
Definition: elements.h:129
Accessor to the data with type Type::Record.
Definition: accessors.hh:327
#define NDEF
Definition: mesh.cc:64
const Ret val(const string &key) const
unsigned int n_sides() const
#define xprintf(...)
Definition: system.hh:100
#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:46
unsigned int n_corners()
Definition: mesh.cc:181
unsigned int index(const T *pointer) const
Definition: sys_vector.hh:383
#define FOR_NEIGHBOURS(_mesh_, it)
Definition: mesh.h:415
ElementVector bc_elements
Definition: mesh.h:213
void close()
Definition: region.cc:220
SideIter side(const unsigned int loc_index)
FullIter begin()
Definition: sys_vector.hh:393
const Record & close() const
Definition: type_record.cc:314
Region implicit_boundary_region()
Definition: region.cc:87
unsigned int edge_idx_
Definition: boundaries.h:87
vector< int > boundary_elements_id_
Definition: mesh.h:322
unsigned int bc_ele_idx_
Definition: boundaries.h:88
void count_element_types()
Definition: mesh.cc:200
Dedicated class for storing path to input and output files.
Definition: file_path.hh:32
Definition: system.hh:72
void read_gmsh_from_stream(istream &in)
Definition: mesh.cc:216
vector< vector< unsigned int > > master_elements
Definition: mesh.h:230
vector< Neighbour > vb_neighbours_
Definition: mesh.h:235
void read_mesh(Mesh *mesh, const RegionDB::MapElementIDToRegionID *el_to_reg_map=NULL)
int n_sides_
Definition: mesh.h:239
static Input::Type::Record region_input_type
Definition: region.hh:294
static Input::Type::Record input_type
Definition: mesh.h:112
RegionDB region_db_
Definition: mesh.h:336
std::vector< Edge > edges
Vector of MH edges, this should not be part of the geometrical mesh.
Definition: mesh.h:216
vector< Intersection > intersections
Definition: mesh.h:224
unsigned int n_neighs_vb
Definition: elements.h:127
int n_insides
Definition: mesh.h:237
void intersect_element_lists(vector< unsigned int > const &nodes_list, vector< unsigned int > &intersection_element_list)
Definition: mesh.cc:298
void make_edge_permutations()
Definition: mesh.cc:525
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc...
Mesh * mesh_
Definition: boundaries.h:89
unsigned int n_nodes() const
Definition: mesh.h:137
Record type proxy class.
Definition: type_record.hh:169
ElementIter element_
Definition: neighbours.h:148
SideIter * side_
Definition: edges.h:49
void element_to_neigh_vb()
Definition: mesh.cc:596
static Input::Type::Record region_set_input_type
Definition: region.hh:298
const Node * node(unsigned int i) const
Definition: side_impl.hh:35
void reinit(Input::Record in_record)
Definition: mesh.cc:106
Reader for (slightly) modified JSON files.
FullIter end()
Returns FullFullIterer of the fictions past the end element.
Definition: sys_vector.hh:397
#define FOR_SIDES(_mesh_, it)
Definition: mesh.h:408
void init_from_input()
Definition: mesh.cc:227
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:326
vector< int > bulk_elements_id_
Definition: mesh.h:322
IntersectionType get_type() const
void make_neighbours_and_edges()
Definition: mesh.cc:363
NodeVector node_vector
Vector of nodes of the mesh.
Definition: mesh.h:203
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:205
unsigned int n_nodes() const
Definition: side_impl.hh:22
~Mesh()
Destructor.
Definition: mesh.cc:151
BoundingBox bounding_box()
Definition: elements.h:113
Record & declare_key(const string &key, const KeyType &type, const Default &default_value, const string &description)
Definition: type_record.cc:430