Flow123d  JS_before_hm-887-g601087d
field_fe.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 field_fe.cc
15  * @brief
16  */
17 
18 
19 #include <limits>
20 
21 #include "fields/field_fe.hh"
22 #include "la/vector_mpi.hh"
23 #include "fields/field_instances.hh" // for instantiation macros
25 #include "input/input_type.hh"
26 #include "fem/fe_p.hh"
27 #include "fem/fe_system.hh"
28 #include "fem/dh_cell_accessor.hh"
29 #include "fem/mapping_p1.hh"
30 #include "io/reader_cache.hh"
31 #include "io/msh_gmshreader.h"
32 #include "mesh/accessors.hh"
33 #include "mesh/range_wrapper.hh"
34 #include "mesh/bc_mesh.hh"
36 
37 #include "system/sys_profiler.hh"
41 
42 
43 
44 
45 /// Implementation.
46 
47 namespace it = Input::Type;
48 
49 
50 
51 
53 
54 
55 
56 /************************************************************************************
57  * Implementation of FieldFE methods
58  */
59 template <int spacedim, class Value>
60 const Input::Type::Record & FieldFE<spacedim, Value>::get_input_type()
61 {
62  return it::Record("FieldFE", FieldAlgorithmBase<spacedim,Value>::template_name()+" Field given by finite element approximation.")
66  "GMSH mesh with data. Can be different from actual computational mesh.")
68  "Section where to find the field.\n Some sections are specific to file format: "
69  "point_data/node_data, cell_data/element_data, -/element_node_data, native/-.\n"
70  "If not given by a user, we try to find the field in all sections, but we report an error "
71  "if it is found in more than one section.")
73  "The values of the Field are read from the ```$ElementData``` section with field name given by this key.")
74  .declare_key("default_value", IT::Double(), IT::Default::optional(),
75  "Default value is set on elements which values have not been listed in the mesh data file.")
76  .declare_key("time_unit", IT::String(), IT::Default::read_time("Common unit of TimeGovernor."),
77  "Definition of the unit of all times defined in the mesh data file.")
78  .declare_key("read_time_shift", TimeGovernor::get_input_time_type(), IT::Default("0.0"),
79  "This key allows reading field data from the mesh data file shifted in time. Considering the time 't', field descriptor with time 'T', "
80  "time shift 'S', then if 't > T', we read the time frame 't + S'.")
82  IT::Default("\"equivalent_mesh\""), "Type of interpolation applied to the input spatial data.\n"
83  "The default value 'equivalent_mesh' assumes the data being constant on elements living on the same mesh "
84  "as the computational mesh, but possibly with different numbering. In the case of the same numbering, "
85  "the user can set 'identical_mesh' to omit algorithm for guessing node and element renumbering. "
86  "Alternatively, in case of different input mesh, several interpolation algorithms are available.")
87  .close();
88 }
89 
90 template <int spacedim, class Value>
92 {
93  return it::Selection("FE_discretization",
94  "Specify the section in mesh input file where field data is listed.\nSome sections are specific to file format.")
95  .add_value(OutputTime::DiscreteSpace::NODE_DATA, "node_data", "point_data (VTK) / node_data (GMSH)")
96  .add_value(OutputTime::DiscreteSpace::ELEM_DATA, "element_data", "cell_data (VTK) / element_data (GMSH)")
97  .add_value(OutputTime::DiscreteSpace::CORNER_DATA, "element_node_data", "element_node_data (only for GMSH)")
98  .add_value(OutputTime::DiscreteSpace::NATIVE_DATA, "native_data", "native_data (only for VTK)")
99  .close();
100 }
101 
102 template <int spacedim, class Value>
104 {
105  return it::Selection("interpolation", "Specify interpolation of the input data from its input mesh to the computational mesh.")
106  .add_value(DataInterpolation::identic_msh, "identic_mesh", "Topology and indices of nodes and elements of"
107  "the input mesh and the computational mesh are identical. "
108  "This interpolation is typically used for GMSH input files containing only the field values without "
109  "explicit mesh specification.")
110  .add_value(DataInterpolation::equivalent_msh, "equivalent_mesh", "Topologies of the input mesh and the computational mesh "
111  "are the same, the node and element numbering may differ. "
112  "This interpolation can be used also for VTK input data.") // default value
113  .add_value(DataInterpolation::gauss_p0, "P0_gauss", "Topologies of the input mesh and the computational mesh may differ. "
114  "Constant values on the elements of the computational mesh are evaluated using the Gaussian quadrature of the fixed order 4, "
115  "where the quadrature points and their values are found in the input mesh and input data using the BIH tree search."
116  )
117  .add_value(DataInterpolation::interp_p0, "P0_intersection", "Topologies of the input mesh and the computational mesh may differ. "
118  "Can be applied only for boundary fields. For every (boundary) element of the computational mesh the "
119  "intersection with the input mesh is computed. Constant values on the elements of the computational mesh "
120  "are evaluated as the weighted average of the (constant) values on the intersecting elements of the input mesh.")
121  .close();
122 }
123 
124 template <int spacedim, class Value>
126  Input::register_class< FieldFE<spacedim, Value>, unsigned int >("FieldFE") +
128 
129 
130 
131 template <int spacedim, class Value>
133 : FieldAlgorithmBase<spacedim, Value>(n_comp),
134  field_name_("")
135 {
136  this->is_constant_in_space_ = false;
137 }
138 
139 
140 template <int spacedim, class Value>
141 VectorMPI FieldFE<spacedim, Value>::set_fe_data(std::shared_ptr<DOFHandlerMultiDim> dh,
142  unsigned int component_index, VectorMPI dof_values)
143 {
144  dh_ = dh;
145  if (dof_values.size()==0) { //create data vector according to dof handler - Warning not tested yet
146  data_vec_ = dh_->create_vector();
148  } else {
149  data_vec_ = dof_values;
150  }
151 
152  unsigned int ndofs = dh_->max_elem_dofs();
153 
154  // initialization data of value handlers
155  FEValueInitData init_data;
156  init_data.dh = dh_;
157  init_data.data_vec = data_vec_;
158  init_data.ndofs = ndofs;
159  init_data.n_comp = this->n_comp();
160  init_data.comp_index = component_index;
161 
162  // initialize value handler objects
163  value_handler0_.initialize(init_data);
164  value_handler1_.initialize(init_data);
165  value_handler2_.initialize(init_data);
166  value_handler3_.initialize(init_data);
167 
168  // set discretization
169  discretization_ = OutputTime::DiscreteSpace::UNDEFINED;
170  interpolation_ = DataInterpolation::equivalent_msh;
171 
172  return data_vec_;
173 }
174 
175 
176 /**
177  * Returns one value in one given point. ResultType can be used to avoid some costly calculation if the result is trivial.
178  */
179 template <int spacedim, class Value>
180 typename Value::return_type const & FieldFE<spacedim, Value>::value(const Point &p, const ElementAccessor<spacedim> &elm)
181 {
182  switch (elm.dim()) {
183  case 0:
184  return value_handler0_.value(p, elm);
185  case 1:
186  return value_handler1_.value(p, elm);
187  case 2:
188  return value_handler2_.value(p, elm);
189  case 3:
190  return value_handler3_.value(p, elm);
191  default:
192  ASSERT(false).error("Invalid element dimension!");
193  }
194 
195  return this->r_value_;
196 }
197 
198 
199 
200 /**
201  * Returns std::vector of scalar values in several points at once.
202  */
203 template <int spacedim, class Value>
206 {
207  ASSERT_EQ( point_list.size(), value_list.size() ).error();
208  ASSERT_DBG( point_list.n_rows() == spacedim && point_list.n_cols() == 1).error("Invalid point size.\n");
209 
210  switch (elm.dim()) {
211  case 0:
212  value_handler0_.value_list(point_list, elm, value_list);
213  break;
214  case 1:
215  value_handler1_.value_list(point_list, elm, value_list);
216  break;
217  case 2:
218  value_handler2_.value_list(point_list, elm, value_list);
219  break;
220  case 3:
221  value_handler3_.value_list(point_list, elm, value_list);
222  break;
223  default:
224  ASSERT(false).error("Invalid element dimension!");
225  }
226 }
227 
228 
229 
230 template <int spacedim, class Value>
232  ElementCacheMap &cache_map, unsigned int region_idx)
233 {
234  ASSERT( !boundary_dofs_ ).error("boundary field NOT supported!!\n");
235  std::shared_ptr<EvalPoints> eval_points = cache_map.eval_points();
237 
238  if (fe_values_.size() == 0) {
239  // initialize FEValues objects (when first using)
240  std::array<Quadrature, 4> quads{QGauss(0, 1), this->init_quad<1>(eval_points), this->init_quad<2>(eval_points), this->init_quad<3>(eval_points)};
241  fe_values_.resize(4);
242  fe_values_[0].initialize(quads[0], *fe_[0_d], update_values);
243  fe_values_[1].initialize(quads[1], *fe_[1_d], update_values);
244  fe_values_[2].initialize(quads[2], *fe_[2_d], update_values);
245  fe_values_[3].initialize(quads[3], *fe_[3_d], update_values);
246  }
247 
248  auto update_cache_data = cache_map.update_cache_data();
249  unsigned int region_in_cache = update_cache_data.region_cache_indices_range_.find(region_idx)->second;
250 
251  for (unsigned int i_elm=update_cache_data.region_element_cache_range_[region_in_cache];
252  i_elm<update_cache_data.region_element_cache_range_[region_in_cache+1]; ++i_elm) {
253  unsigned int elm_idx = cache_map.elm_idx_on_position(i_elm);
254  ElementAccessor<spacedim> elm(dh_->mesh(), elm_idx);
255  fe_values_[elm.dim()].reinit( elm );
256 
257  DHCellAccessor cell = dh_->cell_accessor_from_element( elm_idx );
258  LocDofVec loc_dofs = cell.get_loc_dof_indices();
259 
260  for (unsigned int i_ep=0; i_ep<eval_points->max_size(); ++i_ep) { // i_eval_point
261  //DHCellAccessor cache_cell = cache_map(cell);
262  int field_cache_idx = cache_map.get_field_value_cache_index(cache_map(cell).element_cache_index(), i_ep);
263  if (field_cache_idx < 0) continue; // skip
264  mat_value.fill(0.0);
265  for (unsigned int i_dof=0; i_dof<loc_dofs.n_elem; i_dof++) {
266  mat_value += data_vec_[loc_dofs[i_dof]] * this->handle_fe_shape(elm.dim(), i_dof, i_ep, 0);
267  }
268  data_cache.data().set(field_cache_idx) = mat_value;
269  }
270  }
271 }
272 
273 
274 template <int spacedim, class Value>
275 template <unsigned int dim>
276 Quadrature FieldFE<spacedim, Value>::init_quad(std::shared_ptr<EvalPoints> eval_points)
277 {
278  Quadrature quad(dim, eval_points->size(dim));
279  for (unsigned int k=0; k<eval_points->size(dim); k++)
280  quad.set(k) = eval_points->local_point<dim>(k);
281  return quad;
282 }
283 
284 
285 template <int spacedim, class Value>
287  this->init_unit_conversion_coefficient(rec, init_data);
288  this->in_rec_ = rec;
289  flags_ = init_data.flags_;
290 
291 
292  // read data from input record
293  reader_file_ = FilePath( rec.val<FilePath>("mesh_data_file") );
294  field_name_ = rec.val<std::string>("field_name");
295  if (! rec.opt_val<OutputTime::DiscreteSpace>("input_discretization", discretization_) ) {
296  discretization_ = OutputTime::DiscreteSpace::UNDEFINED;
297  }
298  if (! rec.opt_val<DataInterpolation>("interpolation", interpolation_) ) {
299  interpolation_ = DataInterpolation::equivalent_msh;
300  }
301  if (! rec.opt_val("default_value", default_value_) ) {
302  default_value_ = numeric_limits<double>::signaling_NaN();
303  }
304 }
305 
306 
307 
308 template <int spacedim, class Value>
309 void FieldFE<spacedim, Value>::set_mesh(const Mesh *mesh, bool boundary_domain) {
310  // Mesh can be set only for field initialized from input.
312  ASSERT(field_name_ != "").error("Uninitialized FieldFE, did you call init_from_input()?\n");
313  this->boundary_domain_ = boundary_domain;
314  switch (this->interpolation_) {
315  case DataInterpolation::identic_msh:
317  break;
318  case DataInterpolation::equivalent_msh:
319  if (!ReaderCache::check_compatible_mesh( reader_file_, const_cast<Mesh &>(*mesh) )) {
320  this->interpolation_ = DataInterpolation::gauss_p0;
321  WarningOut().fmt("Source mesh of FieldFE '{}' is not compatible with target mesh.\nInterpolation of input data will be changed to 'P0_gauss'.\n",
322  field_name_);
323  }
324  break;
325  case DataInterpolation::gauss_p0:
326  {
327  auto source_mesh = ReaderCache::get_mesh(reader_file_);
328  ReaderCache::get_element_ids(reader_file_, *(source_mesh.get()) );
329  break;
330  }
331  case DataInterpolation::interp_p0:
332  {
333  auto source_msh = ReaderCache::get_mesh(reader_file_);
334  ReaderCache::get_element_ids(reader_file_, *(source_msh.get()) );
335  if (!boundary_domain) {
336  this->interpolation_ = DataInterpolation::gauss_p0;
337  WarningOut().fmt("Interpolation 'P0_intersection' of FieldFE '{}' can't be used on bulk region.\nIt will be changed to 'P0_gauss'.\n", field_name_);
338  }
339  break;
340  }
341  }
342  this->make_dof_handler(mesh);
343  }
344 }
345 
346 
347 
348 template <int spacedim, class Value>
350  ASSERT(this->boundary_domain_);
351 
352  auto bc_mesh = dh_->mesh()->get_bc_mesh();
353  unsigned int n_comp = this->value_.n_rows() * this->value_.n_cols();
354  boundary_dofs_ = std::make_shared< std::vector<IntIdx> >( n_comp * bc_mesh->n_elements() );
355  std::vector<IntIdx> &in_vec = *( boundary_dofs_.get() );
356  unsigned int j = 0; // actual index to boundary_dofs_ vector
357 
358  for (auto ele : bc_mesh->elements_range()) {
359  IntIdx elm_shift = n_comp * ele.idx();
360  for (unsigned int i=0; i<n_comp; ++i, ++j) {
361  in_vec[j] = elm_shift + i;
362  }
363  }
364 
369 
371 }
372 
373 
374 
375 template <int spacedim, class Value>
377  // temporary solution - these objects will be set through FieldCommon
378  switch (this->value_.n_rows() * this->value_.n_cols()) { // by number of components
379  case 1: { // scalar
381  break;
382  }
383  case 3: { // vector
384  MixedPtr<FE_P_disc> fe_base(0) ;
385  fe_ = mixed_fe_system(fe_base, FEType::FEVector, 3);
386  break;
387  }
388  case 9: { // tensor
389  MixedPtr<FE_P_disc> fe_base(0) ;
390  fe_ = mixed_fe_system(fe_base, FEType::FETensor, 9);
391  break;
392  }
393  default:
394  ASSERT(false).error("Should not happen!\n");
395  }
396 
397  std::shared_ptr<DOFHandlerMultiDim> dh_par = std::make_shared<DOFHandlerMultiDim>( const_cast<Mesh &>(*mesh) );
398  std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( &const_cast<Mesh &>(*mesh), fe_);
399  dh_par->distribute_dofs(ds);
400  dh_ = dh_par;
401  unsigned int ndofs = dh_->max_elem_dofs();
402 
403  if (this->boundary_domain_) fill_boundary_dofs(); // temporary solution for boundary mesh
404  else data_vec_ = VectorMPI::sequential( dh_->lsize() ); // allocate data_vec_
405 
406  // initialization data of value handlers
407  FEValueInitData init_data;
408  init_data.dh = dh_;
409  init_data.data_vec = data_vec_;
410  init_data.ndofs = ndofs;
411  init_data.n_comp = this->n_comp();
412  init_data.comp_index = 0;
413 
414  // initialize value handler objects
415  value_handler0_.initialize(init_data);
416  value_handler1_.initialize(init_data);
417  value_handler2_.initialize(init_data);
418  value_handler3_.initialize(init_data);
419 }
420 
421 
422 
423 template <int spacedim, class Value>
425  // Time can be set only for field initialized from input.
427  ASSERT(field_name_ != "").error("Uninitialized FieldFE, did you call init_from_input()?\n");
428  ASSERT_PTR(dh_)(field_name_).error("Null target mesh pointer of finite element field, did you call set_mesh()?\n");
429  if ( reader_file_ == FilePath() ) return false;
430 
431  unsigned int n_components = this->value_.n_rows() * this->value_.n_cols();
432  double time_unit_coef = time.read_coef(in_rec_.find<string>("time_unit"));
433  double time_shift = time.read_time( in_rec_.find<Input::Tuple>("read_time_shift") );
434  double read_time = (time.end()+time_shift) / time_unit_coef;
435  BaseMeshReader::HeaderQuery header_query(field_name_, read_time, this->discretization_, dh_->hash());
436  ReaderCache::get_reader(reader_file_)->find_header(header_query);
437  // TODO: use default and check NaN values in data_vec
438 
439  unsigned int n_entities;
440  bool is_native = (header_query.discretization == OutputTime::DiscreteSpace::NATIVE_DATA);
441  bool boundary;
442  if (is_native || this->interpolation_==DataInterpolation::identic_msh || this->interpolation_==DataInterpolation::equivalent_msh) {
443  n_entities = dh_->mesh()->n_elements();
444  boundary = this->boundary_domain_;
445  } else {
446  n_entities = ReaderCache::get_mesh(reader_file_)->n_elements();
447  boundary = false;
448  }
449  auto input_data_cache = ReaderCache::get_reader(reader_file_)->template get_element_data<double>(n_entities, n_components,
450  boundary, this->component_idx_);
451  CheckResult checked_data = ReaderCache::get_reader(reader_file_)->scale_and_check_limits(field_name_,
453 
454 
455  if (checked_data == CheckResult::not_a_number) {
456  THROW( ExcUndefElementValue() << EI_Field(field_name_) );
457  }
458 
459  if (is_native) {
460  this->calculate_native_values(input_data_cache);
461  } else if (this->interpolation_==DataInterpolation::identic_msh || this->interpolation_==DataInterpolation::equivalent_msh) {
462  this->calculate_elementwise_values(input_data_cache);
463  } else if (this->interpolation_==DataInterpolation::gauss_p0) {
464  this->interpolate_gauss(input_data_cache);
465  } else { // DataInterpolation::interp_p0
466  this->interpolate_intersection(input_data_cache);
467  }
468 
469  return true;
470  } else return false;
471 
472 }
473 
474 
475 template <int spacedim, class Value>
477 {
478  static const unsigned int quadrature_order = 4; // parameter of quadrature
479  std::shared_ptr<Mesh> source_mesh = ReaderCache::get_mesh(reader_file_);
480  std::vector<unsigned int> searched_elements; // stored suspect elements in calculating the intersection
481  std::vector<arma::vec::fixed<3>> q_points; // real coordinates of quadrature points
482  std::vector<double> q_weights; // weights of quadrature points
483  unsigned int quadrature_size=0; // size of quadrature point and weight vector
484  std::vector<double> sum_val(dh_->max_elem_dofs()); // sum of value of one quadrature point
485  unsigned int elem_count; // count of intersect (source) elements of one quadrature point
486  std::vector<double> elem_value(dh_->max_elem_dofs()); // computed value of one (target) element
487  bool contains; // sign if source element contains quadrature point
488 
489  {
490  // set size of vectors to maximal count of quadrature points
491  QGauss quad(3, quadrature_order);
492  q_points.resize(quad.size());
493  q_weights.resize(quad.size());
494  }
495 
496  for (auto cell : dh_->own_range()) {
497  auto ele = cell.elm();
498  std::fill(elem_value.begin(), elem_value.end(), 0.0);
499  switch (cell.dim()) {
500  case 0:
501  quadrature_size = 1;
502  q_points[0] = *ele.node(0);
503  q_weights[0] = 1.0;
504  break;
505  case 1:
506  quadrature_size = value_handler1_.compute_quadrature(q_points, q_weights, ele, quadrature_order);
507  break;
508  case 2:
509  quadrature_size = value_handler2_.compute_quadrature(q_points, q_weights, ele, quadrature_order);
510  break;
511  case 3:
512  quadrature_size = value_handler3_.compute_quadrature(q_points, q_weights, ele, quadrature_order);
513  break;
514  }
515  searched_elements.clear();
516  source_mesh->get_bih_tree().find_bounding_box(ele.bounding_box(), searched_elements);
517 
518  for (unsigned int i=0; i<quadrature_size; ++i) {
519  std::fill(sum_val.begin(), sum_val.end(), 0.0);
520  elem_count = 0;
521  for (std::vector<unsigned int>::iterator it = searched_elements.begin(); it!=searched_elements.end(); it++) {
522  ElementAccessor<3> elm = source_mesh->element_accessor(*it);
523  contains=false;
524  switch (elm->dim()) {
525  case 0:
526  contains = arma::norm(*elm.node(0) - q_points[i], 2) < 4*std::numeric_limits<double>::epsilon();
527  break;
528  case 1:
529  contains = MappingP1<1,3>::contains_point(q_points[i], elm);
530  break;
531  case 2:
532  contains = MappingP1<2,3>::contains_point(q_points[i], elm);
533  break;
534  case 3:
535  contains = MappingP1<3,3>::contains_point(q_points[i], elm);
536  break;
537  default:
538  ASSERT(false).error("Invalid element dimension!");
539  }
540  if ( contains ) {
541  // projection point in element
542  unsigned int index = sum_val.size() * (*it);
543  for (unsigned int j=0; j < sum_val.size(); j++) {
544  sum_val[j] += (*data_vec)[index+j];
545  }
546  ++elem_count;
547  }
548  }
549 
550  if (elem_count > 0) {
551  for (unsigned int j=0; j < sum_val.size(); j++) {
552  elem_value[j] += (sum_val[j] / elem_count) * q_weights[i];
553  }
554  }
555  }
556 
557  LocDofVec loc_dofs;
558  if (this->boundary_domain_) loc_dofs = value_handler1_.get_loc_dof_indices(cell.elm_idx());
559  else loc_dofs = cell.get_loc_dof_indices();
560 
561  ASSERT_LE_DBG(loc_dofs.n_elem, elem_value.size());
562  for (unsigned int i=0; i < elem_value.size(); i++) {
563  ASSERT_LT_DBG( loc_dofs[i], (int)data_vec_.size());
564  data_vec_[loc_dofs[i]] = elem_value[i] * this->unit_conversion_coefficient_;
565  }
566  }
567 }
568 
569 
570 template <int spacedim, class Value>
572 {
573  std::shared_ptr<Mesh> source_mesh = ReaderCache::get_mesh(reader_file_);
574  std::vector<unsigned int> searched_elements; // stored suspect elements in calculating the intersection
575  std::vector<double> value(dh_->max_elem_dofs());
576  double total_measure;
577  double measure = 0;
578 
579  Mesh *mesh;
580  if (this->boundary_domain_) mesh = dh_->mesh()->get_bc_mesh();
581  else mesh = dh_->mesh();
582  for (auto elm : mesh->elements_range()) {
583  if (elm.dim() == 3) {
584  xprintf(Err, "Dimension of element in target mesh must be 0, 1 or 2! elm.idx() = %d\n", elm.idx());
585  }
586 
587  double epsilon = 4* numeric_limits<double>::epsilon() * elm.measure();
588 
589  // gets suspect elements
590  if (elm.dim() == 0) {
591  searched_elements.clear();
592  source_mesh->get_bih_tree().find_point(*elm.node(0), searched_elements);
593  } else {
594  BoundingBox bb = elm.bounding_box();
595  searched_elements.clear();
596  source_mesh->get_bih_tree().find_bounding_box(bb, searched_elements);
597  }
598 
599  // set zero values of value object
600  std::fill(value.begin(), value.end(), 0.0);
601  total_measure=0.0;
602 
603  START_TIMER("compute_pressure");
604  ADD_CALLS(searched_elements.size());
605 
606 
607  for (std::vector<unsigned int>::iterator it = searched_elements.begin(); it!=searched_elements.end(); it++)
608  {
609  ElementAccessor<3> ele = source_mesh->element_accessor(*it);
610  if (ele->dim() == 3) {
611  // get intersection (set measure = 0 if intersection doesn't exist)
612  switch (elm.dim()) {
613  case 0: {
614  arma::vec::fixed<3> real_point = *elm.node(0);
615  arma::mat::fixed<3, 4> elm_map = MappingP1<3,3>::element_map(ele);
616  arma::vec::fixed<4> unit_point = MappingP1<3,3>::project_real_to_unit(real_point, elm_map);
617 
618  measure = (std::fabs(arma::sum( unit_point )-1) <= 1e-14
619  && arma::min( unit_point ) >= 0)
620  ? 1.0 : 0.0;
621  break;
622  }
623  case 1: {
625  ComputeIntersection<1,3> CI(elm, ele, source_mesh.get());
626  CI.init();
627  CI.compute(is);
628 
629  IntersectionLocal<1,3> ilc(is);
630  measure = ilc.compute_measure() * elm.measure();
631  break;
632  }
633  case 2: {
635  ComputeIntersection<2,3> CI(elm, ele, source_mesh.get());
636  CI.init();
637  CI.compute(is);
638 
639  IntersectionLocal<2,3> ilc(is);
640  measure = 2 * ilc.compute_measure() * elm.measure();
641  break;
642  }
643  }
644 
645  //adds values to value_ object if intersection exists
646  if (measure > epsilon) {
647  unsigned int index = value.size() * (*it);
648  std::vector<double> &vec = *( data_vec.get() );
649  for (unsigned int i=0; i < value.size(); i++) {
650  value[i] += vec[index+i] * measure;
651  }
652  total_measure += measure;
653  }
654  }
655  }
656 
657  // computes weighted average, store it to data vector
658  if (total_measure > epsilon) {
660 
661  LocDofVec loc_dofs;
662  if (this->boundary_domain_) loc_dofs = value_handler1_.get_loc_dof_indices(elm.idx());
663  else{
664  DHCellAccessor cell = dh_->cell_accessor_from_element(elm.idx());
665  loc_dofs = cell.get_loc_dof_indices();
666  }
667 
668  ASSERT_LE_DBG(loc_dofs.n_elem, value.size());
669  for (unsigned int i=0; i < value.size(); i++) {
670  (*data_vector)[ loc_dofs[i] ] = value[i] / total_measure;
671  }
672  } else {
673  WarningOut().fmt("Processed element with idx {} is out of source mesh!\n", elm.idx());
674  }
675  END_TIMER("compute_pressure");
676 
677  }
678 }
679 
680 
681 template <int spacedim, class Value>
683 {
684  // Same algorithm as in output of Node_data. Possibly code reuse.
685  unsigned int dof_size, data_vec_i;
686  std::vector<unsigned int> count_vector(data_vec_.size(), 0);
689  std::vector<LongIdx> global_dof_indices(dh_->max_elem_dofs());
690 
691  // iterate through cells, assembly MPIVector
692  for (auto cell : dh_->own_range()) {
693  dof_size = cell.get_dof_indices(global_dof_indices);
694  LocDofVec loc_dofs = cell.get_loc_dof_indices();
695  data_vec_i = cell.elm_idx() * dof_size;
696  ASSERT_EQ_DBG(dof_size, loc_dofs.n_elem);
697  for (unsigned int i=0; i<dof_size; ++i, ++data_vec_i) {
698  (*data_vector)[ loc_dofs[i] ] += (*data_cache)[ global_dof_indices[i] ];
699  ++count_vector[ loc_dofs[i] ];
700  }
701  }
702 
703  // compute averages of values
704  for (unsigned int i=0; i<data_vec_.size(); ++i) {
705  if (count_vector[i]>0) (*data_vector)[i] /= count_vector[i];
706  }
707 }
708 
709 
710 template <int spacedim, class Value>
712 {
713  // Same algorithm as in output of Node_data. Possibly code reuse.
714  unsigned int data_vec_i;
715  std::vector<unsigned int> count_vector(data_vec_.size(), 0);
717 
718  // iterate through elements, assembly global vector and count number of writes
719  if (this->boundary_domain_) {
720  Mesh *mesh = dh_->mesh()->get_bc_mesh();
721  for (auto ele : mesh->elements_range()) { // remove special case for rank == 0 - necessary for correct output
722  LocDofVec loc_dofs = value_handler1_.get_loc_dof_indices(ele.idx());
723  data_vec_i = ele.idx() * dh_->max_elem_dofs();
724  for (unsigned int i=0; i<loc_dofs.n_elem; ++i, ++data_vec_i) {
725  ASSERT_LT_DBG(loc_dofs[i], (LongIdx)data_vec_.size());
726  data_vec_[ loc_dofs[i] ] += (*data_cache)[data_vec_i];
727  ++count_vector[ loc_dofs[i] ];
728  }
729  }
730  }
731  else {
732  // iterate through cells, assembly global vector and count number of writes - prepared solution for further development
733  for (auto cell : dh_->own_range()) {
734  LocDofVec loc_dofs = cell.get_loc_dof_indices();
735  data_vec_i = cell.elm_idx() * dh_->max_elem_dofs();
736  for (unsigned int i=0; i<loc_dofs.n_elem; ++i, ++data_vec_i) {
737  ASSERT_LT_DBG(loc_dofs[i], (LongIdx)data_vec_.size());
738  data_vec_[ loc_dofs[i] ] += (*data_cache)[data_vec_i];
739  ++count_vector[ loc_dofs[i] ];
740  }
741  }
742  }
743 
744  // compute averages of values
745  for (unsigned int i=0; i<data_vec_.size(); ++i) {
746  if (count_vector[i]>0) data_vec_[i] /= count_vector[i];
747  }
748 }
749 
750 
751 template <int spacedim, class Value>
753  ASSERT_EQ(output_data_cache.n_values() * output_data_cache.n_comp(), dh_->distr()->lsize()).error();
754  double loc_values[output_data_cache.n_comp()];
755  unsigned int i;
756 
758  for (auto dh_cell : dh_->own_range()) {
759  LocDofVec loc_dofs = dh_cell.get_loc_dof_indices();
760  for (i=0; i<loc_dofs.n_elem; ++i) loc_values[i] = (*data_vec)[ loc_dofs[i] ];
761  for ( ; i<output_data_cache.n_comp(); ++i) loc_values[i] = numeric_limits<double>::signaling_NaN();
762  output_data_cache.store_value( dh_cell.local_idx(), loc_values );
763  }
764 
765  output_data_cache.set_dof_handler_hash( dh_->hash() );
766 }
767 
768 
769 
770 template <int spacedim, class Value>
771 inline unsigned int FieldFE<spacedim, Value>::data_size() const {
772  return data_vec_.size();
773 }
774 
775 
776 
777 template <int spacedim, class Value>
780 }
781 
782 
783 
784 template <int spacedim, class Value>
787 }
788 
789 
790 
791 template <int spacedim, class Value>
793  unsigned int i_dof, unsigned int i_qp, unsigned int comp_index)
794 {
796  for (unsigned int c=0; c<Value::NRows_*Value::NCols_; ++c)
797  v(c/spacedim,c%spacedim) = fe_values_[dim].shape_value_component(i_dof, i_qp, comp_index+c);
798  if (Value::NRows_ == Value::NCols_)
799  return v;
800  else
801  return v.t();
802 }
803 
804 
805 
806 template <int spacedim, class Value>
808 {}
809 
810 
811 // Instantiations of FieldFE
unsigned int comp_index
index of component (of vector_value/tensor_value)
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
Shape function values.
Definition: update_flags.hh:87
virtual ~FieldFE()
Destructor.
Definition: field_fe.cc:807
void init_unit_conversion_coefficient(const Input::Record &rec, const struct FieldAlgoBaseInitData &init_data)
Init value of unit_conversion_coefficient_ from input.
unsigned int size() const
Return size of output data.
Definition: vector_mpi.hh:98
Bounding box in 3d ambient space.
Definition: bounding_box.hh:54
void set_dof_handler_hash(std::size_t hash)
void interpolate_gauss(ElementDataCache< double >::ComponentDataPtr data_vec)
Interpolate data (use Gaussian distribution) over all elements of target mesh.
Definition: field_fe.cc:476
unsigned int size() const
Definition: armor.hh:718
FEValueHandler< 0, spacedim, Value > value_handler0_
Value handler that allows get value of 0D elements.
Definition: field_fe.hh:196
Armor::Array< double >::ArrayMatSet set(uint i)
Definition: quadrature.hh:98
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:332
std::shared_ptr< std::vector< T > > ComponentDataPtr
ArmaVec< double, N > vec
Definition: armor.hh:861
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
MixedPtr< FiniteElement > fe_
Definition: field_fe.hh:209
double read_coef(Input::Iterator< std::string > unit_it) const
Classes with algorithms for computation of intersections of meshes.
Some value(s) is set to NaN.
void initialize(FEValueInitData init_data)
Initialize data members.
std::vector< FEValues< spacedim > > fe_values_
List of FEValues objects of dimensions 0,1,2,3 used for value calculation.
Definition: field_fe.hh:243
unsigned int n_values() const
unsigned int component_idx_
Specify if the field is part of a MultiField and which component it is.
LocDofVec get_loc_dof_indices(unsigned int cell_idx) const
TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
static void get_element_ids(const FilePath &file_path, const Mesh &mesh)
Definition: reader_cache.cc:80
const UpdateCacheHelper & update_cache_data() const
Return update cache data helper.
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
std::string field_name_
field name read from input
Definition: field_fe.hh:215
void local_to_ghost_begin()
local_to_ghost_{begin,end} updates the ghost values on neighbouring processors from local values ...
Definition: vector_mpi.hh:123
double compute_measure() const override
Computes the relative measure of intersection object.
int IntIdx
Definition: index_types.hh:25
void zero_entries()
Definition: vector_mpi.hh:82
Abstract linear system class.
Definition: balance.hh:40
void value_list(const Armor::array &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list)
Returns std::vector of scalar values in several points at once.
FieldFlag::Flags flags_
FEValueHandler< 2, spacedim, Value > value_handler2_
Value handler that allows get value of 2D elements.
Definition: field_fe.hh:200
VectorMPI data_vec_
Store data of Field.
Definition: field_fe.hh:193
void store_value(unsigned int idx, const T *value)
Fundamental simplicial intersections.
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
static BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map)
Definition: mapping_p1.cc:69
FilePath reader_file_
mesh reader file
Definition: field_fe.hh:212
#define INSTANCE_ALL(field)
void init()
Initializes lower dimensional objects. Sets correctly the pointers to Plucker coordinates and product...
Directing class of FieldValueCache.
Definition: mesh.h:78
Iterator< Ret > find(const string &key) const
void initialize(FEValueInitData init_data)
Initialize data members.
Cell accessor allow iterate over DOF handler cells.
Helper struct stores data for initizalize descentants of FieldAlgorithmBase.
Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Returns one value in one given point.
Class FESystem for compound finite elements.
VectorDataPtr data_ptr()
Getter for shared pointer of output data.
Definition: vector_mpi.hh:75
uint n_cols() const
Definition: armor.hh:710
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
void set_mesh(const Mesh *mesh, bool boundary_domain) override
Definition: field_fe.cc:309
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:347
#define ADD_CALLS(n_calls)
Increase number of calls in actual timer.
double default_value_
Default value of element if not set in mesh data file.
Definition: field_fe.hh:227
FEValueHandler< 3, spacedim, Value > value_handler3_
Value handler that allows get value of 3D elements.
Definition: field_fe.hh:202
bool set_time(const TimeStep &time) override
Definition: field_fe.cc:424
static const Input::Type::Selection & get_interp_selection_input_type()
Definition: field_fe.cc:103
Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Returns one value in one given point.
Value::return_type r_value_
void fill_boundary_dofs()
Definition: field_fe.cc:349
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
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
Symmetric Gauss-Legendre quadrature formulae on simplices.
unsigned int data_size() const
Definition: field_fe.cc:771
constexpr bool match(Mask mask) const
Definition: flag_array.hh:163
void init()
Initializes lower dimensional objects. Sets correctly the pointers to Plucker coordinates and product...
double unit_conversion_coefficient_
Coeficient of conversion of user-defined unit.
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:196
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
bool boundary_domain_
Is set in set_mesh method. Value true means, that we accept only boundary element accessors in the va...
Definition: field_fe.hh:233
bool opt_val(const string &key, Ret &value) const
OutputTime::DiscreteSpace discretization_
Specify section where to find the field data in input mesh file.
Definition: field_fe.hh:218
std::shared_ptr< DOFHandlerMultiDim > dh_
DOF handler object.
Definition: field_fe.hh:191
Quadrature init_quad(std::shared_ptr< EvalPoints > eval_points)
Initialize FEValues object of given dimension.
Definition: field_fe.cc:276
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:541
Armor::ArmaMat< typename Value::element_type, Value::NRows_, Value::NCols_ > handle_fe_shape(unsigned int dim, unsigned int i_dof, unsigned int i_qp, unsigned int comp_index)
Definition: field_fe.cc:792
static const Input::Type::Tuple & get_input_time_type(double lower_bound=-std::numeric_limits< double >::max(), double upper_bound=std::numeric_limits< double >::max())
const Armor::Array< elm_type > & data() const
Return data vector.
DataInterpolation interpolation_
Specify type of FE data interpolation.
Definition: field_fe.hh:221
Definitions of basic Lagrangean finite elements with polynomial shape functions.
static FileName input()
The factory function for declaring type FileName for input files.
Definition: type_base.cc:526
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
const Ret val(const string &key) const
#define xprintf(...)
Definition: system.hh:93
void resize(unsigned int local_size)
Definition: vector_mpi.cc:50
std::unordered_map< unsigned int, unsigned int > region_cache_indices_range_
#define START_TIMER(tag)
Starts a timer with specified tag.
Selection & add_value(const int value, const std::string &key, const std::string &description="", TypeBase::attribute_map attributes=TypeBase::attribute_map())
Adds one new value with name given by key to the Selection.
void interpolate_intersection(ElementDataCache< double >::ComponentDataPtr data_vec)
Interpolate data (use intersection library) over all elements of target mesh.
Definition: field_fe.cc:571
std::shared_ptr< std::vector< IntIdx > > boundary_dofs_
Definition: field_fe.hh:240
void cache_update(FieldValueCache< typename Value::element_type > &data_cache, ElementCacheMap &cache_map, unsigned int region_idx) override
Definition: field_fe.cc:231
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
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Definition: field_fe.cc:180
ArrayMatSet set(uint index)
Definition: armor.hh:821
uint n_rows() const
Definition: armor.hh:705
virtual Range< ElementAccessor< 3 > > elements_range() const
Returns range of bulk elements.
Definition: mesh.cc:1088
void local_to_ghost_data_scatter_begin()
Call begin scatter functions (local to ghost) on data vector.
Definition: field_fe.cc:778
static bool check_compatible_mesh(const FilePath &file_path, Mesh &mesh)
Definition: reader_cache.cc:73
double end() const
unsigned int ndofs
number of dofs
Space< spacedim >::Point Point
void local_to_ghost_data_scatter_end()
Call end scatter functions (local to ghost) on data vector.
Definition: field_fe.cc:785
void set_boundary_dofs_vector(std::shared_ptr< std::vector< IntIdx > > boundary_dofs)
TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
unsigned int compute_quadrature(std::vector< arma::vec::fixed< 3 >> &q_points, std::vector< double > &q_weights, const ElementAccessor< spacedim > &elm, unsigned int order=3)
Compute real coordinates and weights (use QGauss) for given element.
CheckResult
Return type of method that checked data stored in ElementDataCache (NaN values, limits) ...
std::shared_ptr< DOFHandlerMultiDim > dh
DOF handler object.
unsigned int elm_idx_on_position(unsigned pos) const
Return idx of element stored at given position of ElementCacheMap.
VectorMPI set_fe_data(std::shared_ptr< DOFHandlerMultiDim > dh, unsigned int component_index=0, VectorMPI dof_values=VectorMPI::sequential(0))
Definition: field_fe.cc:141
double read_time(Input::Iterator< Input::Tuple > time_it, double default_time=std::numeric_limits< double >::quiet_NaN()) const
typename arma::Mat< Type >::template fixed< nr, nc > ArmaMat
Definition: armor.hh:502
void local_to_ghost_end()
local_to_ghost_{begin,end} updates the ghost values on neighbouring processors from local values ...
Definition: vector_mpi.hh:127
#define ASSERT_LE_DBG(a, b)
Definition of comparative assert macro (Less or Equal) only for debug mode.
Definition: asserts.hh:308
FEValueHandler< 1, spacedim, Value > value_handler1_
Value handler that allows get value of 1D elements.
Definition: field_fe.hh:198
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:216
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
int get_field_value_cache_index(unsigned int elm_idx, unsigned int loc_point_idx) const
Return index of point in FieldValueCache.
static const Input::Type::Selection & get_disc_selection_input_type()
Definition: field_fe.cc:91
static Default read_time(const std::string &description)
The factory function to make an default value that will be specified at the time when a key will be r...
Definition: type_record.hh:97
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: index_types.hh:24
bool is_constant_in_space_
Flag detects that field is only dependent on time.
void value_list(const Armor::array &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list)
Returns std::vector of scalar values in several points at once.
Value value_
Last value, prevents passing large values (vectors) by value.
Input::Record in_rec_
Accessor to Input::Record.
Definition: field_fe.hh:230
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
Definition: asserts.hh:336
const Selection & close() const
Close the Selection, no more values can be added.
FieldFE(unsigned int n_comp=0)
Definition: field_fe.cc:132
#define ASSERT_DBG(expr)
std::shared_ptr< EvalPoints > eval_points() const
Getter of eval_points object.
MixedPtr< FESystem > mixed_fe_system(MixedPtr< FiniteElement > fe, Args &&...args)
Definition: fe_system.hh:171
Definitions of particular quadrature rules on simplices.
Class for 2D-2D intersections.
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:258
void set_boundary_dofs_vector(std::shared_ptr< std::vector< IntIdx > > boundary_dofs)
TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
#define END_TIMER(tag)
Ends a timer with specified tag.
DataInterpolation
Definition: field_fe.hh:58
Definition: system.hh:65
virtual void value_list(const Armor::array &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list)
Definition: field_fe.cc:204
virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData &init_data)
Definition: field_fe.cc:286
static std::shared_ptr< Mesh > get_mesh(const FilePath &file_path)
Definition: reader_cache.cc:40
Record type proxy class.
Definition: type_record.hh:182
VectorMPI data_vec
Store data of Field.
Accessor to the data with type Type::Tuple.
Definition: accessors.hh:412
static ElementMap element_map(ElementAccessor< 3 > elm)
Definition: mapping_p1.cc:48
unsigned int n_comp() const
static constexpr Mask equation_input
The field is data parameter of the owning equation. (default on)
Definition: field_flag.hh:33
unsigned int dim() const
Definition: accessors.hh:157
void calculate_elementwise_values(ElementDataCache< double >::ComponentDataPtr data_cache)
Calculate elementwise data over all elements of target mesh.
Definition: field_fe.cc:711
static bool contains_point(arma::vec point, ElementAccessor< 3 > elm)
Test if element contains given point.
Definition: mapping_p1.cc:96
Internal auxiliary class representing intersection object of simplex<dimA> and simplex<dimB>.
static std::shared_ptr< BaseMeshReader > get_reader(const FilePath &file_path)
Definition: reader_cache.cc:36
std::shared_ptr< VectorData > VectorDataPtr
Definition: vector_mpi.hh:46
unsigned int n_comp
number of components
unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:87
void calculate_native_values(ElementDataCache< double >::ComponentDataPtr data_cache)
Calculate native data over all elements of target mesh.
Definition: field_fe.cc:682
FieldFlag::Flags flags_
Field flags.
Definition: field_fe.hh:224
Initialization structure of FEValueHandler class.
Class for declaration of the input data that are in string format.
Definition: type_base.hh:589
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
Definition: exceptions.hh:53
Representation of one time step..
Template for classes storing finite set of named values.
void make_dof_handler(const Mesh *mesh)
Create DofHandler object.
Definition: field_fe.cc:376
Implementation of range helper class.
void native_data_to_cache(ElementDataCache< double > &output_data_cache)
Definition: field_fe.cc:752
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
Definition: global_defs.h:180
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual)
Definition: asserts.hh:328
static constexpr Mask declare_input
The field can be set from input. The key in input field descriptor is declared. (default on) ...
Definition: field_flag.hh:35
static VectorMPI sequential(unsigned int size)
Definition: vector_mpi.cc:44
Definition: field.hh:60
unsigned int n_comp() const
Internal class representing intersection object.
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:300