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