Flow123d  jenkins-Flow123d-windows32-release-multijob-28
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_reader.h"
62 #include "mesh/msh_gmshreader.h"
63 #include "mesh/region.hh"
64 
65 #define NDEF -1
66 
67 namespace IT = Input::Type;
68 
69 
71  = IT::Record("Mesh","Record with mesh related data." )
73  "Input file with mesh description.")
75  "List of additional region definitions not contained in the mesh.")
77  "List of region set definitions. There are three region sets implicitly defined:\n"
78  "ALL (all regions of the mesh), BOUNDARY (all boundary regions), and BULK (all bulk regions)")
79  .declare_key("partitioning", Partitioning::input_type, IT::Default("any_neighboring"), "Parameters of mesh partitioning algorithms.\n" )
80  .close();
81 
82 
83 
84 const unsigned int Mesh::undef_idx;
85 
86 Mesh::Mesh(const std::string &input_str, MPI_Comm comm)
87 :comm_(comm)
88 {
89 
90  Input::JSONToStorage reader( input_str, Mesh::input_type );
92 
94 }
95 
96 
97 
99 : in_record_(in_record),
100  comm_(com)
101 {
103 }
104 
105 
106 
108 {
109 
110  n_insides = NDEF;
111  n_exsides = NDEF;
112  n_sides_ = NDEF;
113 
114  // number of element of particular dimension
115  n_lines = 0;
116  n_triangles = 0;
117  n_tetrahedras = 0;
118 
119  for (int d=0; d<3; d++) max_edge_sides_[d] = 0;
120 
121  // Initialize numbering of nodes on sides.
122  // This is temporary solution, until class Element is templated
123  // by dimension. Then we can replace Mesh::side_nodes by
124  // RefElement<dim>::side_nodes.
125 
126  // indices of side nodes in element node array
127  // Currently this is made ad libitum
128  // with some ordering here we can get sides with correct orientation.
129  // This speedup normal calculation.
130 
131  side_nodes.resize(3); // three side dimensions
132  for(int i=0; i < 3; i++) {
133  side_nodes[i].resize(i+2); // number of sides
134  for(int j=0; j < i+2; j++)
135  side_nodes[i][j].resize(i+1);
136  }
137 
138  for (unsigned int sid=0; sid<RefElement<1>::n_sides; sid++)
139  for (unsigned int nid=0; nid<RefElement<1>::n_nodes_per_side; nid++)
140  side_nodes[0][sid][nid] = RefElement<1>::side_nodes[sid][nid];
141 
142  for (unsigned int sid=0; sid<RefElement<2>::n_sides; sid++)
143  for (unsigned int nid=0; nid<RefElement<2>::n_nodes_per_side; nid++)
144  side_nodes[1][sid][nid] = RefElement<2>::side_nodes[sid][nid];
145 
146  for (unsigned int sid=0; sid<RefElement<3>::n_sides; sid++)
147  for (unsigned int nid=0; nid<RefElement<3>::n_nodes_per_side; nid++)
148  side_nodes[2][sid][nid] = RefElement<3>::side_nodes[sid][nid];
149 }
150 
151 
152 unsigned int Mesh::n_sides()
153 {
154  if (n_sides_ == NDEF) {
155  n_sides_=0;
156  FOR_ELEMENTS(this, ele) n_sides_ += ele->n_sides();
157  }
158  return n_sides_;
159 }
160 
161 
162 
164  return part_.get();
165 }
166 
167 
168 //=============================================================================
169 // COUNT ELEMENT TYPES
170 //=============================================================================
171 
173  FOR_ELEMENTS(this, elm)
174  switch (elm->dim()) {
175  case 1:
176  n_lines++;
177  break;
178  case 2:
179  n_triangles++;
180  break;
181  case 3:
182  n_tetrahedras++;
183  break;
184  }
185 }
186 
187 
188 void Mesh::read_gmsh_from_stream(istream &in) {
189 
190  START_TIMER("Reading mesh - from_stream");
191 
192  GmshMeshReader reader(in);
193  reader.read_mesh(this);
194  setup_topology();
195 }
196 
197 
198 
200  START_TIMER("Reading mesh - init_from_input");
201 
202  Input::Array region_list;
203  RegionDB::MapElementIDToRegionID el_to_reg_map;
204 
205  // create regions from our input
206  if (in_record_.opt_val("regions", region_list)) {
207  region_db_.read_regions_from_input(region_list, el_to_reg_map);
208  }
209  // read raw mesh, add regions from GMSH file
210  GmshMeshReader reader( in_record_.val<FilePath>("mesh_file") );
211  reader.read_mesh(this, &el_to_reg_map);
212  // possibly add implicit_boundary region, close region_db_.
213  setup_topology();
214  // create sets
215  Input::Array set_list;
216  if (in_record_.opt_val("sets", set_list)) {
218  }
219 }
220 
221 
222 
223 
225  START_TIMER("MESH - setup topology");
226 
228 
229  // check mesh quality
230  FOR_ELEMENTS(this, ele)
231  if (ele->quality_measure_smooth() < 0.001) xprintf(Warn, "Bad quality (<0.001) of the element %u.\n", ele.id());
232 
237 
238  region_db_.close();
239  part_ = boost::make_shared<Partitioning>(this, in_record_.val<Input::Record>("partitioning") );
240 }
241 
242 
243 //
245 {
246 
247  n_insides = 0;
248  n_exsides = 0;
249  //FOR_SIDES(this, sde ) {
250  // DBGMSG( "ele: %d edge: %d\n", sde->element().index(), sde->edge_idx());
251  //}
252  FOR_SIDES(this, sde ) {
253  if (sde->is_external()) n_exsides++;
254  else n_insides++;
255  }
256 }
257 
258 
259 
261  // for each node we make a list of elements that use this node
262  node_elements.resize(node_vector.size());
263 
264  FOR_ELEMENTS( this, e )
265  for (unsigned int n=0; n<e->n_nodes(); n++)
266  node_elements[node_vector.index(e->node[n])].push_back(e->index());
267 
268  for (vector<vector<unsigned int> >::iterator n=node_elements.begin(); n!=node_elements.end(); n++)
269  stable_sort(n->begin(), n->end());
270 }
271 
272 
273 void Mesh::intersect_element_lists(vector<unsigned int> const &nodes_list, vector<unsigned int> &intersection_element_list)
274 {
275  if (nodes_list.size() == 0) {
276  intersection_element_list.clear();
277  } else if (nodes_list.size() == 1) {
278  intersection_element_list = node_elements[ nodes_list[0] ];
279  } else {
280  vector<unsigned int>::const_iterator it1=nodes_list.begin();
282  intersection_element_list.resize( node_elements[*it1].size() ); // make enough space
283 
284  it1=set_intersection(
285  node_elements[*it1].begin(), node_elements[*it1].end(),
286  node_elements[*it2].begin(), node_elements[*it2].end(),
287  intersection_element_list.begin());
288  intersection_element_list.resize(it1-intersection_element_list.begin()); // resize to true size
289 
290  for(;it2<nodes_list.end();++it2) {
291  it1=set_intersection(
292  intersection_element_list.begin(), intersection_element_list.end(),
293  node_elements[*it2].begin(), node_elements[*it2].end(),
294  intersection_element_list.begin());
295  intersection_element_list.resize(it1-intersection_element_list.begin()); // resize to true size
296  }
297  }
298 }
299 
300 
302  unsigned int dim, unsigned int &element_idx) {
303  bool is_neighbour = false;
304 
305  vector<unsigned int>::iterator e_dest=element_list.begin();
306  for( vector<unsigned int>::iterator ele = element_list.begin(); ele!=element_list.end(); ++ele)
307  if (elements[*ele].dim() == dim) { // keep only indexes of elements of same dimension
308  *e_dest=*ele;
309  ++e_dest;
310  } else if (elements[*ele].dim() == dim-1) { // get only first element of lower dimension
311  if (is_neighbour) xprintf(UsrErr, "Too matching elements id: %d and id: %d in the same mesh.\n",
312  elements(*ele).id(), elements(element_idx).id() );
313 
314  is_neighbour = true;
315  element_idx = *ele;
316  }
317  element_list.resize( e_dest - element_list.begin());
318  return is_neighbour;
319 }
320 
321 bool Mesh::same_sides(const SideIter &si, vector<unsigned int> &side_nodes) {
322  // 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)
323  unsigned int ni=0;
324  while ( ni < si->n_nodes()
325  && find(side_nodes.begin(), side_nodes.end(), node_vector.index( si->node(ni) ) ) != side_nodes.end() ) ni++;
326  return ( ni == si->n_nodes() );
327 }
328 
329 /**
330  * TODO:
331  * - use std::is_any for setting is_neigbour
332  * - possibly make appropriate constructors for Edge and Neighbour
333  * - check side!=-1 when searching neigbouring element
334  * - process bc_elements first, there should be no Neigh, but check it
335  * set Edge and boundary there
336  */
337 
339 {
340  Neighbour neighbour;
341  Edge *edg;
342  unsigned int ngh_element_idx, last_edge_idx;
343 
345 
346  // pointers to created edges
347  //vector<Edge *> tmp_edges;
348  edges.resize(0); // be sure that edges are empty
349 
351  vector<unsigned int> intersection_list; // list of elements in intersection of node element lists
352 
353  for( ElementFullIter bc_ele = bc_elements.begin(); bc_ele != bc_elements.end(); ++bc_ele) {
354  // Find all elements that share this side.
355  side_nodes.resize(bc_ele->n_nodes());
356  for (unsigned n=0; n<bc_ele->n_nodes(); n++) side_nodes[n] = node_vector.index(bc_ele->node[n]);
357  intersect_element_lists(side_nodes, intersection_list);
358  bool is_neighbour = find_lower_dim_element(element, intersection_list, bc_ele->dim() +1, ngh_element_idx);
359  if (is_neighbour) {
360  xprintf(UsrErr, "Boundary element (id: %d) match a regular element (id: %d) of lower dimension.\n",
361  bc_ele.id(), element(ngh_element_idx).id());
362  } else {
363  if (intersection_list.size() == 0) {
364  // no matching dim+1 element found
365  xprintf(Warn, "Lonely boundary element, id: %d, region: %d, dimension %d.\n", bc_ele.id(), bc_ele->region().id(), bc_ele->dim());
366  continue; // skip the boundary element
367  }
368  last_edge_idx=edges.size();
369  edges.resize(last_edge_idx+1);
370  edg = &( edges.back() );
371  edg->n_sides = 0;
372  edg->side_ = new struct SideIter[ intersection_list.size() ];
373 
374  // common boundary object
375  unsigned int bdr_idx=boundary_.size();
376  boundary_.resize(bdr_idx+1);
377  Boundary &bdr=boundary_.back();
378  bdr.bc_ele_idx_ = bc_ele.index();
379  bdr.edge_idx_ = last_edge_idx;
380  bdr.mesh_=this;
381 
382  // for 1d boundaries there can be more then one 1d elements connected to the boundary element
383  // we do not detect this case later in the main search over bulk elements
384  for( vector<unsigned int>::iterator isect = intersection_list.begin(); isect!=intersection_list.end(); ++isect) {
385  Element *elem = &(element[*isect]);
386  for (unsigned int ecs=0; ecs<elem->n_sides(); ecs++) {
387  SideIter si = elem->side(ecs);
388  if ( same_sides( si, side_nodes) ) {
389  edg->side_[ edg->n_sides++ ] = si;
390  elem->edge_idx_[ecs] = last_edge_idx;
391 
392  if (elem->boundary_idx_ == NULL) {
393  elem->boundary_idx_ = new unsigned int [ elem->n_sides() ];
394  std::fill( elem->boundary_idx_, elem->boundary_idx_ + elem->n_sides(), Mesh::undef_idx);
395  }
396  elem->boundary_idx_[ecs] = bdr_idx;
397  break; // next element in intersection list
398  }
399  }
400  }
401 
402  }
403 
404  }
405  // Now we go through all element sides and create edges and neighbours
406  FOR_ELEMENTS( this, e )
407  {
408  for (unsigned int s=0; s<e->n_sides(); s++)
409  {
410  // skip sides that were already found
411  if (e->edge_idx_[s] != Mesh::undef_idx) continue;
412 
413 
414  // Find all elements that share this side.
415  side_nodes.resize(e->side(s)->n_nodes());
416  for (unsigned n=0; n<e->side(s)->n_nodes(); n++) side_nodes[n] = node_vector.index(e->side(s)->node(n));
417  intersect_element_lists(side_nodes, intersection_list);
418 
419  bool is_neighbour = find_lower_dim_element(element, intersection_list, e->dim(), ngh_element_idx);
420 
421  if (is_neighbour) { // edge connects elements of different dimensions
422  neighbour.element_ = &(element[ngh_element_idx]);
423  } else { // edge connects only elements of the same dimension
424  // Allocate the array of sides.
425  last_edge_idx=edges.size();
426  edges.resize(last_edge_idx+1);
427  edg = &( edges.back() );
428  edg->n_sides = 0;
429  edg->side_ = new struct SideIter[ intersection_list.size() ];
430  if (intersection_list.size() > max_edge_sides_[e->dim()-1])
431  max_edge_sides_[e->dim()-1] = intersection_list.size();
432 
433  if (intersection_list.size() == 1) { // outer edge, create boundary object as well
434  edg->n_sides=1;
435  edg->side_[0] = e->side(s);
436  e->edge_idx_[s] = last_edge_idx;
437 
438  if (e->boundary_idx_ == NULL) {
439  e->boundary_idx_ = new unsigned int [ e->n_sides() ];
440  std::fill( e->boundary_idx_, e->boundary_idx_ + e->n_sides(), Mesh::undef_idx);
441  }
442 
443  unsigned int bdr_idx=boundary_.size();
444  boundary_.resize(bdr_idx+1);
445  Boundary &bdr=boundary_.back();
446  e->boundary_idx_[s] = bdr_idx;
447 
448  // fill boundary element
449  ElementFullIter bc_ele = bc_elements.add_item( -bdr_idx ); // use negative bcd index as ID,
450  bc_ele->init(e->dim()-1, this, region_db_.implicit_boundary_region() );
451  for(unsigned int ni = 0; ni< side_nodes.size(); ni++) bc_ele->node[ni] = &( node_vector[side_nodes[ni]] );
452 
453  // fill Boundary object
454  bdr.edge_idx_ = last_edge_idx;
455  bdr.bc_ele_idx_ = bc_ele.index();
456  bdr.mesh_=this;
457 
458  continue; // next side of element e
459  }
460  }
461 
462  // go through the elements connected to the edge or neighbour
463  for( vector<unsigned int>::iterator isect = intersection_list.begin(); isect!=intersection_list.end(); ++isect) {
464  Element *elem = &(element[*isect]);
465  for (unsigned int ecs=0; ecs<elem->n_sides(); ecs++) {
466  if (elem->edge_idx_[ecs] != Mesh::undef_idx) continue;
467  SideIter si = elem->side(ecs);
468  if ( same_sides( si, side_nodes) ) {
469  if (is_neighbour) {
470  // create a new edge and neighbour for this side, and element to the edge
471  last_edge_idx=edges.size();
472  edges.resize(last_edge_idx+1);
473  edg = &( edges.back() );
474  edg->n_sides = 1;
475  edg->side_ = new struct SideIter[1];
476  edg->side_[0] = si;
477  elem->edge_idx_[ecs] = last_edge_idx;
478 
479  neighbour.edge_idx_ = last_edge_idx;
480 
481  vb_neighbours_.push_back(neighbour); // copy neighbour with this edge setting
482  } else {
483  // connect the side to the edge, and side to the edge
484  edg->side_[ edg->n_sides++ ] = si;
485  elem->edge_idx_[ecs] = last_edge_idx;
486  }
487  break; // next element from intersection list
488  }
489  } // search for side of other connected element
490  } // connected elements
491  ASSERT( is_neighbour || ( (unsigned int) edg->n_sides ) == intersection_list.size(), "Some connected sides were not found.\n");
492  } // for element sides
493  } // for elements
494 
495  xprintf( Msg, "Created %d edges and %d neighbours.\n", edges.size(), vb_neighbours_.size() );
496 }
497 
498 
499 
501 {
502  for (EdgeVector::iterator edg=edges.begin(); edg!=edges.end(); edg++)
503  {
504  // side 0 is reference, so its permutation is 0
505  edg->side(0)->element()->permutation_idx_[edg->side(0)->el_idx()] = 0;
506 
507  if (edg->n_sides > 1)
508  {
509  map<const Node*,unsigned int> node_numbers;
510  unsigned int permutation[edg->side(0)->n_nodes()];
511 
512  for (unsigned int i=0; i<edg->side(0)->n_nodes(); i++)
513  node_numbers[edg->side(0)->node(i)] = i;
514 
515  for (int sid=1; sid<edg->n_sides; sid++)
516  {
517  for (unsigned int i=0; i<edg->side(0)->n_nodes(); i++)
518  permutation[node_numbers[edg->side(sid)->node(i)]] = i;
519 
520  switch (edg->side(0)->dim())
521  {
522  case 0:
523  edg->side(sid)->element()->permutation_idx_[edg->side(sid)->el_idx()] = RefElement<1>::permutation_index(permutation);
524  break;
525  case 1:
526  edg->side(sid)->element()->permutation_idx_[edg->side(sid)->el_idx()] = RefElement<2>::permutation_index(permutation);
527  break;
528  case 2:
529  edg->side(sid)->element()->permutation_idx_[edg->side(sid)->el_idx()] = RefElement<3>::permutation_index(permutation);
530  break;
531  }
532  }
533  }
534  }
535 
536  for (vector<Neighbour>::iterator nb=vb_neighbours_.begin(); nb!=vb_neighbours_.end(); nb++)
537  {
538  map<const Node*,unsigned int> node_numbers;
539  unsigned int permutation[nb->element()->n_nodes()];
540 
541  // element of lower dimension is reference, so
542  // we calculate permutation for the adjacent side
543  for (unsigned int i=0; i<nb->element()->n_nodes(); i++)
544  node_numbers[nb->element()->node[i]] = i;
545 
546  for (unsigned int i=0; i<nb->side()->n_nodes(); i++)
547  permutation[node_numbers[nb->side()->node(i)]] = i;
548 
549  switch (nb->side()->dim())
550  {
551  case 0:
552  nb->side()->element()->permutation_idx_[nb->side()->el_idx()] = RefElement<1>::permutation_index(permutation);
553  break;
554  case 1:
555  nb->side()->element()->permutation_idx_[nb->side()->el_idx()] = RefElement<2>::permutation_index(permutation);
556  break;
557  case 2:
558  nb->side()->element()->permutation_idx_[nb->side()->el_idx()] = RefElement<3>::permutation_index(permutation);
559  break;
560  }
561  }
562 }
563 
564 
565 
566 
567 
568 /**
569  * Set Element->boundaries_ for all external sides. (temporary solution)
570  */
572 {
573  /*
574  // set to non zero all pointers including boundary connected to lower dim elements
575  // these have only one side per edge
576  Boundary empty_boundary;
577 
578 
579  FOR_ELEMENTS(this, ele) {
580  // is there any outer side
581  bool outer=false;
582  FOR_ELEMENT_SIDES(ele, si)
583  {
584  if ( ele->side(si)->edge()->n_sides == 1) {
585  outer=true;
586  break;
587  }
588  }
589  if (outer) {
590  // for elements on the boundary set boundaries_
591  FOR_ELEMENT_SIDES(ele,si)
592  if ( ele->side(si)->edge()->n_sides == 1)
593  ele->boundaries_[si] = &empty_boundary;
594  else
595  ele->boundaries_[si] = NULL;
596 
597  } else {
598  // can delete boundaries on internal elements !!
599  delete ele->boundaries_;
600  ele->boundaries_=NULL;
601  }
602  }
603 
604  int count=0;
605  // pass through neighbours and set to NULL internal interfaces
606  FOR_NEIGHBOURS(this, ngh ) {
607  SideIter s = ngh->side();
608  if (s->element()->boundaries_ == NULL) continue;
609  s->element()->boundaries_[ s->el_idx() ] = NULL;
610  }
611 
612  // count remaining
613  unsigned int n_boundaries=0;
614  FOR_ELEMENTS(this, ele) {
615  if (ele->boundaries_ == NULL) continue;
616  FOR_ELEMENT_SIDES(ele, si)
617  if (ele->boundaries_[si]) n_boundaries ++;
618  }
619 
620  // fill boundaries
621  BoundaryFullIter bcd(boundary);
622  unsigned int ni;
623 
624  boundary.reserve(n_boundaries);
625  DBGMSG("bc_elements size after read: %d\n", bc_elements.size());
626  bc_elements.reserve(n_boundaries);
627  FOR_ELEMENTS(this, ele) {
628  if (ele->boundaries_ == NULL) continue;
629  FOR_ELEMENT_SIDES(ele, si)
630  if (ele->boundaries_[si]) {
631  // add boundary object
632  bcd = boundary.add_item();
633 
634 
635  // fill boundary element
636  Element * bc_ele = bc_elements.add_item( -bcd.index() ); // use negative bcd index as ID,
637  bc_ele->dim_ = ele->dim()-1;
638  bc_ele->node = new Node * [bc_ele->n_nodes()];
639  FOR_ELEMENT_NODES(bc_ele, ni) {
640  bc_ele->node[ni] = (Node *)ele->side(si)->node(ni);
641  }
642 
643  // fill Boudary object
644  bcd->side = ele->side(si);
645  bcd->bc_element_ = bc_ele;
646  ele->boundaries_[si] = bcd;
647 
648  }
649  }*/
650 }
651 
652 
653 //=============================================================================
654 //
655 //=============================================================================
657 {
658 
659  xprintf( MsgVerb, " Element to neighbours of vb2 type... ")/*orig verb 5*/;
660 
661  FOR_ELEMENTS(this,ele) ele->n_neighs_vb =0;
662 
663  // count vb neighs per element
664  FOR_NEIGHBOURS(this, ngh ) ngh->element_->n_neighs_vb++;
665 
666  // Allocation of the array per element
667  FOR_ELEMENTS(this, ele )
668  if( ele->n_neighs_vb > 0 ) {
669  ele->neigh_vb = new struct Neighbour* [ele->n_neighs_vb];
670  ele->n_neighs_vb=0;
671  }
672 
673  // fill
674  ElementIter ele;
675  FOR_NEIGHBOURS(this, ngh ) {
676  ele = ngh->element();
677  ele->neigh_vb[ ele->n_neighs_vb++ ] = &( *ngh );
678  }
679 
680  xprintf( MsgVerb, "O.K.\n")/*orig verb 6*/;
681 }
682 
683 
684 
685 
689 
690 
692  /* Algorithm:
693  *
694  * 1) create BIH tree
695  * 2) for every 1D, find list of candidates
696  * 3) compute intersections for 1d, store it to master_elements
697  *
698  */
699  BIHTree bih_tree( this );
700  master_elements.resize(n_elements());
701 
702  for(unsigned int i_ele=0; i_ele<n_elements(); i_ele++) {
703  Element &ele = this->element[i_ele];
704 
705  if (ele.dim() == 1) {
706  vector<unsigned int> candidate_list;
707  //bih_tree.find_bounding_box(ele.bounding_box(), candidate_list);
708  //DBGMSG("1d el: %d n_cand: %d\n", ele.index(), candidate_list.size() );
709  //for(int i_elm: candidate_list) {
710  for(unsigned int i_elm=0; i_elm<n_elements(); i_elm++) {
711  ElementFullIter elm = this->element( i_elm );
712  if (elm->dim() == 2) {
713  IntersectionLocal *intersection;
714  GetIntersection( TAbscissa(ele), TTriangle(*elm), intersection);
715  //DBGMSG("2d el: %d %p\n", elm->index(), intersection);
716  if (intersection && intersection->get_type() == IntersectionLocal::line) {
717 
718  master_elements[i_ele].push_back( intersections.size() );
719  intersections.push_back( Intersection(this->element(i_ele), elm, intersection) );
720  //DBGMSG("isec: %d %d %g\n" , ele.index(), elm->index(),
721  // intersections.back().intersection_true_size() );
722  }
723  }
724 
725  }
726  }
727  }
728 
729 }
730 
731 
732 
733 ElementAccessor<3> Mesh::element_accessor(unsigned int idx, bool boundary) {
734  return ElementAccessor<3>(this, idx, boundary);
735 }
736 
737 
738 
739 vector<int> const & Mesh::elements_id_maps( bool boundary_domain) const
740 {
741  if (bulk_elements_id_.size() ==0) {
743  int last_id;
744 
745  bulk_elements_id_.resize(n_elements());
746  map_it = bulk_elements_id_.begin();
747  last_id = -1;
748  for(unsigned int idx=0; idx < element.size(); idx++, ++map_it) {
749  int id = element.get_id(idx);
750  if (last_id >= id) xprintf(UsrErr, "Element IDs in non-increasing order, ID: %d\n", id);
751  last_id=*map_it = id;
752 // DBGMSG("bulk map: %d\n", *map_it);
753  }
754 
756  map_it = boundary_elements_id_.begin();
757  last_id = -1;
758  for(unsigned int idx=0; idx < bc_elements.size(); idx++, ++map_it) {
759  int id = bc_elements.get_id(idx);
760  // We set ID for boundary elements created by the mesh itself to "-1"
761  // this force gmsh reader to skip all remaining entries in boundary_elements_id_
762  // and thus report error for any remaining data lines
763  if (id < 0) last_id=*map_it=-1;
764  else {
765  if (last_id >= id) xprintf(UsrErr, "Element IDs in non-increasing order, ID: %d\n", id);
766  last_id=*map_it = id;
767  }
768 // DBGMSG("bc map: %d\n", *map_it);
769  }
770  }
771 
772  if (boundary_domain) return boundary_elements_id_;
773  return bulk_elements_id_;
774 }
775 
776 //-----------------------------------------------------------------------------
777 // vim: set cindent:
int n_triangles
Definition: mesh.h:236
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:86
void read_regions_from_input(Input::Array region_list, MapElementIDToRegionID &map)
Definition: region.cc:474
vector< vector< unsigned int > > node_elements
Definition: mesh.h:322
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:479
Definition: system.hh:75
void read_sets_from_input(Input::Array arr)
Definition: region.cc:374
void count_side_types()
Definition: mesh.cc:244
int MPI_Comm
Definition: mpi.h:141
void make_intersec_elements()
Definition: mesh.cc:691
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:359
int n_lines
Definition: mesh.h:235
static const unsigned int undef_idx
Definition: mesh.h:110
void create_node_element_lists()
Definition: mesh.cc:260
static unsigned int permutation_index(unsigned int p[n_nodes_per_side])
Definition: ref_element.cc:158
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:325
boost::shared_ptr< Partitioning > part_
Definition: mesh.h:335
int index() const
Definition: sys_vector.hh:88
Input::Record in_record_
Definition: mesh.h:339
void create_external_boundary()
Definition: mesh.cc:571
FullIter add_item(int id)
Definition: sys_vector.hh:379
int n_sides
Definition: edges.h:48
vector< vector< vector< unsigned int > > > side_nodes
Definition: mesh.h:248
Definition: edges.h:38
int n_tetrahedras
Definition: mesh.h:237
ElementAccessor< 3 > element_accessor(unsigned int idx, bool boundary=false)
Definition: mesh.cc:733
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:163
unsigned int n_sides()
Definition: mesh.cc:152
Class for declaration of inputs sequences.
Definition: type_base.hh:230
int n_exsides
Definition: mesh.h:232
unsigned int dim() const
unsigned int size() const
Returns size of the container. This is independent of the allocated space.
Definition: sys_vector.hh:412
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:321
#define ASSERT(...)
Definition: global_defs.h:120
Definition: system.hh:75
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:235
vector< int > const & elements_id_maps(bool boundary_domain) const
Definition: mesh.cc:739
void setup_topology()
Definition: mesh.cc:224
Neighbour ** neigh_vb
Definition: elements.h:129
Accessor to the data with type Type::Record.
Definition: accessors.hh:308
#define NDEF
Definition: mesh.cc:65
const Ret val(const string &key) const
unsigned int n_sides() const
#define xprintf(...)
Definition: system.hh:104
#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:394
#define FOR_NEIGHBOURS(_mesh_, it)
Definition: mesh.h:410
ElementVector bc_elements
Definition: mesh.h:206
void close()
Definition: region.cc:224
SideIter side(const unsigned int loc_index)
FullIter begin()
Definition: sys_vector.hh:404
const Record & close() const
Definition: type_record.cc:290
Region implicit_boundary_region()
Definition: region.cc:85
unsigned int edge_idx_
Definition: boundaries.h:94
vector< int > boundary_elements_id_
Definition: mesh.h:317
unsigned int bc_ele_idx_
Definition: boundaries.h:95
void count_element_types()
Definition: mesh.cc:172
Dedicated class for storing path to input and output files.
Definition: file_path.hh:34
Definition: system.hh:75
void read_gmsh_from_stream(istream &in)
Definition: mesh.cc:188
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:233
static Input::Type::Record region_input_type
Definition: region.hh:306
static Input::Type::Record input_type
Definition: mesh.h:111
RegionDB region_db_
Definition: mesh.h:331
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:231
void intersect_element_lists(vector< unsigned int > const &nodes_list, vector< unsigned int > &intersection_element_list)
Definition: mesh.cc:273
void make_edge_permutations()
Definition: mesh.cc:500
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc...
Mesh * mesh_
Definition: boundaries.h:96
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:656
static Input::Type::Record region_set_input_type
Definition: region.hh:310
const Node * node(unsigned int i) const
Definition: side_impl.hh:35
void reinit(Input::Record in_record)
Definition: mesh.cc:107
Reader for (slightly) modified JSON files.
FullIter end()
Returns FullFullIterer of the fictions past the end element.
Definition: sys_vector.hh:408
#define FOR_SIDES(_mesh_, it)
Definition: mesh.h:403
void init_from_input()
Definition: mesh.cc:199
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:301
vector< int > bulk_elements_id_
Definition: mesh.h:317
IntersectionType get_type() const
void make_neighbours_and_edges()
Definition: mesh.cc:338
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:390