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