Flow123d  jenkins-Flow123d-windows32-release-multijob-51
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 
151 unsigned int Mesh::n_sides()
152 {
153  if (n_sides_ == NDEF) {
154  n_sides_=0;
155  FOR_ELEMENTS(this, ele) n_sides_ += ele->n_sides();
156  }
157  return n_sides_;
158 }
159 
160 
161 
163  return part_.get();
164 }
165 
166 
167 //=============================================================================
168 // COUNT ELEMENT TYPES
169 //=============================================================================
170 
172  FOR_ELEMENTS(this, elm)
173  switch (elm->dim()) {
174  case 1:
175  n_lines++;
176  break;
177  case 2:
178  n_triangles++;
179  break;
180  case 3:
181  n_tetrahedras++;
182  break;
183  }
184 }
185 
186 
187 void Mesh::read_gmsh_from_stream(istream &in) {
188 
189  START_TIMER("Reading mesh - from_stream");
190 
191  GmshMeshReader reader(in);
192  reader.read_mesh(this);
193  setup_topology();
194 }
195 
196 
197 
199  START_TIMER("Reading mesh - init_from_input");
200 
201  Input::Array region_list;
202  RegionDB::MapElementIDToRegionID el_to_reg_map;
203 
204  // create regions from our input
205  if (in_record_.opt_val("regions", region_list)) {
206  region_db_.read_regions_from_input(region_list, el_to_reg_map);
207  }
208  // read raw mesh, add regions from GMSH file
209  GmshMeshReader reader( in_record_.val<FilePath>("mesh_file") );
210  reader.read_mesh(this, &el_to_reg_map);
211  // possibly add implicit_boundary region, close region_db_.
212  setup_topology();
213  // create sets
214  Input::Array set_list;
215  if (in_record_.opt_val("sets", set_list)) {
217  }
218 }
219 
220 
221 
222 
224  START_TIMER("MESH - setup topology");
225 
227 
228  // check mesh quality
229  FOR_ELEMENTS(this, ele)
230  if (ele->quality_measure_smooth() < 0.001) xprintf(Warn, "Bad quality (<0.001) of the element %u.\n", ele.id());
231 
236 
237  region_db_.close();
238  part_ = boost::make_shared<Partitioning>(this, in_record_.val<Input::Record>("partitioning") );
239 }
240 
241 
242 //
244 {
245 
246  n_insides = 0;
247  n_exsides = 0;
248  FOR_SIDES(this, sde ) {
249  if (sde->is_external()) n_exsides++;
250  else n_insides++;
251  }
252 }
253 
254 
255 
257  // for each node we make a list of elements that use this node
258  node_elements.resize(node_vector.size());
259 
260  FOR_ELEMENTS( this, e )
261  for (unsigned int n=0; n<e->n_nodes(); n++)
262  node_elements[node_vector.index(e->node[n])].push_back(e->index());
263 
264  for (vector<vector<unsigned int> >::iterator n=node_elements.begin(); n!=node_elements.end(); n++)
265  stable_sort(n->begin(), n->end());
266 }
267 
268 
269 void Mesh::intersect_element_lists(vector<unsigned int> const &nodes_list, vector<unsigned int> &intersection_element_list)
270 {
271  if (nodes_list.size() == 0) {
272  intersection_element_list.clear();
273  } else if (nodes_list.size() == 1) {
274  intersection_element_list = node_elements[ nodes_list[0] ];
275  } else {
276  vector<unsigned int>::const_iterator it1=nodes_list.begin();
278  intersection_element_list.resize( node_elements[*it1].size() ); // make enough space
279 
280  it1=set_intersection(
281  node_elements[*it1].begin(), node_elements[*it1].end(),
282  node_elements[*it2].begin(), node_elements[*it2].end(),
283  intersection_element_list.begin());
284  intersection_element_list.resize(it1-intersection_element_list.begin()); // resize to true size
285 
286  for(;it2<nodes_list.end();++it2) {
287  it1=set_intersection(
288  intersection_element_list.begin(), intersection_element_list.end(),
289  node_elements[*it2].begin(), node_elements[*it2].end(),
290  intersection_element_list.begin());
291  intersection_element_list.resize(it1-intersection_element_list.begin()); // resize to true size
292  }
293  }
294 }
295 
296 
298  unsigned int dim, unsigned int &element_idx) {
299  bool is_neighbour = false;
300 
301  vector<unsigned int>::iterator e_dest=element_list.begin();
302  for( vector<unsigned int>::iterator ele = element_list.begin(); ele!=element_list.end(); ++ele)
303  if (elements[*ele].dim() == dim) { // keep only indexes of elements of same dimension
304  *e_dest=*ele;
305  ++e_dest;
306  } else if (elements[*ele].dim() == dim-1) { // get only first element of lower dimension
307  if (is_neighbour) xprintf(UsrErr, "Too matching elements id: %d and id: %d in the same mesh.\n",
308  elements(*ele).id(), elements(element_idx).id() );
309 
310  is_neighbour = true;
311  element_idx = *ele;
312  }
313  element_list.resize( e_dest - element_list.begin());
314  return is_neighbour;
315 }
316 
317 bool Mesh::same_sides(const SideIter &si, vector<unsigned int> &side_nodes) {
318  // 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)
319  unsigned int ni=0;
320  while ( ni < si->n_nodes()
321  && find(side_nodes.begin(), side_nodes.end(), node_vector.index( si->node(ni) ) ) != side_nodes.end() ) ni++;
322  return ( ni == si->n_nodes() );
323 }
324 
325 /**
326  * TODO:
327  * - use std::is_any for setting is_neigbour
328  * - possibly make appropriate constructors for Edge and Neighbour
329  * - check side!=-1 when searching neigbouring element
330  * - process bc_elements first, there should be no Neigh, but check it
331  * set Edge and boundary there
332  */
333 
335 {
336  Neighbour neighbour;
337  Edge *edg;
338  unsigned int ngh_element_idx, last_edge_idx;
339 
341 
342  // pointers to created edges
343  //vector<Edge *> tmp_edges;
344  edges.resize(0); // be sure that edges are empty
345 
347  vector<unsigned int> intersection_list; // list of elements in intersection of node element lists
348 
349  for( ElementFullIter bc_ele = bc_elements.begin(); bc_ele != bc_elements.end(); ++bc_ele) {
350  // Find all elements that share this side.
351  side_nodes.resize(bc_ele->n_nodes());
352  for (unsigned n=0; n<bc_ele->n_nodes(); n++) side_nodes[n] = node_vector.index(bc_ele->node[n]);
353  intersect_element_lists(side_nodes, intersection_list);
354  bool is_neighbour = find_lower_dim_element(element, intersection_list, bc_ele->dim() +1, ngh_element_idx);
355  if (is_neighbour) {
356  xprintf(UsrErr, "Boundary element (id: %d) match a regular element (id: %d) of lower dimension.\n",
357  bc_ele.id(), element(ngh_element_idx).id());
358  } else {
359  if (intersection_list.size() == 0) {
360  // no matching dim+1 element found
361  xprintf(Warn, "Lonely boundary element, id: %d, region: %d, dimension %d.\n", bc_ele.id(), bc_ele->region().id(), bc_ele->dim());
362  continue; // skip the boundary element
363  }
364  last_edge_idx=edges.size();
365  edges.resize(last_edge_idx+1);
366  edg = &( edges.back() );
367  edg->n_sides = 0;
368  edg->side_ = new struct SideIter[ intersection_list.size() ];
369 
370  // common boundary object
371  unsigned int bdr_idx=boundary_.size();
372  boundary_.resize(bdr_idx+1);
373  Boundary &bdr=boundary_.back();
374  bdr.bc_ele_idx_ = bc_ele.index();
375  bdr.edge_idx_ = last_edge_idx;
376  bdr.mesh_=this;
377 
378  // for 1d boundaries there can be more then one 1d elements connected to the boundary element
379  // we do not detect this case later in the main search over bulk elements
380  for( vector<unsigned int>::iterator isect = intersection_list.begin(); isect!=intersection_list.end(); ++isect) {
381  Element *elem = &(element[*isect]);
382  for (unsigned int ecs=0; ecs<elem->n_sides(); ecs++) {
383  SideIter si = elem->side(ecs);
384  if ( same_sides( si, side_nodes) ) {
385  edg->side_[ edg->n_sides++ ] = si;
386  elem->edge_idx_[ecs] = last_edge_idx;
387 
388  if (elem->boundary_idx_ == NULL) {
389  elem->boundary_idx_ = new unsigned int [ elem->n_sides() ];
390  std::fill( elem->boundary_idx_, elem->boundary_idx_ + elem->n_sides(), Mesh::undef_idx);
391  }
392  elem->boundary_idx_[ecs] = bdr_idx;
393  break; // next element in intersection list
394  }
395  }
396  }
397 
398  }
399 
400  }
401  // Now we go through all element sides and create edges and neighbours
402  FOR_ELEMENTS( this, e )
403  {
404  for (unsigned int s=0; s<e->n_sides(); s++)
405  {
406  // skip sides that were already found
407  if (e->edge_idx_[s] != Mesh::undef_idx) continue;
408 
409 
410  // Find all elements that share this side.
411  side_nodes.resize(e->side(s)->n_nodes());
412  for (unsigned n=0; n<e->side(s)->n_nodes(); n++) side_nodes[n] = node_vector.index(e->side(s)->node(n));
413  intersect_element_lists(side_nodes, intersection_list);
414 
415  bool is_neighbour = find_lower_dim_element(element, intersection_list, e->dim(), ngh_element_idx);
416 
417  if (is_neighbour) { // edge connects elements of different dimensions
418  neighbour.element_ = &(element[ngh_element_idx]);
419  } else { // edge connects only elements of the same dimension
420  // Allocate the array of sides.
421  last_edge_idx=edges.size();
422  edges.resize(last_edge_idx+1);
423  edg = &( edges.back() );
424  edg->n_sides = 0;
425  edg->side_ = new struct SideIter[ intersection_list.size() ];
426  if (intersection_list.size() > max_edge_sides_[e->dim()-1])
427  max_edge_sides_[e->dim()-1] = intersection_list.size();
428 
429  if (intersection_list.size() == 1) { // outer edge, create boundary object as well
430  edg->n_sides=1;
431  edg->side_[0] = e->side(s);
432  e->edge_idx_[s] = last_edge_idx;
433 
434  if (e->boundary_idx_ == NULL) {
435  e->boundary_idx_ = new unsigned int [ e->n_sides() ];
436  std::fill( e->boundary_idx_, e->boundary_idx_ + e->n_sides(), Mesh::undef_idx);
437  }
438 
439  unsigned int bdr_idx=boundary_.size();
440  boundary_.resize(bdr_idx+1);
441  Boundary &bdr=boundary_.back();
442  e->boundary_idx_[s] = bdr_idx;
443 
444  // fill boundary element
445  ElementFullIter bc_ele = bc_elements.add_item( -bdr_idx ); // use negative bcd index as ID,
446  bc_ele->init(e->dim()-1, this, region_db_.implicit_boundary_region() );
447  for(unsigned int ni = 0; ni< side_nodes.size(); ni++) bc_ele->node[ni] = &( node_vector[side_nodes[ni]] );
448 
449  // fill Boundary object
450  bdr.edge_idx_ = last_edge_idx;
451  bdr.bc_ele_idx_ = bc_ele.index();
452  bdr.mesh_=this;
453 
454  continue; // next side of element e
455  }
456  }
457 
458  // go through the elements connected to the edge or neighbour
459  for( vector<unsigned int>::iterator isect = intersection_list.begin(); isect!=intersection_list.end(); ++isect) {
460  Element *elem = &(element[*isect]);
461  for (unsigned int ecs=0; ecs<elem->n_sides(); ecs++) {
462  if (elem->edge_idx_[ecs] != Mesh::undef_idx) continue;
463  SideIter si = elem->side(ecs);
464  if ( same_sides( si, side_nodes) ) {
465  if (is_neighbour) {
466  // create a new edge and neighbour for this side, and element to the edge
467  last_edge_idx=edges.size();
468  edges.resize(last_edge_idx+1);
469  edg = &( edges.back() );
470  edg->n_sides = 1;
471  edg->side_ = new struct SideIter[1];
472  edg->side_[0] = si;
473  elem->edge_idx_[ecs] = last_edge_idx;
474 
475  neighbour.edge_idx_ = last_edge_idx;
476 
477  vb_neighbours_.push_back(neighbour); // copy neighbour with this edge setting
478  } else {
479  // connect the side to the edge, and side to the edge
480  edg->side_[ edg->n_sides++ ] = si;
481  elem->edge_idx_[ecs] = last_edge_idx;
482  }
483  break; // next element from intersection list
484  }
485  } // search for side of other connected element
486  } // connected elements
487  ASSERT( is_neighbour || ( (unsigned int) edg->n_sides ) == intersection_list.size(), "Some connected sides were not found.\n");
488  } // for element sides
489  } // for elements
490 
491  xprintf( Msg, "Created %d edges and %d neighbours.\n", edges.size(), vb_neighbours_.size() );
492 }
493 
494 
495 
497 {
498  for (EdgeVector::iterator edg=edges.begin(); edg!=edges.end(); edg++)
499  {
500  // side 0 is reference, so its permutation is 0
501  edg->side(0)->element()->permutation_idx_[edg->side(0)->el_idx()] = 0;
502 
503  if (edg->n_sides > 1)
504  {
505  map<const Node*,unsigned int> node_numbers;
506  unsigned int permutation[edg->side(0)->n_nodes()];
507 
508  for (unsigned int i=0; i<edg->side(0)->n_nodes(); i++)
509  node_numbers[edg->side(0)->node(i)] = i;
510 
511  for (int sid=1; sid<edg->n_sides; sid++)
512  {
513  for (unsigned int i=0; i<edg->side(0)->n_nodes(); i++)
514  permutation[node_numbers[edg->side(sid)->node(i)]] = i;
515 
516  switch (edg->side(0)->dim())
517  {
518  case 0:
519  edg->side(sid)->element()->permutation_idx_[edg->side(sid)->el_idx()] = RefElement<1>::permutation_index(permutation);
520  break;
521  case 1:
522  edg->side(sid)->element()->permutation_idx_[edg->side(sid)->el_idx()] = RefElement<2>::permutation_index(permutation);
523  break;
524  case 2:
525  edg->side(sid)->element()->permutation_idx_[edg->side(sid)->el_idx()] = RefElement<3>::permutation_index(permutation);
526  break;
527  }
528  }
529  }
530  }
531 
532  for (vector<Neighbour>::iterator nb=vb_neighbours_.begin(); nb!=vb_neighbours_.end(); nb++)
533  {
534  map<const Node*,unsigned int> node_numbers;
535  unsigned int permutation[nb->element()->n_nodes()];
536 
537  // element of lower dimension is reference, so
538  // we calculate permutation for the adjacent side
539  for (unsigned int i=0; i<nb->element()->n_nodes(); i++)
540  node_numbers[nb->element()->node[i]] = i;
541 
542  for (unsigned int i=0; i<nb->side()->n_nodes(); i++)
543  permutation[node_numbers[nb->side()->node(i)]] = i;
544 
545  switch (nb->side()->dim())
546  {
547  case 0:
548  nb->side()->element()->permutation_idx_[nb->side()->el_idx()] = RefElement<1>::permutation_index(permutation);
549  break;
550  case 1:
551  nb->side()->element()->permutation_idx_[nb->side()->el_idx()] = RefElement<2>::permutation_index(permutation);
552  break;
553  case 2:
554  nb->side()->element()->permutation_idx_[nb->side()->el_idx()] = RefElement<3>::permutation_index(permutation);
555  break;
556  }
557  }
558 }
559 
560 
561 
562 
563 
564 //=============================================================================
565 //
566 //=============================================================================
568 {
569 
570  xprintf( MsgVerb, " Element to neighbours of vb2 type... ")/*orig verb 5*/;
571 
572  FOR_ELEMENTS(this,ele) ele->n_neighs_vb =0;
573 
574  // count vb neighs per element
575  FOR_NEIGHBOURS(this, ngh ) ngh->element_->n_neighs_vb++;
576 
577  // Allocation of the array per element
578  FOR_ELEMENTS(this, ele )
579  if( ele->n_neighs_vb > 0 ) {
580  ele->neigh_vb = new struct Neighbour* [ele->n_neighs_vb];
581  ele->n_neighs_vb=0;
582  }
583 
584  // fill
585  ElementIter ele;
586  FOR_NEIGHBOURS(this, ngh ) {
587  ele = ngh->element();
588  ele->neigh_vb[ ele->n_neighs_vb++ ] = &( *ngh );
589  }
590 
591  xprintf( MsgVerb, "O.K.\n")/*orig verb 6*/;
592 }
593 
594 
595 
596 
600 
601 
603  /* Algorithm:
604  *
605  * 1) create BIH tree
606  * 2) for every 1D, find list of candidates
607  * 3) compute intersections for 1d, store it to master_elements
608  *
609  */
610  BIHTree bih_tree( this );
611  master_elements.resize(n_elements());
612 
613  for(unsigned int i_ele=0; i_ele<n_elements(); i_ele++) {
614  Element &ele = this->element[i_ele];
615 
616  if (ele.dim() == 1) {
617  vector<unsigned int> candidate_list;
618  for(unsigned int i_elm=0; i_elm<n_elements(); i_elm++) {
619  ElementFullIter elm = this->element( i_elm );
620  if (elm->dim() == 2) {
621  IntersectionLocal *intersection;
622  GetIntersection( TAbscissa(ele), TTriangle(*elm), intersection);
623  if (intersection && intersection->get_type() == IntersectionLocal::line) {
624 
625  master_elements[i_ele].push_back( intersections.size() );
626  intersections.push_back( Intersection(this->element(i_ele), elm, intersection) );
627  }
628  }
629 
630  }
631  }
632  }
633 
634 }
635 
636 
637 
638 ElementAccessor<3> Mesh::element_accessor(unsigned int idx, bool boundary) {
639  return ElementAccessor<3>(this, idx, boundary);
640 }
641 
642 
643 
644 vector<int> const & Mesh::elements_id_maps( bool boundary_domain) const
645 {
646  if (bulk_elements_id_.size() ==0) {
648  int last_id;
649 
650  bulk_elements_id_.resize(n_elements());
651  map_it = bulk_elements_id_.begin();
652  last_id = -1;
653  for(unsigned int idx=0; idx < element.size(); idx++, ++map_it) {
654  int id = element.get_id(idx);
655  if (last_id >= id) xprintf(UsrErr, "Element IDs in non-increasing order, ID: %d\n", id);
656  last_id=*map_it = id;
657  }
658 
660  map_it = boundary_elements_id_.begin();
661  last_id = -1;
662  for(unsigned int idx=0; idx < bc_elements.size(); idx++, ++map_it) {
663  int id = bc_elements.get_id(idx);
664  // We set ID for boundary elements created by the mesh itself to "-1"
665  // this force gmsh reader to skip all remaining entries in boundary_elements_id_
666  // and thus report error for any remaining data lines
667  if (id < 0) last_id=*map_it=-1;
668  else {
669  if (last_id >= id) xprintf(UsrErr, "Element IDs in non-increasing order, ID: %d\n", id);
670  last_id=*map_it = id;
671  }
672  }
673  }
674 
675  if (boundary_domain) return boundary_elements_id_;
676  return bulk_elements_id_;
677 }
678 
679 //-----------------------------------------------------------------------------
680 // vim: set cindent:
int n_triangles
Definition: mesh.h:235
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
vector< vector< unsigned int > > node_elements
Definition: mesh.h:320
Accessor to input data conforming to declared Array.
Definition: accessors.hh:521
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:243
int MPI_Comm
Definition: mpi.h:141
void make_intersec_elements()
Definition: mesh.cc:602
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:357
int n_lines
Definition: mesh.h:234
static const unsigned int undef_idx
Definition: mesh.h:110
void create_node_element_lists()
Definition: mesh.cc:256
static unsigned int permutation_index(unsigned int p[n_nodes_per_side])
Definition: ref_element.cc:146
static Default obligatory()
Definition: type_record.hh:87
???
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:323
boost::shared_ptr< Partitioning > part_
Definition: mesh.h:333
int index() const
Definition: sys_vector.hh:88
Input::Record in_record_
Definition: mesh.h:337
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:247
Definition: edges.h:38
int n_tetrahedras
Definition: mesh.h:236
ElementAccessor< 3 > element_accessor(unsigned int idx, bool boundary=false)
Definition: mesh.cc:638
vector< Boundary > boundary_
Definition: mesh.h:202
I/O functions with filename storing, able to track current line in opened file. All standard stdio fu...
Partitioning * get_part()
Definition: mesh.cc:162
unsigned int n_sides()
Definition: mesh.cc:151
Class for declaration of inputs sequences.
Definition: type_base.hh:230
int n_exsides
Definition: mesh.h:231
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:100
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:137
bool same_sides(const SideIter &si, vector< unsigned int > &side_nodes)
Definition: mesh.cc:317
#define ASSERT(...)
Definition: global_defs.h:121
Definition: system.hh:72
static FileName input()
Definition: type_base.hh:443
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:644
void setup_topology()
Definition: mesh.cc:223
Neighbour ** neigh_vb
Definition: elements.h:129
Accessor to the data with type Type::Record.
Definition: accessors.hh:308
#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 index(const T *pointer) const
Definition: sys_vector.hh:383
#define FOR_NEIGHBOURS(_mesh_, it)
Definition: mesh.h:408
ElementVector bc_elements
Definition: mesh.h:206
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:286
Region implicit_boundary_region()
Definition: region.cc:87
unsigned int edge_idx_
Definition: boundaries.h:87
vector< int > boundary_elements_id_
Definition: mesh.h:315
unsigned int bc_ele_idx_
Definition: boundaries.h:88
void count_element_types()
Definition: mesh.cc:171
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:187
vector< vector< unsigned int > > master_elements
Definition: mesh.h:223
vector< Neighbour > vb_neighbours_
Definition: mesh.h:228
void read_mesh(Mesh *mesh, const RegionDB::MapElementIDToRegionID *el_to_reg_map=NULL)
int n_sides_
Definition: mesh.h:232
static Input::Type::Record region_input_type
Definition: region.hh:294
static Input::Type::Record input_type
Definition: mesh.h:111
RegionDB region_db_
Definition: mesh.h:329
std::vector< Edge > edges
Vector of MH edges, this should not be part of the geometrical mesh.
Definition: mesh.h:209
vector< Intersection > intersections
Definition: mesh.h:217
unsigned int n_neighs_vb
Definition: elements.h:127
int n_insides
Definition: mesh.h:230
void intersect_element_lists(vector< unsigned int > const &nodes_list, vector< unsigned int > &intersection_element_list)
Definition: mesh.cc:269
void make_edge_permutations()
Definition: mesh.cc:496
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:133
Record type proxy class.
Definition: type_record.hh:161
ElementIter element_
Definition: neighbours.h:148
SideIter * side_
Definition: edges.h:49
void element_to_neigh_vb()
Definition: mesh.cc:567
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:401
void init_from_input()
Definition: mesh.cc:198
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:297
vector< int > bulk_elements_id_
Definition: mesh.h:315
IntersectionType get_type() const
void make_neighbours_and_edges()
Definition: mesh.cc:334
NodeVector node_vector
Vector of nodes of the mesh.
Definition: mesh.h:196
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:198
unsigned int n_nodes() const
Definition: side_impl.hh:22
Record & declare_key(const string &key, const KeyType &type, const Default &default_value, const string &description)
Definition: type_record.cc:386