Flow123d  JS_before_hm-2127-g0296217bc
output_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 output_mesh.cc
15  * @brief Classes for auxiliary output mesh.
16  */
17 
18 #include "system/index_types.hh"
19 #include "output_mesh.hh"
20 #include "output_element.hh"
21 #include "mesh/mesh.h"
22 #include "mesh/ref_element.hh"
23 #include "mesh/accessors.hh"
24 #include "mesh/node_accessor.hh"
25 #include "mesh/range_wrapper.hh"
26 #include "la/distribution.hh"
27 
28 
29 namespace IT=Input::Type;
30 
32  return IT::Record("OutputMesh", "Parameters of the refined output mesh. [Not impemented]")
33  .declare_key("max_level", IT::Integer(1,20),IT::Default("3"),
34  "Maximal level of refinement of the output mesh.")
35  .declare_key("refine_by_error", IT::Bool(), IT::Default("false"),
36  "Set true for using ``error_control_field``. Set false for global uniform refinement to max_level.")
37  .declare_key("error_control_field",IT::String(), IT::Default::optional(),
38  "Name of an output field, according to which the output mesh will be refined. The field must be a SCALAR one.")
39  .declare_key("refinement_error_tolerance",IT::Double(0.0), IT::Default("0.01"),
40  "Tolerance for element refinement by error. If tolerance is reached, refinement is stopped."
41  "Relative difference between error control field and its linear approximation on element is computed"
42  "and compared with tolerance.")
43  .close();
44 }
45 
47 :
48  orig_mesh_(&mesh),
49  max_level_(0),
50  refine_by_error_(false),
51  refinement_error_tolerance_(0.0),
52  el_ds_(nullptr),
53  node_ds_(nullptr)
54 {
55 }
56 
57 
59 :
60  input_record_(in_rec),
61  orig_mesh_(&mesh),
62  max_level_(input_record_.val<int>("max_level")),
63  refine_by_error_(input_record_.val<bool>("refine_by_error")),
64  refinement_error_tolerance_(input_record_.val<double>("refinement_error_tolerance")),
65  el_ds_(nullptr),
66  node_ds_(nullptr)
67 {
68 }
69 
71 {
72  // Refined mesh creates own special distributions and local to global maps and needs destroy these objects.
73  if ( (el_ds_!=nullptr) && (mesh_type_ == MeshType::refined)) {
74  delete[] el_4_loc_;
75  delete[] node_4_loc_;
76  delete el_ds_;
77  delete node_ds_;
78  }
79 }
80 
82 {
84 // ASSERT_DBG(offsets_->n_values() > 0);
85  return OutputElementIterator(OutputElement(0, shared_from_this()));
86 }
87 
89 {
91 // ASSERT_DBG(offsets_->n_values() > 0);
92  return OutputElementIterator(OutputElement(offsets_->n_values()-1, shared_from_this()));
93 }
94 
95 
97 {
98  error_control_field_func_ = error_control_field_func;
99 }
100 
102 {
104  return offsets_->n_values()-1;
105 }
106 
108 {
110  return nodes_->n_values();
111 }
112 
114 {
115  unsigned int elm_idx[1];
116  unsigned int node_idx[1];
117  unsigned int region_idx[1];
118  int partition[1];
119  elem_ids_ = std::make_shared< ElementDataCache<unsigned int> >("elements_ids", (unsigned int)1, this->n_elements());
120  node_ids_ = std::make_shared< ElementDataCache<unsigned int> >("node_ids", (unsigned int)1, this->n_nodes());
121  region_ids_ = std::make_shared< ElementDataCache<unsigned int> >("region_ids", (unsigned int)1, this->n_elements());
122  partitions_ = std::make_shared< ElementDataCache<int> >("partitions", (unsigned int)1, this->n_elements());
123  OutputElementIterator it = this->begin();
124  for (unsigned int i = 0; i < this->n_elements(); ++i, ++it) {
125  if (mesh_type_ == MeshType::orig) elm_idx[0] = orig_mesh_->find_elem_id(it->idx());
126  else elm_idx[0] = it->idx();
127  elem_ids_->store_value( i, elm_idx );
128 
129  region_idx[0] = orig_mesh_->element_accessor(it->idx()).region().id();
130  region_ids_->store_value( i, region_idx );
131 
132  partition[0] = orig_mesh_->element_accessor(it->idx())->pid();
133  partitions_->store_value( i, partition );
134 
135  std::vector< unsigned int > node_list = it->node_list();
136  for (unsigned int j = 0; j < it->n_nodes(); ++j) {
137  if (mesh_type_ == MeshType::orig) node_idx[0] = orig_mesh_->find_node_id(node_list[j]);
138  else node_idx[0] = node_list[j];
139  node_ids_->store_value( node_list[j], node_idx );
140  }
141  }
142 }
143 
144 
146 {
147  return (nodes_ && connectivity_ && offsets_);
148 }
149 
150 
152 {
153  ASSERT( !is_created() ).error("Multiple initialization of OutputMesh!\n");
154 
155  DebugOut() << "Create output submesh containing only local elements.";
156 
157  unsigned int ele_id = 0,
158  offset = 0, // offset of node indices of element in node vector
159  coord_id = 0, // coordinate id in node vector
160  conn_id = 0; // index to connectivity vector
161  ElementAccessor<3> elm;
162 
168 
169  const unsigned int n_local_elements = el_ds_->lsize();
170  unsigned int n_nodes = node_ds_->end( node_ds_->np()-1 );
171  std::vector<unsigned int> local_nodes_map(n_nodes, undef_idx); // map global to local ids of nodes
172  for (unsigned int i=0; i<n_local_nodes_; ++i) local_nodes_map[ node_4_loc_[i] ] = i;
173 
174  orig_element_indices_ = std::make_shared<std::vector<unsigned int>>(n_local_elements);
175  offsets_ = std::make_shared<ElementDataCache<unsigned int>>("offsets", (unsigned int)ElementDataCacheBase::N_SCALAR, n_local_elements+1);
176  auto &offset_vec = *( offsets_->get_component_data(0).get() );
177 
178  offset_vec[0] = 0;
179  loc_4_el_ = new LongIdx [ el_ds_->size() ];
180  for (unsigned int loc_el = 0; loc_el < el_ds_->size(); loc_el++) loc_4_el_[loc_el] = -1;
181  for (unsigned int loc_el = 0; loc_el < n_local_elements; loc_el++) {
182  loc_4_el_[ el_4_loc_[loc_el] ] = loc_el;
183  elm = orig_mesh_->element_accessor( el_4_loc_[loc_el] );
184  // increase offset by number of nodes of the simplicial element
185  offset += elm->dim() + 1;
186  offset_vec[ele_id+1] = offset;
187  (*orig_element_indices_)[ele_id] = el_4_loc_[loc_el];
188  ele_id++;
189  }
190 
191  connectivity_ = std::make_shared<ElementDataCache<unsigned int>>("connectivity", (unsigned int)ElementDataCacheBase::N_SCALAR,
192  offset_vec[offset_vec.size()-1]);
193  auto &connectivity_vec = *( connectivity_->get_component_data(0).get() );
194  for (unsigned int loc_el = 0; loc_el < n_local_elements; loc_el++) {
195  elm = orig_mesh_->element_accessor( el_4_loc_[loc_el] );
196  for (unsigned int li=0; li<elm->n_nodes(); li++) {
197  ASSERT_DBG(local_nodes_map[ elm.node(li).idx() ] != undef_idx)(elm.node(li).idx()).error("Undefined global to local node index!");
198  connectivity_vec[conn_id++] = local_nodes_map[ elm.node(li).idx() ];
199  }
200  }
201 
202  // set coords of nodes
203  nodes_ = std::make_shared<ElementDataCache<double>>("", (unsigned int)ElementDataCacheBase::N_VECTOR, n_local_nodes_);
204  auto &node_vec = *( nodes_->get_component_data(0) );
205  for(unsigned int i_node=0; i_node<local_nodes_map.size(); ++i_node) {
206  if (local_nodes_map[i_node]==undef_idx) continue; // skip element if it is not local
207  auto node = *orig_mesh_->node(i_node);
208  coord_id = 3*local_nodes_map[i_node]; // id of first coordinates in node_vec
209  node_vec[coord_id++] = node[0];
210  node_vec[coord_id++] = node[1];
211  node_vec[coord_id] = node[2];
212  }
213 }
214 
215 
216 
218 {
219  std::shared_ptr<ElementDataCache<unsigned int>> global_offsets; // needs for creating serial nodes and connectivity caches on zero process
220  auto elems_n_nodes = get_elems_n_nodes(); // collects number of nodes on each elements (for fill master_mesh_->offsets_)
221  int rank = el_ds_->myp();
222 
223  if (rank==0) {
224  // create serial output mesh, fill offsets cache and orig_element_indices vector
225  unsigned int n_elems = el_ds_->end( el_ds_->np()-1 );
226  master_mesh_ = this->construct_mesh();
227  master_mesh_->orig_element_indices_ = std::make_shared<std::vector<unsigned int>>(n_elems);
228  master_mesh_->offsets_ = std::make_shared<ElementDataCache<unsigned int>>("offsets", ElementDataCacheBase::N_SCALAR, n_elems+1);
229  auto &offsets_vec = *( master_mesh_->offsets_->get_component_data(0).get() );
230  auto &elems_n_nodes_vec = *( elems_n_nodes->get_component_data(0).get() );
231  unsigned int offset=0;
232  offsets_vec[0] = 0;
233  for (unsigned int i=0; i<n_elems; ++i) {
234  offset += elems_n_nodes_vec[i];
235  offsets_vec[i+1] = offset;
236  (*master_mesh_->orig_element_indices_)[i] = i;
237  }
238  global_offsets = master_mesh_->offsets_;
239  }
240 
241  // collects serial caches
242  std::shared_ptr<ElementDataCache<double>> serial_nodes_cache = make_serial_nodes_cache(global_offsets);
243  std::shared_ptr<ElementDataCache<unsigned int>> serial_connectivity_cache = make_serial_connectivity_cache(global_offsets);
244 
245  if (rank==0) {
246  // set serial output mesh caches
247  master_mesh_->connectivity_ = serial_connectivity_cache;
248  master_mesh_->nodes_ = serial_nodes_cache;
249 
250  master_mesh_->mesh_type_ = this->mesh_type_;
251  }
252 }
253 
254 
255 std::shared_ptr<ElementDataCache<unsigned int>> OutputMeshBase::get_elems_n_nodes()
256 {
257  // Compute (locally) number of nodes of each elements
258  ElementDataCache<unsigned int> local_elems_n_nodes("elems_n_nodes", ElementDataCacheBase::N_SCALAR, offsets_->n_values()-1);
259  auto &local_elems_n_nodes_vec = *( local_elems_n_nodes.get_component_data(0).get() );
260  auto &offset_vec = *( offsets_->get_component_data(0).get() );
261  for (unsigned int i=0; i<local_elems_n_nodes.n_values(); ++i) local_elems_n_nodes_vec[i] = offset_vec[i+1] - offset_vec[i];
262 
263  // Collect data, set on zero process
264  std::shared_ptr<ElementDataCache<unsigned int>> global_elems_n_nodes;
265  auto gather_cache = local_elems_n_nodes.gather(el_ds_, el_4_loc_);
266  if (el_ds_->myp()==0) global_elems_n_nodes = std::dynamic_pointer_cast< ElementDataCache<unsigned int> >(gather_cache);
267  return global_elems_n_nodes;
268 }
269 
270 
271 
272 //////////////////////////////////////////////////////////////////////////////////////////////////////////////
273 
274 
276 : OutputMeshBase(mesh)
277 {
278  this->mesh_type_ = MeshType::orig;
279 }
280 
282 : OutputMeshBase(mesh, in_rec)
283 {
284  this->mesh_type_ = MeshType::orig;
285 }
286 
287 
289 {
290 }
291 
292 
294 {
295  ASSERT(0).error("Not implemented yet.");
296 }
297 
299 {
300  ASSERT(0).error("Not implemented yet.");
301  return false;
302 }
303 
304 
305 std::shared_ptr<OutputMeshBase> OutputMesh::construct_mesh()
306 {
307  return std::make_shared<OutputMesh>(*orig_mesh_);
308 }
309 
310 
311 std::shared_ptr<ElementDataCache<double>> OutputMesh::make_serial_nodes_cache(FMT_UNUSED std::shared_ptr<ElementDataCache<unsigned int>> global_offsets)
312 {
313  std::shared_ptr<ElementDataCache<double>> serial_nodes_cache;
314 
315  // collects nodes_ data (coordinates)
316  auto serial_nodes = nodes_->gather(node_ds_, node_4_loc_);
317 
318  if (el_ds_->myp()==0) serial_nodes_cache = std::dynamic_pointer_cast< ElementDataCache<double> >(serial_nodes);
319  return serial_nodes_cache;
320 }
321 
322 
323 std::shared_ptr<ElementDataCache<unsigned int>> OutputMesh::make_serial_connectivity_cache(std::shared_ptr<ElementDataCache<unsigned int>> global_offsets)
324 {
325  std::shared_ptr<ElementDataCache<unsigned int>> serial_connectivity_cache;
326 
327  // re-number connectivity indices from local to global
328  auto &conn_vec = *( connectivity_->get_component_data(0).get() );
329  ElementDataCache<unsigned int> global_conn("connectivity", (unsigned int)1, conn_vec.size()); // holds global indices of nodes
330  auto &global_conn_vec = *( global_conn.get_component_data(0).get() );
331  for(unsigned int i=0; i<conn_vec.size(); i++) {
332  global_conn_vec[i] = node_4_loc_[ conn_vec[i] ];
333  }
334 
335  // collects global connectivities
336  auto &local_offset_vec = *( offsets_->get_component_data(0).get() );
337  auto global_fix_size_conn = global_conn.element_node_cache_fixed_size(local_offset_vec);
338  auto collective_conn = global_fix_size_conn->gather(el_ds_, el_4_loc_);
339 
340  if (el_ds_->myp()==0) {
341  auto &offset_vec = *( global_offsets->get_component_data(0).get() );
342  serial_connectivity_cache = std::dynamic_pointer_cast< ElementDataCache<unsigned int> >( collective_conn->element_node_cache_optimize_size(offset_vec) );
343  }
344  return serial_connectivity_cache;
345 }
346 
347 
348 //////////////////////////////////////////////////////////////////////////////////////////////////////////////
349 
351 : OutputMeshBase(mesh)
352 {
353  this->mesh_type_ = MeshType::discont;
354 }
355 
357 : OutputMeshBase(mesh, in_rec)
358 {
359  this->mesh_type_ = MeshType::discont;
360 }
361 
362 
364 {
365 }
366 
367 
368 template<int dim>
371  const ElementAccessor<spacedim> &ele_acc)
372 {
373  static const unsigned int n_subelements = 1 << dim; //2^dim
374 
375 // The refinement of elements for the output mesh is done using edge splitting
376 // technique (so called red refinement). Since we use this only for better output
377 // visualization of non-polynomial solutions, we do not care for existence of hanging
378 // nodes.
379 // In 2D case, it is straightforward process: find the midpoints of all sides,
380 // connect them and generate 4 triangles. These triangles are congruent and have
381 // equal surface areas.
382 // On the other hand, the 3D case is more complicated. After splitting the
383 // edges, we obtain 4 tetrahedra at the vertices of the original one. The octahedron
384 // that remains in the middle can be subdivided according to one of its three
385 // diagonals. Only the choice of the shortest octahedron diagonal leads to a regular
386 // tetrahedra decomposition. This algorithm originally comes from Bey.
387 // Bey's algorithm (red refinement of tetrahedron):
388 // p.29 https://www5.in.tum.de/pub/Joshi2016_Thesis.pdf
389 // p.108 http://www.bcamath.org/documentos_public/archivos/publicaciones/sergey_book.pdf
390 // https://www.math.uci.edu/~chenlong/iFEM/doc/html/uniformrefine3doc.html#1
391 // J. Bey. Simplicial grid refinement: on Freudenthal's algorithm and the optimal number of congruence classes.
392 // Numer. Math. 85(1):1--29, 2000. p11 Algorithm: RedRefinement3D.
393 // p.4 http://www.vis.uni-stuttgart.de/uploads/tx_vispublications/vis97-grosso.pdf
394 
395  // connectivity of refined element
396  // these arrays are hardwired to the current reference element
397  static const std::vector<unsigned int> conn[] = {
398  {}, //0D
399 
400  //1D:
401  // 0,1 original nodes, 2 is in the middle
402  // get 2 elements
403  {0, 2,
404  2, 1},
405 
406  //2D:
407  // 0,1,2 original nodes
408  // 3,4,5 nodes are in the middle of sides 0,1,2 in the respective order
409  // get 4 elements
410  {0, 3, 4,
411  3, 1, 5,
412  4, 5, 2,
413  3, 5, 4},
414 
415  //3D:
416  // 0,1,2,3 original nodes
417  // 4,5,6,7,8,9 are nodes in the middle of edges 0,1,2,3,4,5 in the respective order
418  // first 4 tetrahedra are at the original nodes
419  // next 4 tetrahedra are from the inner octahedron - 4 choices according to the diagonal
420  {1, 7, 4, 8,
421  7, 2, 5, 9,
422  4, 5, 0, 6,
423  8, 9, 6, 3,
424  7, 4, 8, 9, // 4-9 octahedron diagonal
425  7, 4, 5, 9,
426  4, 8, 9, 6,
427  4, 5, 9, 6},
428 
429  {1, 7, 4, 8,
430  7, 2, 5, 9,
431  4, 5, 0, 6,
432  8, 9, 6, 3,
433  8, 5, 4, 6, // 5-8 octahedron diagonal
434  5, 8, 9, 6,
435  5, 8, 4, 7,
436  8, 5, 9, 7},
437 
438  {1, 7, 4, 8,
439  7, 2, 5, 9,
440  4, 5, 0, 6,
441  8, 9, 6, 3,
442  6, 8, 7, 9, // 6-7 octahedron diagonal
443  8, 6, 7, 4,
444  6, 5, 7, 4,
445  5, 6, 7, 9}
446  };
447 // DBGMSG("level = %d, %d\n", aux_element.level, max_refinement_level_);
448 
449  ASSERT_DBG(dim == aux_element.nodes.size()-1);
450 
451  // if not refining any further, push into final vector
452  if( ! refinement_criterion(aux_element, ele_acc) ) {
453  refinement.push_back(aux_element);
454  return;
455  }
456 
457  std::vector<AuxElement> subelements(n_subelements);
458 
459  const unsigned int n_old_nodes = RefElement<dim>::n_nodes,
460  n_new_nodes = RefElement<dim>::n_lines; // new points are in the center of lines
461 
462  // auxiliary vectors
463  std::vector<Space<spacedim>::Point> nodes = aux_element.nodes;
464 
465  nodes.reserve(n_old_nodes+n_new_nodes);
466 
467  // create new points in the element
468  for(unsigned int e=0; e < n_new_nodes; e++)
469  {
472  nodes.push_back( p / 2.0);
473  //nodes.back().print();
474  }
475 
476  unsigned int diagonal = 0;
477  // find shortest diagonal: [0]:4-9, [1]:5-8 or [2]:6-7
478  if(dim == 3){
479  double min_diagonal = arma::norm(nodes[4]-nodes[9],2);
480  double d = arma::norm(nodes[5]-nodes[8],2);
481  if(d < min_diagonal){
482  min_diagonal = d;
483  diagonal = 1;
484  }
485  d = arma::norm(nodes[6]-nodes[7],2);
486  if(d < min_diagonal){
487  min_diagonal = d;
488  diagonal = 2;
489  }
490  }
491 
492  for(unsigned int i=0; i < n_subelements; i++)
493  {
494  AuxElement& sub_ele = subelements[i];
495  sub_ele.nodes.resize(n_old_nodes);
496  sub_ele.level = aux_element.level+1;
497 
498  // over nodes
499  for(unsigned int j=0; j < n_old_nodes; j++)
500  {
501  unsigned int conn_id = (n_old_nodes)*i + j;
502  sub_ele.nodes[j] = nodes[conn[dim+diagonal][conn_id]];
503  }
504  refine_aux_element<dim>(sub_ele, refinement, ele_acc);
505  }
506 }
507 
508 
509 
510 template void OutputMeshDiscontinuous::refine_aux_element<1>(const OutputMeshDiscontinuous::AuxElement&,std::vector< OutputMeshDiscontinuous::AuxElement >&, const ElementAccessor<spacedim> &);
511 template void OutputMeshDiscontinuous::refine_aux_element<2>(const OutputMeshDiscontinuous::AuxElement&,std::vector< OutputMeshDiscontinuous::AuxElement >&, const ElementAccessor<spacedim> &);
512 template void OutputMeshDiscontinuous::refine_aux_element<3>(const OutputMeshDiscontinuous::AuxElement&,std::vector< OutputMeshDiscontinuous::AuxElement >&, const ElementAccessor<spacedim> &);
513 
514 
516  const ElementAccessor<spacedim> &ele_acc)
517 {
518  // check refinement criteria:
519 
520  //first check max. level
521  bool refine = refinement_criterion_uniform(aux_ele);
522 
523  //if max. level not reached and refinement by error is set
524  if(refine && refine_by_error_)
525  {
526  // compute centre of aux element
527  Space<spacedim>::Point centre({0,0,0});
528  for(auto& v : aux_ele.nodes ) centre += v;
529  centre = centre/aux_ele.nodes.size();
530  return refinement_criterion_error(aux_ele, centre, ele_acc);
531  }
532 
533  return refine;
534 }
535 
537 {
538  return (ele.level < max_level_);
539 }
540 
542  const Space<spacedim>::Point &centre,
543  const ElementAccessor<spacedim> &ele_acc
544  )
545 {
546  ASSERT_DBG(error_control_field_func_).error("Error control field not set!");
547 
548  // evaluate at nodes and center in a single call
549  std::vector<double> val_list(ele.nodes.size()+1);
550  Armor::array point_list(spacedim,1,1+ele.nodes.size());
551  point_list.set(0) = centre;
552  unsigned int i=0;
553  for (auto node : ele.nodes) point_list.set(++i) = node;
554  error_control_field_func_(point_list, ele_acc, val_list);
555 
556  //TODO: compute L1 or L2 error using standard quadrature
557 
558  //compare average value at nodes with value at center
559 
560  double average_val = 0.0;
561  for(unsigned int i=1; i<ele.nodes.size()+1; ++i)//(double& v: nodes_val)
562  average_val += val_list[i];
563  average_val = average_val / ele.nodes.size();
564 
565  double diff = std::abs((average_val - val_list[0])/val_list[0]);
566 // DebugOut().fmt("diff: {} {} {}\n", diff, average_val, val_list[0]);
567  return ( diff > refinement_error_tolerance_);
568 
569 }
570 
571 
572 std::shared_ptr<OutputMeshBase> OutputMeshDiscontinuous::construct_mesh()
573 {
574  return std::make_shared<OutputMeshDiscontinuous>(*orig_mesh_);
575 }
576 
577 
578 std::shared_ptr<ElementDataCache<double>> OutputMeshDiscontinuous::make_serial_nodes_cache(std::shared_ptr<ElementDataCache<unsigned int>> global_offsets)
579 {
580  std::shared_ptr<ElementDataCache<double>> serial_nodes_cache;
581 
582  // Create helper cache of discontinuous node data ordering by elements
583  std::shared_ptr< ElementDataCache<double> > discont_node_cache = std::make_shared<ElementDataCache<double>>("",
584  ElementDataCacheBase::N_VECTOR, this->connectivity_->n_values());
585  auto &discont_node_vec = *( discont_node_cache->get_component_data(0).get() );
586  auto &local_nodes_vec = *( this->nodes_->get_component_data(0).get() );
587  auto &local_conn_vec = *( this->connectivity_->get_component_data(0).get() );
588  auto &local_offset_vec = *( this->offsets_->get_component_data(0).get() );
589  unsigned int i_old, i_new;
590  for (unsigned int i_conn=0; i_conn<this->connectivity_->n_values(); ++i_conn) {
591  i_old = local_conn_vec[i_conn] * ElementDataCacheBase::N_VECTOR;
592  i_new = i_conn * ElementDataCacheBase::N_VECTOR;
593  for(unsigned int i = 0; i < ElementDataCacheBase::N_VECTOR; i++) {
594  discont_node_vec[i_new+i] = local_nodes_vec[i_old+i];
595  }
596  }
597  // Collects node data
598  auto fix_size_node_cache = discont_node_cache->element_node_cache_fixed_size(local_offset_vec);
599  auto collect_fix_size_node_cache = fix_size_node_cache->gather(el_ds_, el_4_loc_);
600 
601  if (el_ds_->myp()==0) {
602  auto &offset_vec = *( global_offsets->get_component_data(0).get() );
603  serial_nodes_cache = std::dynamic_pointer_cast< ElementDataCache<double> >(collect_fix_size_node_cache->element_node_cache_optimize_size(offset_vec));
604  }
605  return serial_nodes_cache;
606 }
607 
608 
609 std::shared_ptr<ElementDataCache<unsigned int>> OutputMeshDiscontinuous::make_serial_connectivity_cache(std::shared_ptr<ElementDataCache<unsigned int>> global_offsets)
610 {
611  std::shared_ptr<ElementDataCache<unsigned int>> serial_connectivity_cache;
612 
613  if (el_ds_->myp()==0) {
614  auto &offset_vec = *( global_offsets->get_component_data(0).get() );
615  serial_connectivity_cache = std::make_shared<ElementDataCache<unsigned int>>("connectivity", (unsigned int)ElementDataCacheBase::N_SCALAR,
616  offset_vec[offset_vec.size()-1]);
617  auto &conn_vec = *( serial_connectivity_cache->get_component_data(0).get() );
618  for (unsigned int i=0; i<conn_vec.size(); ++i) conn_vec[i] = i;
619  }
620  return serial_connectivity_cache;
621 }
622 
623 
625 {
626  ASSERT( !is_created() ).error("Multiple initialization of OutputMesh!\n");
627 
628  DebugOut() << "Create refined discontinuous submesh containing only local elements.";
629  // initial guess of size: n_elements
630  nodes_ = std::make_shared<ElementDataCache<double>>("",(unsigned int)ElementDataCacheBase::N_VECTOR,0);
631  connectivity_ = std::make_shared<ElementDataCache<unsigned int>>("connectivity",(unsigned int)ElementDataCacheBase::N_SCALAR,0);
632  offsets_ = std::make_shared<ElementDataCache<unsigned int>>("offsets",(unsigned int)ElementDataCacheBase::N_SCALAR,0);
633  orig_element_indices_ = std::make_shared<std::vector<unsigned int>>();
634 
635  // index of last node added; set at the end of original ones
636  unsigned int last_offset = 0;
637 
638  auto &node_vec = *( nodes_->get_component_data(0).get() );
639  auto &conn_vec = *( connectivity_->get_component_data(0).get() );
640  auto &offset_vec = *( offsets_->get_component_data(0).get() );
641 
642  node_vec.reserve(4*orig_mesh_->n_nodes());
643  conn_vec.reserve(4*4*orig_mesh_->n_elements());
644  offset_vec.reserve(4*orig_mesh_->n_elements()+4);
645  offset_vec.push_back(0);
646 
647  LongIdx *el_4_loc = orig_mesh_->get_el_4_loc();
648  const unsigned int n_local_elements = orig_mesh_->get_el_ds()->lsize();
649 
650  for (unsigned int loc_el = 0; loc_el < n_local_elements; loc_el++) {
651  auto ele = orig_mesh_->element_accessor( el_4_loc[loc_el] );
652  const unsigned int
653  dim = ele->dim(),
654  ele_idx = ele.idx();
655 
656  AuxElement aux_ele;
657  aux_ele.nodes.resize(ele->n_nodes());
658  aux_ele.level = 0;
659 
660  unsigned int li;
661  for (li=0; li<ele->n_nodes(); li++) {
662  aux_ele.nodes[li] = *ele.node(li);
663  }
664 
665  std::vector<AuxElement> refinement;
666 
667  switch(dim){
668  case 1: this->refine_aux_element<1>(aux_ele, refinement, ele); break;
669  case 2: this->refine_aux_element<2>(aux_ele, refinement, ele); break;
670  case 3: this->refine_aux_element<3>(aux_ele, refinement, ele); break;
671  default: ASSERT(0 < dim && dim < 4);
672  }
673 
674  //skip unrefined element
675 // if(refinement.size() < 2) continue;
676  unsigned int node_offset = node_vec.size(),
677  con_offset = conn_vec.size();
678  node_vec.resize(node_vec.size() + (refinement.size() * (dim+1))*spacedim);
679  conn_vec.resize(conn_vec.size() + refinement.size()*(dim+1));
680 // orig_element_indices_->resize(orig_element_indices_->size() + refinement.size()*(dim+1));
681 
682 // DebugOut() << "ref size = " << refinement.size() << "\n";
683  //gather coords and connectivity (in a continous way inside element)
684  for(unsigned int i=0; i < refinement.size(); i++)
685  {
686  last_offset += dim+1;
687  offset_vec.push_back(last_offset);
688  (*orig_element_indices_).push_back(ele_idx);
689  for(unsigned int j=0; j < dim+1; j++)
690  {
691  unsigned int con = i*(dim+1) + j;
692  conn_vec[con_offset + con] = con_offset + con;
693 
694  for(unsigned int k=0; k < spacedim; k++) {
695  node_vec[node_offset + con*spacedim + k] = refinement[i].nodes[j][k];
696  }
697  }
698  }
699  }
700 
701  conn_vec.shrink_to_fit();
702  node_vec.shrink_to_fit();
703  offset_vec.shrink_to_fit();
704 
705  connectivity_->set_n_values(conn_vec.size());
706  nodes_->set_n_values(node_vec.size() / spacedim);
707  offsets_->set_n_values(offset_vec.size());
708 
709  // Create special distributions and arrays of local to global indexes of refined mesh
710  el_ds_ = new Distribution(offset_vec.size()-1, PETSC_COMM_WORLD);
711  node_ds_ = new Distribution(offset_vec[offset_vec.size()-1], PETSC_COMM_WORLD);
713  el_4_loc_ = new LongIdx [ el_ds_->lsize() ];
714  LongIdx global_el_idx = el_ds_->begin();
715  for (unsigned int i=0; i<el_ds_->lsize(); ++i, ++global_el_idx) {
716  el_4_loc_[i] = global_el_idx;
717  }
718  node_4_loc_ = new LongIdx [ node_ds_->lsize() ];
719  LongIdx global_node_idx = node_ds_->begin();
720  for (unsigned int i=0; i<node_ds_->lsize(); ++i, ++global_node_idx) {
721  node_4_loc_[i] = global_node_idx;
722  }
723 
724  mesh_type_ = MeshType::refined;
725 }
726 
727 
729 {
730  master_mesh_ = this->construct_mesh();
731  master_mesh_->offsets_ = this->offsets_;
732  master_mesh_->orig_element_indices_ = this->orig_element_indices_;
733 
734  auto &conn_vec = *( this->connectivity_->get_component_data(0).get() );
735  master_mesh_->connectivity_ = std::make_shared<ElementDataCache<unsigned int>>("connectivity", ElementDataCacheBase::N_SCALAR,
736  conn_vec.size());
737  auto &master_conn_vec = *( master_mesh_->connectivity_->get_component_data(0).get() );
738  for (unsigned int i=0; i<master_conn_vec.size(); ++i) master_conn_vec[i] = i;
739 
740  master_mesh_->nodes_ = std::make_shared<ElementDataCache<double>>("", ElementDataCacheBase::N_VECTOR, conn_vec.size());
741  auto &node_vec = *( this->nodes_->get_component_data(0).get() );
742  auto &master_node_vec = *( master_mesh_->nodes_->get_component_data(0).get() );
743  unsigned int i_own, i_master, j;
744  for (unsigned int i=0; i<conn_vec.size(); ++i) {
745  i_own = conn_vec[i]*ElementDataCacheBase::N_VECTOR;
746  i_master = i*ElementDataCacheBase::N_VECTOR;
747  for (j=0; j<ElementDataCacheBase::N_VECTOR; ++j) {
748  master_node_vec[i_master+j] = node_vec[i_own+j];
749  }
750  }
751 
752  master_mesh_->mesh_type_ = this->mesh_type_;
753 }
Space::Point
Armor::ArmaVec< double, spacedim > Point
Definition: point.hh:42
OutputMeshDiscontinuous::refinement_criterion_error
bool refinement_criterion_error(const AuxElement &ele, const Space< spacedim >::Point &centre, const ElementAccessor< spacedim > &ele_acc)
Refinement flag - measures discretisation error according to error control field.
Definition: output_mesh.cc:541
Input::Type::Bool
Class for declaration of the input of type Bool.
Definition: type_base.hh:452
Distribution::np
unsigned int np() const
get num of processors
Definition: distribution.hh:105
OutputMeshBase::spacedim
static const unsigned int spacedim
Shortcut instead of spacedim template. We suppose only spacedim=3 at the moment.
Definition: output_mesh.hh:74
OutputMeshBase::node_4_loc_
LongIdx * node_4_loc_
Index set assigning to local node index its global index.
Definition: output_mesh.hh:237
ref_element.hh
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc.
OutputElementIterator
Iter< OutputElement > OutputElementIterator
Definition: output_mesh.hh:37
ElementDataCacheBase::N_SCALAR
@ N_SCALAR
Definition: element_data_cache_base.hh:40
Distribution::lsize
unsigned int lsize(int proc) const
get local size
Definition: distribution.hh:115
OutputMeshDiscontinuous::AuxElement
Auxiliary structure defining element of refined output mesh.
Definition: output_mesh.hh:295
RefElement
Definition: ref_element.hh:339
OutputMesh::refinement_criterion
bool refinement_criterion()
Definition: output_mesh.cc:298
OutputMeshDiscontinuous::~OutputMeshDiscontinuous
~OutputMeshDiscontinuous()
Definition: output_mesh.cc:363
Iter
General iterator template. Provides iterator over objects of type Object in some container.
Definition: general_iterator.hh:32
Distribution::myp
unsigned int myp() const
get my processor
Definition: distribution.hh:107
ElementDataCache
Definition: element_data_cache.hh:44
ElementDataCache::get_component_data
ComponentDataPtr get_component_data(unsigned int component_idx)
Return vector of element data for get component.
Definition: element_data_cache.cc:71
ASSERT
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
Definition: asserts.hh:347
MeshBase::get_el_ds
Distribution * get_el_ds() const
Definition: mesh.h:119
RefElement::interact
static const IdxVector<(InDim >OutDim ? InDim+1 :dim-InDim) > interact(TInteraction< OutDim, InDim > interaction, bool inv=false)
Element::dim
unsigned int dim() const
Definition: elements.h:120
OutputMeshBase::OutputElement
friend class OutputElement
Friend provides access to vectors for element accessor class.
Definition: output_mesh.hh:242
distribution.hh
Support classes for parallel programing.
OutputMeshBase::el_ds_
Distribution * el_ds_
Parallel distribution of elements.
Definition: output_mesh.hh:236
Input::Type::Integer
Class for declaration of the integral input data.
Definition: type_base.hh:483
OutputMeshDiscontinuous::refine_aux_element
void refine_aux_element(const AuxElement &aux_element, std::vector< AuxElement > &refinement, const ElementAccessor< spacedim > &ele_acc)
Performs the actual refinement of AuxElement. Recurrent.
Definition: output_mesh.cc:369
ElementDataCache::gather
std::shared_ptr< ElementDataCacheBase > gather(Distribution *distr, LongIdx *local_to_global) override
Implements ElementDataCacheBase::gather.
Definition: element_data_cache.cc:309
ASSERT_DBG
#define ASSERT_DBG(expr)
Definition: include_fadbad.hh:28
Input::Type::Double
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:534
OutputMeshBase::OutputMeshDiscontinuous
friend class OutputMeshDiscontinuous
Definition: output_mesh.hh:247
OutputMeshBase::el_4_loc_
LongIdx * el_4_loc_
Index set assigning to local element index its global index.
Definition: output_mesh.hh:234
OutputMeshBase::error_control_field_func_
ErrorControlFieldFunc error_control_field_func_
Refinement error control field function (hold value_list function of field).
Definition: output_mesh.hh:196
std::vector< unsigned int >
OutputMeshBase::offsets_
std::shared_ptr< ElementDataCache< unsigned int > > offsets_
Vector of offsets of node indices of elements. Maps elements to their nodes in connectivity_.
Definition: output_mesh.hh:210
ElementAccessor< 3 >
Distribution::size
unsigned int size() const
get global size
Definition: distribution.hh:118
OutputMeshBase::n_nodes
unsigned int n_nodes()
Returns number of nodes.
Definition: output_mesh.cc:107
OutputMeshBase::elem_ids_
std::shared_ptr< ElementDataCache< unsigned int > > elem_ids_
Vector gets ids of elements. Data is used in GMSH output.
Definition: output_mesh.hh:215
OutputMeshBase::partitions_
std::shared_ptr< ElementDataCache< int > > partitions_
Vector gets partitions of elements. Data is used in GMSH output.
Definition: output_mesh.hh:219
OutputMeshBase::orig_mesh_
Mesh * orig_mesh_
Pointer to the computational mesh.
Definition: output_mesh.hh:190
OutputMeshDiscontinuous::refinement_criterion
bool refinement_criterion(const AuxElement &ele, const ElementAccessor< spacedim > &ele_acc)
Collects different refinement criteria results.
Definition: output_mesh.cc:515
OutputMesh::construct_mesh
std::shared_ptr< OutputMeshBase > construct_mesh() override
Implements OutputMeshBase::construct_mesh.
Definition: output_mesh.cc:305
OutputMeshDiscontinuous::make_serial_connectivity_cache
std::shared_ptr< ElementDataCache< unsigned int > > make_serial_connectivity_cache(std::shared_ptr< ElementDataCache< unsigned int >> global_offsets) override
Implements OutputMeshBase::make_serial_connectivity_cache.
Definition: output_mesh.cc:609
index_types.hh
OutputMeshBase::OutputMeshBase
OutputMeshBase(Mesh &mesh)
Constructor. Takes computational mesh as a parameter.
Definition: output_mesh.cc:46
Mesh::get_node_ds
Distribution * get_node_ds() const
Definition: mesh.h:434
Distribution::end
unsigned int end(int proc) const
get last local index +1
Definition: distribution.hh:112
ASSERT_PTR_DBG
#define ASSERT_PTR_DBG(ptr)
Definition of assert macro checking non-null pointer (PTR) only for debug mode.
Definition: asserts.hh:340
OutputMeshBase::mesh_type_
MeshType mesh_type_
Type of OutputMesh.
Definition: output_mesh.hh:198
OutputMeshBase::ErrorControlFieldFunc
std::function< void(const Armor::array &, const ElementAccessor< spacedim > &, std::vector< double > &)> ErrorControlFieldFunc
Definition: output_mesh.hh:77
OutputMeshBase::node_ids_
std::shared_ptr< ElementDataCache< unsigned int > > node_ids_
Vector gets ids of nodes. Data is used in GMSH output.
Definition: output_mesh.hh:213
OutputMeshBase::make_serial_connectivity_cache
virtual std::shared_ptr< ElementDataCache< unsigned int > > make_serial_connectivity_cache(std::shared_ptr< ElementDataCache< unsigned int >> global_offsets)=0
OutputMesh::create_refined_sub_mesh
void create_refined_sub_mesh() override
Implements OutputMeshBase::create_refined_sub_mesh.
Definition: output_mesh.cc:293
Input::Type::Default
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
OutputMeshBase::make_serial_nodes_cache
virtual std::shared_ptr< ElementDataCache< double > > make_serial_nodes_cache(std::shared_ptr< ElementDataCache< unsigned int >> global_offsets)=0
OutputMeshBase::~OutputMeshBase
virtual ~OutputMeshBase()
Definition: output_mesh.cc:70
Distribution
Definition: distribution.hh:50
accessors.hh
Input::Record
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
Mesh::get_node_4_loc
LongIdx * get_node_4_loc() const
Definition: mesh.h:437
undef_idx
const unsigned int undef_idx
Definition: index_types.hh:32
OutputMeshBase::n_elements
unsigned int n_elements()
Returns number of element.
Definition: output_mesh.cc:101
ElementDataCache::element_node_cache_fixed_size
std::shared_ptr< ElementDataCacheBase > element_node_cache_fixed_size(std::vector< unsigned int > &offset_vec) override
Implements ElementDataCacheBase::element_node_cache_fixed_size.
Definition: element_data_cache.cc:362
OutputMesh::make_serial_nodes_cache
std::shared_ptr< ElementDataCache< double > > make_serial_nodes_cache(std::shared_ptr< ElementDataCache< unsigned int >> global_offsets) override
Implements OutputMeshBase::make_serial_nodes_cache.
Definition: output_mesh.cc:311
MeshBase::element_accessor
ElementAccessor< 3 > element_accessor(unsigned int idx) const
Create and return ElementAccessor to element of given idx.
Definition: mesh.cc:866
MeshBase::node
NodeAccessor< 3 > node(unsigned int idx) const
Create and return NodeAccessor to node of given idx.
Definition: mesh.cc:872
OutputMeshDiscontinuous::construct_mesh
std::shared_ptr< OutputMeshBase > construct_mesh() override
Implements OutputMeshBase::construct_mesh.
Definition: output_mesh.cc:572
Input::Type::Record::declare_key
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:503
OutputMeshBase::n_local_nodes_
unsigned int n_local_nodes_
Hold number of local nodes (own + ghost), value is equal with size of node_4_loc array.
Definition: output_mesh.hh:239
OutputMeshBase::construct_mesh
virtual std::shared_ptr< OutputMeshBase > construct_mesh()=0
OutputMeshDiscontinuous::make_serial_nodes_cache
std::shared_ptr< ElementDataCache< double > > make_serial_nodes_cache(std::shared_ptr< ElementDataCache< unsigned int >> global_offsets) override
Implements OutputMeshBase::make_serial_nodes_cache.
Definition: output_mesh.cc:578
mesh.h
OutputMeshBase
Base class for Output mesh.
Definition: output_mesh.hh:70
OutputMeshBase::set_error_control_field
void set_error_control_field(ErrorControlFieldFunc error_control_field_func)
Selects the error control field computing function of output field set according to input record.
Definition: output_mesh.cc:96
ElementAccessor::node
NodeAccessor< 3 > node(unsigned int ni) const
Definition: accessors.hh:233
Armor::Array::set
ArrayMatSet set(uint index)
Definition: armor.hh:838
ElementDataCacheBase::N_VECTOR
@ N_VECTOR
Definition: element_data_cache_base.hh:41
Input::Type::Record::close
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
Input::Type
Definition: balance.hh:41
Input::Type::Record
Record type proxy class.
Definition: type_record.hh:182
Interaction
Definition: ref_element.hh:286
OutputMeshDiscontinuous::create_refined_sub_mesh
void create_refined_sub_mesh() override
Implements OutputMeshBase::create_refined_sub_mesh.
Definition: output_mesh.cc:624
MeshBase::find_elem_id
int find_elem_id(unsigned int pos) const
Return element id (in GMSH file) of element of given position in element vector.
Definition: mesh.h:138
LongIdx
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: index_types.hh:24
OutputMeshBase::end
OutputElementIterator end()
Gives iterator to the LAST element of the output mesh.
Definition: output_mesh.cc:88
OutputMeshBase::master_mesh_
std::shared_ptr< OutputMeshBase > master_mesh_
Definition: output_mesh.hh:227
Mesh
Definition: mesh.h:355
OutputMeshBase::max_level_
const unsigned int max_level_
Maximal level of refinement.
Definition: output_mesh.hh:193
OutputMeshBase::node_ds_
Distribution * node_ds_
Parallel distribution of nodes. Depends on elements distribution.
Definition: output_mesh.hh:238
OutputMeshBase::create_id_caches
void create_id_caches()
Create nodes and elements data caches.
Definition: output_mesh.cc:113
Input::Type::String
Class for declaration of the input data that are in string format.
Definition: type_base.hh:582
OutputMeshBase::orig_element_indices_
std::shared_ptr< std::vector< unsigned int > > orig_element_indices_
Vector of element indices in the computational mesh. (Important when refining.)
Definition: output_mesh.hh:203
OutputMeshDiscontinuous::refinement_criterion_uniform
bool refinement_criterion_uniform(const AuxElement &ele)
Refinement flag - checks only maximal level of refinement.
Definition: output_mesh.cc:536
MeshBase::get_el_4_loc
LongIdx * get_el_4_loc() const
Definition: mesh.h:122
OutputMeshBase::create_sub_mesh
void create_sub_mesh()
Definition: output_mesh.cc:151
OutputMesh::~OutputMesh
~OutputMesh()
Definition: output_mesh.cc:288
output_mesh.hh
Classes for auxiliary output mesh.
Armor::Array< double >
OutputMeshDiscontinuous::make_parallel_master_mesh
void make_parallel_master_mesh() override
Overrides OutputMeshBase::make_parallel_master_mesh.
Definition: output_mesh.cc:728
OutputMeshBase::refine_by_error_
bool refine_by_error_
True, if output mesh is to be refined by error criterion.
Definition: output_mesh.hh:199
OutputMesh::make_serial_connectivity_cache
std::shared_ptr< ElementDataCache< unsigned int > > make_serial_connectivity_cache(std::shared_ptr< ElementDataCache< unsigned int >> global_offsets) override
Implements OutputMeshBase::make_serial_connectivity_cache.
Definition: output_mesh.cc:323
OutputMeshBase::get_elems_n_nodes
std::shared_ptr< ElementDataCache< unsigned int > > get_elems_n_nodes()
Compute and return number of nodes for each elements (compute from offsets)
Definition: output_mesh.cc:255
OutputMeshDiscontinuous::AuxElement::nodes
std::vector< Space< spacedim >::Point > nodes
Definition: output_mesh.hh:296
OutputMeshBase::is_created
bool is_created()
Check if nodes_, connectivity_ and offsets_ data caches are created.
Definition: output_mesh.cc:145
MeshBase::n_elements
unsigned int n_elements() const
Definition: mesh.h:111
OutputMeshBase::loc_4_el_
LongIdx * loc_4_el_
Index set assigning to global index its local element index.
Definition: output_mesh.hh:235
Mesh::n_local_nodes
unsigned int n_local_nodes() const
Definition: mesh.h:440
OutputMeshBase::OutputMesh
friend class OutputMesh
Definition: output_mesh.hh:246
DebugOut
#define DebugOut()
Macro defining 'debug' record of log.
Definition: logger.hh:284
ASSERT_PTR
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
Definition: asserts.hh:336
ElementDataCacheBase::n_values
unsigned int n_values() const
Definition: element_data_cache_base.hh:131
MeshBase::find_node_id
int find_node_id(unsigned int pos) const
Return node id (in GMSH file) of node of given position in node vector.
Definition: mesh.h:156
Input::Type::Default::optional
static Default optional()
The factory function to make an empty default value which is optional.
Definition: type_record.hh:124
output_element.hh
Class OutputElement and its iterator OutputElementIterator on the output mesh.
OutputMeshBase::refinement_error_tolerance_
double refinement_error_tolerance_
Tolerance for error criterion refinement.
Definition: output_mesh.hh:200
OutputMeshBase::get_input_type
static const Input::Type::Record & get_input_type()
The specification of output mesh.
Definition: output_mesh.cc:31
MeshBase::n_nodes
unsigned int n_nodes() const
Definition: mesh.h:167
OutputMeshBase::make_serial_master_mesh
void make_serial_master_mesh()
Synchronize parallel data and create serial COLECTIVE output mesh on zero process.
Definition: output_mesh.cc:217
OutputMeshBase::region_ids_
std::shared_ptr< ElementDataCache< unsigned int > > region_ids_
Vector gets ids of regions. Data is used in GMSH output.
Definition: output_mesh.hh:217
Element::n_nodes
unsigned int n_nodes() const
Definition: elements.h:125
OutputMeshBase::begin
OutputElementIterator begin()
Gives iterator to the FIRST element of the output mesh.
Definition: output_mesh.cc:81
OutputMeshBase::nodes_
std::shared_ptr< ElementDataCache< double > > nodes_
Vector of node coordinates. [spacedim x n_nodes].
Definition: output_mesh.hh:206
range_wrapper.hh
Implementation of range helper class.
Distribution::begin
unsigned int begin(int proc) const
get starting local index
Definition: distribution.hh:109
FMT_UNUSED
#define FMT_UNUSED
Definition: posix.h:75
node_accessor.hh
OutputMeshBase::connectivity_
std::shared_ptr< ElementDataCache< unsigned int > > connectivity_
Vector maps the nodes to their coordinates in vector nodes_.
Definition: output_mesh.hh:208
OutputMeshDiscontinuous::AuxElement::level
unsigned int level
Definition: output_mesh.hh:297