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