Flow123d  release_3.0.0-680-gbed5aba
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 "io/reader_cache.hh"
30 #include "io/msh_gmshreader.h"
31 #include "mesh/accessors.hh"
32 #include "mesh/range_wrapper.hh"
33 #include "mesh/bc_mesh.hh"
35 
36 #include "system/sys_profiler.hh"
40 
41 
42 
43 
44 /// Implementation.
45 
46 namespace it = Input::Type;
47 
48 
49 
50 
52 
53 
54 template <int spacedim, class Value>
55 const Input::Type::Record & FieldFE<spacedim, Value>::get_input_type()
56 {
57  return it::Record("FieldFE", FieldAlgorithmBase<spacedim,Value>::template_name()+" Field given by finite element approximation.")
61  "GMSH mesh with data. Can be different from actual computational mesh.")
63  "Section where to find the field.\n Some sections are specific to file format: "
64  "point_data/node_data, cell_data/element_data, -/element_node_data, native/-.\n"
65  "If not given by a user, we try to find the field in all sections, but we report an error "
66  "if it is found in more than one section.")
68  "The values of the Field are read from the ```$ElementData``` section with field name given by this key.")
69  .declare_key("default_value", IT::Double(), IT::Default::optional(),
70  "Default value is set on elements which values have not been listed in the mesh data file.")
71  .declare_key("time_unit", IT::String(), IT::Default::read_time("Common unit of TimeGovernor."),
72  "Definition of the unit of all times defined in the mesh data file.")
73  .declare_key("read_time_shift", TimeGovernor::get_input_time_type(), IT::Default("0.0"),
74  "This key allows reading field data from the mesh data file shifted in time. Considering the time 't', field descriptor with time 'T', "
75  "time shift 'S', then if 't > T', we read the time frame 't + S'.")
77  IT::Default("\"equivalent_mesh\""), "Type of interpolation applied to the input spatial data.\n"
78  "The default value 'equivalent_mesh' assumes the data being constant on elements living on the same mesh "
79  "as the computational mesh, but possibly with different numbering. In the case of the same numbering, "
80  "the user can set 'identical_mesh' to omit algorithm for guessing node and element renumbering. "
81  "Alternatively, in case of different input mesh, several interpolation algorithms are available.")
82  .close();
83 }
84 
85 template <int spacedim, class Value>
87 {
88  return it::Selection("FE_discretization",
89  "Specify the section in mesh input file where field data is listed.\nSome sections are specific to file format.")
90  .add_value(OutputTime::DiscreteSpace::NODE_DATA, "node_data", "point_data (VTK) / node_data (GMSH)")
91  .add_value(OutputTime::DiscreteSpace::ELEM_DATA, "element_data", "cell_data (VTK) / element_data (GMSH)")
92  .add_value(OutputTime::DiscreteSpace::CORNER_DATA, "element_node_data", "element_node_data (only for GMSH)")
93  .add_value(OutputTime::DiscreteSpace::NATIVE_DATA, "native_data", "native_data (only for VTK)")
94  .close();
95 }
96 
97 template <int spacedim, class Value>
99 {
100  return it::Selection("interpolation", "Specify interpolation of the input data from its input mesh to the computational mesh.")
101  .add_value(DataInterpolation::identic_msh, "identic_mesh", "Topology and indices of nodes and elements of"
102  "the input mesh and the computational mesh are identical. "
103  "This interpolation is typically used for GMSH input files containing only the field values without "
104  "explicit mesh specification.")
105  .add_value(DataInterpolation::equivalent_msh, "equivalent_mesh", "Topologies of the input mesh and the computational mesh "
106  "are the same, the node and element numbering may differ. "
107  "This interpolation can be used also for VTK input data.") // default value
108  .add_value(DataInterpolation::gauss_p0, "P0_gauss", "Topologies of the input mesh and the computational mesh may differ. "
109  "Constant values on the elements of the computational mesh are evaluated using the Gaussian quadrature of the fixed order 4, "
110  "where the quadrature points and their values are found in the input mesh and input data using the BIH tree search."
111  )
112  .add_value(DataInterpolation::interp_p0, "P0_intersection", "Topologies of the input mesh and the computational mesh may differ. "
113  "Can be applied only for boundary fields. For every (boundary) element of the computational mesh the "
114  "intersection with the input mesh is computed. Constant values on the elements of the computational mesh "
115  "are evaluated as the weighted average of the (constant) values on the intersecting elements of the input mesh.")
116  .close();
117 }
118 
119 template <int spacedim, class Value>
121  Input::register_class< FieldFE<spacedim, Value>, unsigned int >("FieldFE") +
123 
124 
125 
126 template <int spacedim, class Value>
128 : FieldAlgorithmBase<spacedim, Value>(n_comp),
129  data_vec_(nullptr),
130  field_name_("")
131 {
132  this->is_constant_in_space_ = false;
133 }
134 
135 
136 template <int spacedim, class Value>
137 void FieldFE<spacedim, Value>::set_fe_data(std::shared_ptr<DOFHandlerMultiDim> dh,
138  MappingP1<1,3> *map1,
139  MappingP1<2,3> *map2,
140  MappingP1<3,3> *map3,
141  VectorMPI *data)
142 {
143  dh_ = dh;
144  data_vec_ = data;
145  reinit_fe_data(map1, map2, map3);
146 }
147 
148 
149 
150 template <int spacedim, class Value>
151 void FieldFE<spacedim, Value>::set_native_dh(std::shared_ptr<DOFHandlerMultiDim> dh)
152 {
153  dh_ = dh;
154  reinit_fe_data(nullptr, nullptr, nullptr);
155 }
156 
157 
158 
159 template <int spacedim, class Value>
161  MappingP1<2,3> *map2,
162  MappingP1<3,3> *map3)
163 {
164  //ASSERT_EQ(dh_->n_global_dofs(), data_vec_->size());
165 
166  unsigned int ndofs = dh_->max_elem_dofs();
167  dof_indices_.resize(ndofs);
168 
169  // initialization data of value handlers
170  FEValueInitData init_data;
171  init_data.dh = dh_;
172  init_data.data_vec = data_vec_;
173  init_data.ndofs = ndofs;
174  init_data.n_comp = this->n_comp();
175 
176  // initialize value handler objects
177  value_handler0_.initialize(init_data);
178  value_handler1_.initialize(init_data, map1);
179  value_handler2_.initialize(init_data, map2);
180  value_handler3_.initialize(init_data, map3);
181 
182  // set discretization
183  discretization_ = OutputTime::DiscreteSpace::UNDEFINED;
184  interpolation_ = DataInterpolation::equivalent_msh;
185 }
186 
187 
188 
189 /**
190  * Returns one value in one given point. ResultType can be used to avoid some costly calculation if the result is trivial.
191  */
192 template <int spacedim, class Value>
193 typename Value::return_type const & FieldFE<spacedim, Value>::value(const Point &p, const ElementAccessor<spacedim> &elm)
194 {
195  switch (elm.dim()) {
196  case 0:
197  return value_handler0_.value(p, elm);
198  case 1:
199  return value_handler1_.value(p, elm);
200  case 2:
201  return value_handler2_.value(p, elm);
202  case 3:
203  return value_handler3_.value(p, elm);
204  default:
205  ASSERT(false).error("Invalid element dimension!");
206  }
207 
208  return this->r_value_;
209 }
210 
211 
212 
213 /**
214  * Returns std::vector of scalar values in several points at once.
215  */
216 template <int spacedim, class Value>
219 {
220  ASSERT_EQ( point_list.size(), value_list.size() ).error();
221 
222  switch (elm.dim()) {
223  case 0:
224  value_handler0_.value_list(point_list, elm, value_list);
225  break;
226  case 1:
227  value_handler1_.value_list(point_list, elm, value_list);
228  break;
229  case 2:
230  value_handler2_.value_list(point_list, elm, value_list);
231  break;
232  case 3:
233  value_handler3_.value_list(point_list, elm, value_list);
234  break;
235  default:
236  ASSERT(false).error("Invalid element dimension!");
237  }
238 }
239 
240 
241 
242 template <int spacedim, class Value>
244  this->init_unit_conversion_coefficient(rec, init_data);
245  this->in_rec_ = rec;
246  flags_ = init_data.flags_;
247 
248 
249  // read data from input record
250  reader_file_ = FilePath( rec.val<FilePath>("mesh_data_file") );
251  field_name_ = rec.val<std::string>("field_name");
252  if (! rec.opt_val<OutputTime::DiscreteSpace>("input_discretization", discretization_) ) {
253  discretization_ = OutputTime::DiscreteSpace::UNDEFINED;
254  }
255  if (! rec.opt_val<DataInterpolation>("interpolation", interpolation_) ) {
256  interpolation_ = DataInterpolation::equivalent_msh;
257  }
258  if (! rec.opt_val("default_value", default_value_) ) {
259  default_value_ = numeric_limits<double>::signaling_NaN();
260  }
261 }
262 
263 
264 
265 template <int spacedim, class Value>
266 void FieldFE<spacedim, Value>::set_mesh(const Mesh *mesh, bool boundary_domain) {
267  // Mesh can be set only for field initialized from input.
269  ASSERT(field_name_ != "").error("Uninitialized FieldFE, did you call init_from_input()?\n");
270  this->boundary_domain_ = boundary_domain;
271  this->make_dof_handler(mesh);
272  switch (this->interpolation_) {
273  case DataInterpolation::identic_msh:
275  break;
276  case DataInterpolation::equivalent_msh:
277  if (!ReaderCache::check_compatible_mesh( reader_file_, const_cast<Mesh &>(*mesh) )) {
278  this->interpolation_ = DataInterpolation::gauss_p0;
279  WarningOut().fmt("Source mesh of FieldFE '{}' is not compatible with target mesh.\nInterpolation of input data will be changed to 'P0_gauss'.\n",
280  field_name_);
281  }
282  break;
283  case DataInterpolation::gauss_p0:
284  {
285  auto source_mesh = ReaderCache::get_mesh(reader_file_);
286  ReaderCache::get_element_ids(reader_file_, *(source_mesh.get()) );
287  break;
288  }
289  case DataInterpolation::interp_p0:
290  {
291  auto source_msh = ReaderCache::get_mesh(reader_file_);
292  ReaderCache::get_element_ids(reader_file_, *(source_msh.get()) );
293  if (!boundary_domain) {
294  this->interpolation_ = DataInterpolation::gauss_p0;
295  WarningOut().fmt("Interpolation 'P0_intersection' of FieldFE '{}' can't be used on bulk region.\nIt will be changed to 'P0_gauss'.\n", field_name_);
296  }
297  break;
298  }
299  }
300  if (boundary_domain) fill_boundary_dofs(); // temporary solution for boundary mesh
301  }
302 }
303 
304 
305 
306 template <int spacedim, class Value>
308  ASSERT(this->boundary_domain_);
309 
310  auto bc_mesh = dh_->mesh()->get_bc_mesh();
311  unsigned int n_comp = this->value_.n_rows() * this->value_.n_cols();
312  boundary_dofs_ = std::make_shared< std::vector<LongIdx> >( n_comp * bc_mesh->n_elements() );
313  std::vector<LongIdx> &in_vec = *( boundary_dofs_.get() );
314  unsigned int j = 0; // actual index to boundary_dofs_ vector
315 
316  for (auto ele : bc_mesh->elements_range()) {
317  LongIdx elm_shift = n_comp * ele.idx();
318  for (unsigned int i=0; i<n_comp; ++i, ++j) {
319  in_vec[j] = elm_shift + i;
320  }
321  }
322 
327 
328  data_vec_->resize(boundary_dofs_->size());
329 }
330 
331 
332 
333 template <int spacedim, class Value>
335  // temporary solution - these objects will be set through FieldCommon
336  switch (this->value_.n_rows() * this->value_.n_cols()) { // by number of components
337  case 1: { // scalar
338  fe0_ = new FE_P_disc<0>(0);
339  fe1_ = new FE_P_disc<1>(0);
340  fe2_ = new FE_P_disc<2>(0);
341  fe3_ = new FE_P_disc<3>(0);
342  break;
343  }
344  case 3: { // vector
345  std::shared_ptr< FiniteElement<0> > fe0_ptr = std::make_shared< FE_P_disc<0> >(0);
346  std::shared_ptr< FiniteElement<1> > fe1_ptr = std::make_shared< FE_P_disc<1> >(0);
347  std::shared_ptr< FiniteElement<2> > fe2_ptr = std::make_shared< FE_P_disc<2> >(0);
348  std::shared_ptr< FiniteElement<3> > fe3_ptr = std::make_shared< FE_P_disc<3> >(0);
349  fe0_ = new FESystem<0>(fe0_ptr, FEType::FEVector, 3);
350  fe1_ = new FESystem<1>(fe1_ptr, FEType::FEVector, 3);
351  fe2_ = new FESystem<2>(fe2_ptr, FEType::FEVector, 3);
352  fe3_ = new FESystem<3>(fe3_ptr, FEType::FEVector, 3);
353  break;
354  }
355  case 9: { // tensor
356  std::shared_ptr< FiniteElement<0> > fe0_ptr = std::make_shared< FE_P_disc<0> >(0);
357  std::shared_ptr< FiniteElement<1> > fe1_ptr = std::make_shared< FE_P_disc<1> >(0);
358  std::shared_ptr< FiniteElement<2> > fe2_ptr = std::make_shared< FE_P_disc<2> >(0);
359  std::shared_ptr< FiniteElement<3> > fe3_ptr = std::make_shared< FE_P_disc<3> >(0);
360  fe0_ = new FESystem<0>(fe0_ptr, FEType::FETensor, 9);
361  fe1_ = new FESystem<1>(fe1_ptr, FEType::FETensor, 9);
362  fe2_ = new FESystem<2>(fe2_ptr, FEType::FETensor, 9);
363  fe3_ = new FESystem<3>(fe3_ptr, FEType::FETensor, 9);
364  break;
365  }
366  default:
367  ASSERT(false).error("Should not happen!\n");
368  }
369 
370  DOFHandlerMultiDim dh_par( const_cast<Mesh &>(*mesh) );
371  std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( &const_cast<Mesh &>(*mesh), fe0_, fe1_, fe2_, fe3_);
372  dh_par.distribute_dofs(ds);
373  dh_ = dh_par.sequential();
374  unsigned int ndofs = dh_->max_elem_dofs();
375  dof_indices_.resize(ndofs);
376 
377  // allocate data_vec_
378  data_vec_ = VectorMPI::sequential( dh_->n_global_dofs() );
379 
380  // initialization data of value handlers
381  FEValueInitData init_data;
382  init_data.dh = dh_;
383  init_data.data_vec = data_vec_;
384  init_data.ndofs = ndofs;
385  init_data.n_comp = this->n_comp();
386 
387  // initialize value handler objects
388  value_handler0_.initialize(init_data);
389  value_handler1_.initialize(init_data);
390  value_handler2_.initialize(init_data);
391  value_handler3_.initialize(init_data);
392 }
393 
394 
395 
396 template <int spacedim, class Value>
398  // Time can be set only for field initialized from input.
400  ASSERT(field_name_ != "").error("Uninitialized FieldFE, did you call init_from_input()?\n");
401  ASSERT_PTR(dh_)(field_name_).error("Null target mesh pointer of finite element field, did you call set_mesh()?\n");
402  if ( reader_file_ == FilePath() ) return false;
403 
404  unsigned int n_components = this->value_.n_rows() * this->value_.n_cols();
405  double time_unit_coef = time.read_coef(in_rec_.find<string>("time_unit"));
406  double time_shift = time.read_time( in_rec_.find<Input::Tuple>("read_time_shift") );
407  double read_time = (time.end()+time_shift) / time_unit_coef;
408  BaseMeshReader::HeaderQuery header_query(field_name_, read_time, this->discretization_, dh_->hash());
409  ReaderCache::get_reader(reader_file_)->find_header(header_query);
410  // TODO: use default and check NaN values in data_vec
411 
412  unsigned int n_entities;
413  bool is_native = (header_query.discretization == OutputTime::DiscreteSpace::NATIVE_DATA);
414  bool boundary;
415  if (is_native || this->interpolation_==DataInterpolation::identic_msh || this->interpolation_==DataInterpolation::equivalent_msh) {
416  n_entities = dh_->mesh()->n_elements();
417  boundary = this->boundary_domain_;
418  } else {
419  n_entities = ReaderCache::get_mesh(reader_file_)->n_elements();
420  boundary = false;
421  }
422  auto data_vec = ReaderCache::get_reader(reader_file_)->template get_element_data<double>(n_entities, n_components,
423  boundary, this->component_idx_);
424  CheckResult checked_data = ReaderCache::get_reader(reader_file_)->scale_and_check_limits(field_name_,
426 
427 
428  if (checked_data == CheckResult::not_a_number) {
429  THROW( ExcUndefElementValue() << EI_Field(field_name_) );
430  }
431 
432  if (is_native || this->interpolation_==DataInterpolation::identic_msh || this->interpolation_==DataInterpolation::equivalent_msh) {
433  this->calculate_native_values(data_vec);
434  } else if (this->interpolation_==DataInterpolation::gauss_p0) {
435  this->interpolate_gauss(data_vec);
436  } else { // DataInterpolation::interp_p0
437  this->interpolate_intersection(data_vec);
438  }
439 
440  return true;
441  } else return false;
442 
443 }
444 
445 
446 template <int spacedim, class Value>
448 {
449  static const unsigned int quadrature_order = 4; // parameter of quadrature
450  std::shared_ptr<Mesh> source_mesh = ReaderCache::get_mesh(reader_file_);
451  std::vector<unsigned int> searched_elements; // stored suspect elements in calculating the intersection
452  std::vector<arma::vec::fixed<3>> q_points; // real coordinates of quadrature points
453  std::vector<double> q_weights; // weights of quadrature points
454  unsigned int quadrature_size; // size of quadrature point and weight vector
455  std::vector<double> sum_val(dh_->max_elem_dofs()); // sum of value of one quadrature point
456  unsigned int elem_count; // count of intersect (source) elements of one quadrature point
457  std::vector<double> elem_value(dh_->max_elem_dofs()); // computed value of one (target) element
458  bool contains; // sign if source element contains quadrature point
459 
460  {
461  // set size of vectors to maximal count of quadrature points
462  QGauss<3> quad(quadrature_order);
463  q_points.resize(quad.size());
464  q_weights.resize(quad.size());
465  }
466 
467  Mesh *mesh;
468  if (this->boundary_domain_) mesh = dh_->mesh()->get_bc_mesh();
469  else mesh = dh_->mesh();
470  for (auto ele : mesh->elements_range()) {
471  std::fill(elem_value.begin(), elem_value.end(), 0.0);
472  switch (ele->dim()) {
473  case 0:
474  quadrature_size = 1;
475  q_points[0] = ele.node(0)->point();
476  q_weights[0] = 1.0;
477  break;
478  case 1:
479  quadrature_size = value_handler1_.compute_quadrature(q_points, q_weights, ele, quadrature_order);
480  break;
481  case 2:
482  quadrature_size = value_handler2_.compute_quadrature(q_points, q_weights, ele, quadrature_order);
483  break;
484  case 3:
485  quadrature_size = value_handler3_.compute_quadrature(q_points, q_weights, ele, quadrature_order);
486  break;
487  }
488  searched_elements.clear();
489  source_mesh->get_bih_tree().find_bounding_box(ele.bounding_box(), searched_elements);
490 
491  for (unsigned int i=0; i<quadrature_size; ++i) {
492  std::fill(sum_val.begin(), sum_val.end(), 0.0);
493  elem_count = 0;
494  for (std::vector<unsigned int>::iterator it = searched_elements.begin(); it!=searched_elements.end(); it++) {
495  ElementAccessor<3> elm = source_mesh->element_accessor(*it);
496  contains=false;
497  switch (elm->dim()) {
498  case 0:
499  contains = arma::norm(elm.node(0)->point()-q_points[i], 2) < 4*std::numeric_limits<double>::epsilon();
500  break;
501  case 1:
502  contains = value_handler1_.get_mapping()->contains_point(q_points[i], elm);
503  break;
504  case 2:
505  contains = value_handler2_.get_mapping()->contains_point(q_points[i], elm);
506  break;
507  case 3:
508  contains = value_handler3_.get_mapping()->contains_point(q_points[i], elm);
509  break;
510  default:
511  ASSERT(false).error("Invalid element dimension!");
512  }
513  if ( contains ) {
514  // projection point in element
515  unsigned int index = sum_val.size() * (*it);
516  for (unsigned int j=0; j < sum_val.size(); j++) {
517  sum_val[j] += (*data_vec)[index+j];
518  }
519  ++elem_count;
520  }
521  }
522 
523  if (elem_count > 0) {
524  for (unsigned int j=0; j < sum_val.size(); j++) {
525  elem_value[j] += (sum_val[j] / elem_count) * q_weights[i];
526  }
527  }
528  }
529 
531  else dh_->get_dof_indices( ele, dof_indices_);
532  for (unsigned int i=0; i < elem_value.size(); i++) {
534  (*data_vec_)[dof_indices_[i]] = elem_value[i] * this->unit_conversion_coefficient_;
535  }
536  }
537 }
538 
539 
540 template <int spacedim, class Value>
542 {
543  std::shared_ptr<Mesh> source_mesh = ReaderCache::get_mesh(reader_file_);
544  std::vector<unsigned int> searched_elements; // stored suspect elements in calculating the intersection
545  std::vector<double> value(dh_->max_elem_dofs());
546  double total_measure, measure;
547 
548  Mesh *mesh;
549  if (this->boundary_domain_) mesh = dh_->mesh()->get_bc_mesh();
550  else mesh = dh_->mesh();
551  for (auto elm : mesh->elements_range()) {
552  if (elm.dim() == 3) {
553  xprintf(Err, "Dimension of element in target mesh must be 0, 1 or 2! elm.idx() = %d\n", elm.idx());
554  }
555 
556  double epsilon = 4* numeric_limits<double>::epsilon() * elm.measure();
557 
558  // gets suspect elements
559  if (elm.dim() == 0) {
560  searched_elements.clear();
561  source_mesh->get_bih_tree().find_point(elm.node(0)->point(), searched_elements);
562  } else {
563  BoundingBox bb = elm.bounding_box();
564  searched_elements.clear();
565  source_mesh->get_bih_tree().find_bounding_box(bb, searched_elements);
566  }
567 
568  // set zero values of value object
569  std::fill(value.begin(), value.end(), 0.0);
570  total_measure=0.0;
571 
572  START_TIMER("compute_pressure");
573  ADD_CALLS(searched_elements.size());
574 
575 
576  MappingP1<3,3> mapping;
577 
578  for (std::vector<unsigned int>::iterator it = searched_elements.begin(); it!=searched_elements.end(); it++)
579  {
580  ElementAccessor<3> ele = source_mesh->element_accessor(*it);
581  if (ele->dim() == 3) {
582  // get intersection (set measure = 0 if intersection doesn't exist)
583  switch (elm.dim()) {
584  case 0: {
585  arma::vec::fixed<3> real_point = elm.node(0)->point();
586  arma::mat::fixed<3, 4> elm_map = mapping.element_map(ele);
587  arma::vec::fixed<4> unit_point = mapping.project_real_to_unit(real_point, elm_map);
588 
589  measure = (std::fabs(arma::sum( unit_point )-1) <= 1e-14
590  && arma::min( unit_point ) >= 0)
591  ? 1.0 : 0.0;
592  break;
593  }
594  case 1: {
596  ComputeIntersection<1,3> CI(elm, ele, source_mesh.get());
597  CI.init();
598  CI.compute(is);
599 
600  IntersectionLocal<1,3> ilc(is);
601  measure = ilc.compute_measure() * elm.measure();
602  break;
603  }
604  case 2: {
606  ComputeIntersection<2,3> CI(elm, ele, source_mesh.get());
607  CI.init();
608  CI.compute(is);
609 
610  IntersectionLocal<2,3> ilc(is);
611  measure = 2 * ilc.compute_measure() * elm.measure();
612  break;
613  }
614  }
615 
616  //adds values to value_ object if intersection exists
617  if (measure > epsilon) {
618  unsigned int index = value.size() * (*it);
619  std::vector<double> &vec = *( data_vec.get() );
620  for (unsigned int i=0; i < value.size(); i++) {
621  value[i] += vec[index+i] * measure;
622  }
623  total_measure += measure;
624  }
625  }
626  }
627 
628  // computes weighted average, store it to data vector
629  if (total_measure > epsilon) {
632  else dh_->get_dof_indices( elm, dof_indices_ );
633  for (unsigned int i=0; i < value.size(); i++) {
634  (*data_vector)[ dof_indices_[i] ] = value[i] / total_measure;
635  }
636  } else {
637  WarningOut().fmt("Processed element with idx {} is out of source mesh!\n", elm.idx());
638  }
639  END_TIMER("compute_pressure");
640 
641  }
642 }
643 
644 
645 template <int spacedim, class Value>
647 {
648  // Same algorithm as in output of Node_data. Possibly code reuse.
649  unsigned int dof_size, data_vec_i;
650  std::vector<unsigned int> count_vector(data_vec_->size(), 0);
653 
654  // iterate through elements, assembly global vector and count number of writes
655  Mesh *mesh;
656  if (this->boundary_domain_) mesh = dh_->mesh()->get_bc_mesh();
657  else mesh = dh_->mesh();
658  for (auto ele : mesh->elements_range()) { // remove special case for rank == 0 - necessary for correct output
659  if (this->boundary_domain_) dof_size = value_handler1_.get_dof_indices( ele, dof_indices_ );
660  else dof_size = dh_->get_dof_indices( ele, dof_indices_ );
661  data_vec_i = ele.idx() * dof_indices_.size();
662  for (unsigned int i=0; i<dof_size; ++i, ++data_vec_i) {
663  (*data_vector)[ dof_indices_[i] ] += (*data_cache)[data_vec_i];
664  ++count_vector[ dof_indices_[i] ];
665  }
666  }
667 
668  // iterate through cells, assembly global vector and count number of writes - prepared solution for further development
669  /*for (auto cell : dh_->own_range()) {
670  dof_size = cell.get_dof_indices(dof_indices_);
671  data_vec_i = cell.elm_idx() * dof_indices_.size();
672  for (unsigned int i=0; i<dof_size; ++i, ++data_vec_i) {
673  (*data_vector)[ dof_indices_[i] ] += (*data_cache)[data_vec_i];
674  ++count_vector[ dof_indices_[i] ];
675  }
676  }*/
677 
678  // compute averages of values
679  for (unsigned int i=0; i<data_vec_->size(); ++i) {
680  if (count_vector[i]>0) (*data_vector)[i] /= count_vector[i];
681  }
682 }
683 
684 
685 template <int spacedim, class Value>
687  ASSERT_EQ(output_data_cache.n_values() * output_data_cache.n_elem(), dh_->n_global_dofs()).error();
688  ASSERT_EQ(output_data_cache.n_elem(), dof_indices_.size()).error();
689  double loc_values[output_data_cache.n_elem()];
690  unsigned int i, dof_filled_size;
691 
693  for (auto ele : dh_->mesh()->elements_range()) {
694  dof_filled_size = dh_->get_dof_indices( ele, dof_indices_);
695  for (i=0; i<dof_filled_size; ++i) loc_values[i] = (*data_vec)[ dof_indices_[0] ];
696  for ( ; i<output_data_cache.n_elem(); ++i) loc_values[i] = numeric_limits<double>::signaling_NaN();
697  output_data_cache.store_value( ele.idx(), loc_values );
698  }
699 
700  output_data_cache.set_dof_handler_hash( dh_->hash() );
701 }
702 
703 
704 
705 template <int spacedim, class Value>
706 inline unsigned int FieldFE<spacedim, Value>::data_size() const {
707  return data_vec_->size();
708 }
709 
710 
711 
712 template <int spacedim, class Value>
714 {}
715 
716 
717 // Instantiations of FieldFE
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: long_idx.hh:22
std::shared_ptr< DOFHandlerMultiDim > sequential()
Returns sequential version of the current dof handler.
Definition: dofhandler.cc:88
virtual ~FieldFE()
Destructor.
Definition: field_fe.cc:713
void init_unit_conversion_coefficient(const Input::Record &rec, const struct FieldAlgoBaseInitData &init_data)
Init value of unit_conversion_coefficient_ from input.
Bounding box in 3d ambient space.
Definition: bounding_box.hh:54
virtual void value_list(const std::vector< Point > &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list)
Definition: field_fe.cc:217
void set_boundary_dofs_vector(std::shared_ptr< std::vector< LongIdx > > boundary_dofs)
TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
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:447
FEValueHandler< 0, spacedim, Value > value_handler0_
Value handler that allows get value of 0D elements.
Definition: field_fe.hh:195
std::shared_ptr< std::vector< T > > ComponentDataPtr
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.
VectorMPI * data_vec
Store data of Field.
unsigned int component_idx_
Specify if the field is part of a MultiField and which component it is.
static void get_element_ids(const FilePath &file_path, const Mesh &mesh)
Definition: reader_cache.cc:80
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
void initialize(FEValueInitData init_data, MappingP1< elemdim, 3 > *map=nullptr)
Initialize data members.
std::string field_name_
field name read from input
Definition: field_fe.hh:220
double compute_measure() const override
Computes the relative measure of intersection object.
void zero_entries()
Definition: vector_mpi.hh:102
void fill_data_to_cache(ElementDataCache< double > &output_data_cache)
Definition: field_fe.cc:686
Abstract linear system class.
Definition: balance.hh:35
FieldFlag::Flags flags_
static VectorMPI * sequential(unsigned int size)
Definition: vector_mpi.hh:60
FEValueHandler< 2, spacedim, Value > value_handler2_
Value handler that allows get value of 2D elements.
Definition: field_fe.hh:199
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
void distribute_dofs(std::shared_ptr< DiscreteSpace > ds)
Distributes degrees of freedom on the mesh needed for the given discrete space.
Definition: dofhandler.cc:263
FilePath reader_file_
mesh reader file
Definition: field_fe.hh:217
void set_fe_data(std::shared_ptr< DOFHandlerMultiDim > dh, MappingP1< 1, 3 > *map1, MappingP1< 2, 3 > *map2, MappingP1< 3, 3 > *map3, VectorMPI *data)
Definition: field_fe.cc:137
#define INSTANCE_ALL(field)
void init()
Initializes lower dimensional objects. Sets correctly the pointers to Plucker coordinates and product...
Definition: mesh.h:80
Iterator< Ret > find(const string &key) const
void initialize(FEValueInitData init_data)
Initialize data members.
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:91
void set_mesh(const Mesh *mesh, bool boundary_domain) override
Definition: field_fe.cc:266
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:346
#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:232
FEValueHandler< 3, spacedim, Value > value_handler3_
Value handler that allows get value of 3D elements.
Definition: field_fe.hh:201
bool set_time(const TimeStep &time) override
Definition: field_fe.cc:397
FiniteElement< 3 > * fe3_
Same as previous, but represents 3D element.
Definition: field_fe.hh:214
static const Input::Type::Selection & get_interp_selection_input_type()
Definition: field_fe.cc:98
Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Returns one value in one given point.
std::shared_ptr< std::vector< LongIdx > > boundary_dofs_
Definition: field_fe.hh:245
Value::return_type r_value_
void fill_boundary_dofs()
Definition: field_fe.cc:307
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:303
Symmetric Gauss-Legendre quadrature formulae on simplices.
unsigned int data_size() const
Definition: field_fe.cc:706
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:195
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
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:238
void set_native_dh(std::shared_ptr< DOFHandlerMultiDim > dh) override
Definition: field_fe.cc:151
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:223
std::shared_ptr< DOFHandlerMultiDim > dh_
DOF handler object.
Definition: field_fe.hh:188
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:541
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())
void value_list(const std::vector< Point > &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.
DataInterpolation interpolation_
Specify type of FE data interpolation.
Definition: field_fe.hh:226
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:99
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:92
void resize(unsigned int local_size)
Definition: vector_mpi.hh:69
#define START_TIMER(tag)
Starts a timer with specified tag.
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Definition: dofhandler.hh:155
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:541
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
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Definition: field_fe.cc:193
virtual Range< ElementAccessor< 3 > > elements_range() const
Returns range of bulk elements.
Definition: mesh.cc:1033
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 reinit_fe_data(MappingP1< 1, 3 > *map1, MappingP1< 2, 3 > *map2, MappingP1< 3, 3 > *map3)
Ensure data setting of methods set_fe_data and set_native_dh.
Definition: field_fe.cc:160
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 get_dof_indices(const ElementAccessor< 3 > &cell, std::vector< LongIdx > &indices) const
TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
double read_time(Input::Iterator< Input::Tuple > time_it, double default_time=std::numeric_limits< double >::quiet_NaN()) const
FEValueHandler< 1, spacedim, Value > value_handler1_
Value handler that allows get value of 1D elements.
Definition: field_fe.hh:197
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:215
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
MappingP1< elemdim, 3 > * get_mapping()
Return mapping object.
FiniteElement< 1 > * fe1_
Same as previous, but represents 1D element.
Definition: field_fe.hh:210
static const Input::Type::Selection & get_disc_selection_input_type()
Definition: field_fe.cc:86
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
FiniteElement< 2 > * fe2_
Same as previous, but represents 2D element.
Definition: field_fe.hh:212
bool is_constant_in_space_
Flag detects that field is only dependent on time.
Value value_
Last value, prevents passing large values (vectors) by value.
VectorMPI * data_vec_
Store data of Field.
Definition: field_fe.hh:190
Input::Record in_rec_
Accessor to Input::Record.
Definition: field_fe.hh:235
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
Definition: asserts.hh:335
const Selection & close() const
Close the Selection, no more values can be added.
FieldFE(unsigned int n_comp=0)
Definition: field_fe.cc:127
const unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:139
void value_list(const std::vector< Point > &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.
std::vector< LongIdx > dof_indices_
Array of indexes to data_vec_, used for get/set values.
Definition: field_fe.hh:192
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:246
void set_boundary_dofs_vector(std::shared_ptr< std::vector< LongIdx > > 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.
FiniteElement< 0 > * fe0_
Definition: field_fe.hh:208
DataInterpolation
Definition: field_fe.hh:59
Definition: system.hh:64
virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData &init_data)
Definition: field_fe.cc:243
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
Accessor to the data with type Type::Tuple.
Definition: accessors.hh:412
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:87
bool contains_point(arma::vec point, ElementAccessor< 3 > elm)
Test if element contains given point.
Definition: mapping_p1.cc:302
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:44
unsigned int n_comp
number of components
void calculate_native_values(ElementDataCache< double >::ComponentDataPtr data_cache)
Calculate native data over all elements of target mesh.
Definition: field_fe.cc:646
FieldFlag::Flags flags_
Field flags.
Definition: field_fe.hh:229
arma::vec3 & point()
Definition: nodes.hh:67
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:334
Implementation of range helper class.
#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:327
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
unsigned int size()
Return size of output data.
Definition: vector_mpi.hh:136
Definition: field.hh:56
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:299
const Node * node(unsigned int ni) const
Definition: accessors.hh:145