Flow123d  JS_before_hm-1623-gd361259
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>
61 const Input::Type::Record & FieldFE<spacedim, Value>::get_input_type()
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  field_name_(""), fe_values_(4)
136 {
137  this->is_constant_in_space_ = false;
138 }
139 
140 
141 template <int spacedim, class Value>
142 VectorMPI FieldFE<spacedim, Value>::set_fe_data(std::shared_ptr<DOFHandlerMultiDim> dh, VectorMPI dof_values, unsigned int block_index)
143 {
144  dh_ = dh;
145  if (dof_values.size()==0) { //create data vector according to dof handler - Warning not tested yet
146  data_vec_ = dh_->create_vector();
148  } else {
149  data_vec_ = dof_values;
150  }
151 
152  if ( block_index == FieldFE<spacedim, Value>::undef_uint ) {
153  this->fill_fe_item<0>();
154  this->fill_fe_item<1>();
155  this->fill_fe_item<2>();
156  this->fill_fe_item<3>();
157  this->fe_ = dh_->ds()->fe();
158  } else {
159  this->fill_fe_system_data<0>(block_index);
160  this->fill_fe_system_data<1>(block_index);
161  this->fill_fe_system_data<2>(block_index);
162  this->fill_fe_system_data<3>(block_index);
164  std::dynamic_pointer_cast<FESystem<0>>( dh_->ds()->fe()[Dim<0>{}] )->fe()[block_index],
165  std::dynamic_pointer_cast<FESystem<1>>( dh_->ds()->fe()[Dim<1>{}] )->fe()[block_index],
166  std::dynamic_pointer_cast<FESystem<2>>( dh_->ds()->fe()[Dim<2>{}] )->fe()[block_index],
167  std::dynamic_pointer_cast<FESystem<3>>( dh_->ds()->fe()[Dim<3>{}] )->fe()[block_index]
168  );
169  }
170 
171  unsigned int ndofs = dh_->max_elem_dofs();
172 
173  // initialization data of value handlers
174  FEValueInitData init_data;
175  init_data.dh = dh_;
176  init_data.data_vec = data_vec_;
177  init_data.ndofs = ndofs;
178  init_data.n_comp = this->n_comp();
179  init_data.mixed_fe = this->fe_;
180 
181  // initialize value handler objects
182  init_data.range_begin = this->fe_item_[0].range_begin_;
183  init_data.range_end = this->fe_item_[0].range_end_;
184  value_handler0_.initialize(init_data);
185  init_data.range_begin = this->fe_item_[1].range_begin_;
186  init_data.range_end = this->fe_item_[1].range_end_;
187  value_handler1_.initialize(init_data);
188  init_data.range_begin = this->fe_item_[2].range_begin_;
189  init_data.range_end = this->fe_item_[2].range_end_;
190  value_handler2_.initialize(init_data);
191  init_data.range_begin = this->fe_item_[3].range_begin_;
192  init_data.range_end = this->fe_item_[3].range_end_;
193  value_handler3_.initialize(init_data);
194 
195  // set discretization
196  discretization_ = OutputTime::DiscreteSpace::UNDEFINED;
197  interpolation_ = DataInterpolation::equivalent_msh;
198 
199  return data_vec_;
200 }
201 
202 
203 /**
204  * Returns one value in one given point. ResultType can be used to avoid some costly calculation if the result is trivial.
205  */
206 template <int spacedim, class Value>
207 typename Value::return_type const & FieldFE<spacedim, Value>::value(const Point &p, const ElementAccessor<spacedim> &elm)
208 {
209  switch (elm.dim()) {
210  case 0:
211  return value_handler0_.value(p, elm);
212  case 1:
213  return value_handler1_.value(p, elm);
214  case 2:
215  return value_handler2_.value(p, elm);
216  case 3:
217  return value_handler3_.value(p, elm);
218  default:
219  ASSERT(false).error("Invalid element dimension!");
220  }
221 
222  return this->r_value_;
223 }
224 
225 
226 
227 /**
228  * Returns std::vector of scalar values in several points at once.
229  */
230 template <int spacedim, class Value>
233 {
234  ASSERT_EQ( point_list.size(), value_list.size() ).error();
235  ASSERT_DBG( point_list.n_rows() == spacedim && point_list.n_cols() == 1).error("Invalid point size.\n");
236 
237  switch (elm.dim()) {
238  case 0:
239  value_handler0_.value_list(point_list, elm, value_list);
240  break;
241  case 1:
242  value_handler1_.value_list(point_list, elm, value_list);
243  break;
244  case 2:
245  value_handler2_.value_list(point_list, elm, value_list);
246  break;
247  case 3:
248  value_handler3_.value_list(point_list, elm, value_list);
249  break;
250  default:
251  ASSERT(false).error("Invalid element dimension!");
252  }
253 }
254 
255 
256 
257 template <int spacedim, class Value>
259  ElementCacheMap &cache_map, unsigned int region_patch_idx)
260 {
261  ASSERT( !boundary_dofs_ ).error("boundary field NOT supported!!\n");
263 
264  unsigned int reg_chunk_begin = cache_map.region_chunk_begin(region_patch_idx);
265  unsigned int reg_chunk_end = cache_map.region_chunk_end(region_patch_idx);
266  unsigned int last_element_idx = -1;
267  DHCellAccessor cell = *( dh_->local_range().begin() ); //needs set variable for correct compiling
268  LocDofVec loc_dofs;
269  unsigned int range_bgn=0, range_end=0;
270 
271  for (unsigned int i_data = reg_chunk_begin; i_data < reg_chunk_end; ++i_data) { // i_eval_point_data
272  unsigned int elm_idx = cache_map.eval_point_data(i_data).i_element_;
273  if (elm_idx != last_element_idx) {
274  ElementAccessor<spacedim> elm(dh_->mesh(), elm_idx);
275  fe_values_[elm.dim()].reinit( elm );
276  cell = dh_->cell_accessor_from_element( elm_idx );
277  loc_dofs = cell.get_loc_dof_indices();
278  last_element_idx = elm_idx;
279  range_bgn = this->fe_item_[elm.dim()].range_begin_;
280  range_end = this->fe_item_[elm.dim()].range_end_;
281  }
282 
283  unsigned int i_ep=cache_map.eval_point_data(i_data).i_eval_point_;
284  //DHCellAccessor cache_cell = cache_map(cell);
285  mat_value.fill(0.0);
286  for (unsigned int i_dof=range_bgn, i_cdof=0; i_dof<range_end; i_dof++, i_cdof++) {
287  mat_value += data_vec_.get(loc_dofs[i_dof]) * this->handle_fe_shape(cell.dim(), i_cdof, i_ep);
288  }
289  data_cache.set(i_data) = mat_value;
290  }
291 }
292 
293 
294 template <int spacedim, class Value>
296 {
297  std::shared_ptr<EvalPoints> eval_points = cache_map.eval_points();
298  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)};
299  fe_values_[0].initialize(quads[0], *this->fe_[0_d], update_values);
300  fe_values_[1].initialize(quads[1], *this->fe_[1_d], update_values);
301  fe_values_[2].initialize(quads[2], *this->fe_[2_d], update_values);
302  fe_values_[3].initialize(quads[3], *this->fe_[3_d], update_values);
303 }
304 
305 
306 template <int spacedim, class Value>
307 template <unsigned int dim>
308 Quadrature FieldFE<spacedim, Value>::init_quad(std::shared_ptr<EvalPoints> eval_points)
309 {
310  Quadrature quad(dim, eval_points->size(dim));
311  for (unsigned int k=0; k<eval_points->size(dim); k++)
312  quad.set(k) = eval_points->local_point<dim>(k);
313  return quad;
314 }
315 
316 
317 template <int spacedim, class Value>
319  this->init_unit_conversion_coefficient(rec, init_data);
320  this->in_rec_ = rec;
321  flags_ = init_data.flags_;
322 
323 
324  // read data from input record
325  reader_file_ = FilePath( rec.val<FilePath>("mesh_data_file") );
326  field_name_ = rec.val<std::string>("field_name");
327  if (! rec.opt_val<OutputTime::DiscreteSpace>("input_discretization", discretization_) ) {
328  discretization_ = OutputTime::DiscreteSpace::UNDEFINED;
329  }
330  if (! rec.opt_val<DataInterpolation>("interpolation", interpolation_) ) {
331  interpolation_ = DataInterpolation::equivalent_msh;
332  }
333  if (! rec.opt_val("default_value", default_value_) ) {
334  default_value_ = numeric_limits<double>::signaling_NaN();
335  }
336 }
337 
338 
339 
340 template <int spacedim, class Value>
341 void FieldFE<spacedim, Value>::set_mesh(const Mesh *mesh, bool boundary_domain) {
342  // Mesh can be set only for field initialized from input.
344  ASSERT(field_name_ != "").error("Uninitialized FieldFE, did you call init_from_input()?\n");
345  this->boundary_domain_ = boundary_domain;
346  if (this->interpolation_ == DataInterpolation::identic_msh) {
348  } else {
349  auto source_mesh = ReaderCache::get_mesh(reader_file_);
350  ReaderCache::get_element_ids(reader_file_, *(source_mesh.get()));
351  if (this->interpolation_ == DataInterpolation::equivalent_msh) {
353  if (source_target_mesh_elm_map_->size() == 0) { // incompatible meshes
354  this->interpolation_ = DataInterpolation::gauss_p0;
355  WarningOut().fmt("Source mesh of FieldFE '{}' is not compatible with target mesh.\nInterpolation of input data will be changed to 'P0_gauss'.\n",
356  field_name_);
357  }
358  } else if (this->interpolation_ == DataInterpolation::interp_p0) {
359  if (!boundary_domain) {
360  this->interpolation_ = DataInterpolation::gauss_p0;
361  WarningOut().fmt("Interpolation 'P0_intersection' of FieldFE '{}' can't be used on bulk region.\nIt will be changed to 'P0_gauss'.\n",
362  field_name_);
363  }
364  }
365  }
366  this->make_dof_handler(mesh);
367  }
368 }
369 
370 
371 
372 template <int spacedim, class Value>
374  ASSERT(this->boundary_domain_);
375 
376  auto bc_mesh = dh_->mesh()->get_bc_mesh();
377  unsigned int n_comp = this->value_.n_rows() * this->value_.n_cols();
378  boundary_dofs_ = std::make_shared< std::vector<IntIdx> >( n_comp * bc_mesh->n_elements() );
379  std::vector<IntIdx> &in_vec = *( boundary_dofs_.get() );
380  unsigned int j = 0; // actual index to boundary_dofs_ vector
381 
382  for (auto ele : bc_mesh->elements_range()) {
383  IntIdx elm_shift = n_comp * ele.idx();
384  for (unsigned int i=0; i<n_comp; ++i, ++j) {
385  in_vec[j] = elm_shift + i;
386  }
387  }
388 
393 
395 }
396 
397 
398 
399 template <int spacedim, class Value>
401 
402  // temporary solution - these objects will be set through FieldCommon
404  switch (this->value_.n_rows() * this->value_.n_cols()) { // by number of components
405  case 1: { // scalar
406  fe = MixedPtr<FE_P_disc>(0);
407  break;
408  }
409  case 3: { // vector
410  MixedPtr<FE_P_disc> fe_base(0) ;
411  fe = mixed_fe_system(fe_base, FEType::FEVector, 3);
412  break;
413  }
414  case 9: { // tensor
415  MixedPtr<FE_P_disc> fe_base(0) ;
416  fe = mixed_fe_system(fe_base, FEType::FETensor, 9);
417  break;
418  }
419  default:
420  ASSERT(false).error("Should not happen!\n");
421  }
422 
423  std::shared_ptr<DOFHandlerMultiDim> dh_par = std::make_shared<DOFHandlerMultiDim>( const_cast<Mesh &>(*mesh) );
424  std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( &const_cast<Mesh &>(*mesh), fe);
425  dh_par->distribute_dofs(ds);
426  dh_ = dh_par;
427  unsigned int ndofs = dh_->max_elem_dofs();
428 
429  this->fill_fe_item<0>();
430  this->fill_fe_item<1>();
431  this->fill_fe_item<2>();
432  this->fill_fe_item<3>();
433  this->fe_ = dh_->ds()->fe();
434 
435  if (this->boundary_domain_) fill_boundary_dofs(); // temporary solution for boundary mesh
436  else data_vec_ = VectorMPI::sequential( dh_->lsize() ); // allocate data_vec_
437 
438  // initialization data of value handlers
439  FEValueInitData init_data;
440  init_data.dh = dh_;
441  init_data.data_vec = data_vec_;
442  init_data.ndofs = ndofs;
443  init_data.n_comp = this->n_comp();
444  init_data.mixed_fe = this->fe_;
445 
446  // initialize value handler objects
447  init_data.range_begin = this->fe_item_[0].range_begin_;
448  init_data.range_end = this->fe_item_[0].range_end_;
449  value_handler0_.initialize(init_data);
450  init_data.range_begin = this->fe_item_[1].range_begin_;
451  init_data.range_end = this->fe_item_[1].range_end_;
452  value_handler1_.initialize(init_data);
453  init_data.range_begin = this->fe_item_[2].range_begin_;
454  init_data.range_end = this->fe_item_[2].range_end_;
455  value_handler2_.initialize(init_data);
456  init_data.range_begin = this->fe_item_[3].range_begin_;
457  init_data.range_end = this->fe_item_[3].range_end_;
458  value_handler3_.initialize(init_data);
459 }
460 
461 
462 
463 template <int spacedim, class Value>
465  // Time can be set only for field initialized from input.
467  ASSERT(field_name_ != "").error("Uninitialized FieldFE, did you call init_from_input()?\n");
468  ASSERT_PTR(dh_)(field_name_).error("Null target mesh pointer of finite element field, did you call set_mesh()?\n");
469  if ( reader_file_ == FilePath() ) return false;
470 
471  unsigned int n_components = this->value_.n_rows() * this->value_.n_cols();
472  double time_unit_coef = time.read_coef(in_rec_.find<Input::Record>("time_unit"));
473  double time_shift = time.read_time( in_rec_.find<Input::Tuple>("read_time_shift") );
474  double read_time = (time.end()+time_shift) / time_unit_coef;
475  BaseMeshReader::HeaderQuery header_query(field_name_, read_time, this->discretization_, dh_->hash());
476  ReaderCache::get_reader(reader_file_)->find_header(header_query);
477  // TODO: use default and check NaN values in data_vec
478 
479  unsigned int n_entities;
480  bool is_native = (header_query.discretization == OutputTime::DiscreteSpace::NATIVE_DATA);
481  bool boundary;
482  if (is_native || this->interpolation_==DataInterpolation::identic_msh || this->interpolation_==DataInterpolation::equivalent_msh) {
483  boundary = this->boundary_domain_;
484  } else {
485  boundary = false;
486  }
487  if (this->interpolation_==DataInterpolation::identic_msh) {
488  n_entities = boundary ? dh_->mesh()->get_bc_mesh()->n_elements() : dh_->mesh()->n_elements();
489  } else {
490  n_entities = boundary ? ReaderCache::get_mesh(reader_file_)->get_bc_mesh()->n_elements() : ReaderCache::get_mesh(reader_file_)->n_elements();
491  }
492  auto input_data_cache = ReaderCache::get_reader(reader_file_)->template get_element_data<double>(n_entities, n_components,
493  boundary, this->component_idx_);
494  CheckResult checked_data = ReaderCache::get_reader(reader_file_)->scale_and_check_limits(field_name_,
496 
497 
498  if (checked_data == CheckResult::not_a_number) {
499  THROW( ExcUndefElementValue() << EI_Field(field_name_) );
500  }
501 
502  if (is_native) {
503  this->calculate_native_values(input_data_cache);
504  } else if (this->interpolation_==DataInterpolation::identic_msh) {
505  this->calculate_identic_values(input_data_cache);
506  } else if (this->interpolation_==DataInterpolation::equivalent_msh) {
507  this->calculate_equivalent_values(input_data_cache);
508  } else if (this->interpolation_==DataInterpolation::gauss_p0) {
509  this->interpolate_gauss(input_data_cache);
510  } else { // DataInterpolation::interp_p0
511  this->interpolate_intersection(input_data_cache);
512  }
513 
514  return true;
515  } else return false;
516 
517 }
518 
519 
520 template <int spacedim, class Value>
522 {
523  static const unsigned int quadrature_order = 4; // parameter of quadrature
524  std::shared_ptr<Mesh> source_mesh = ReaderCache::get_mesh(reader_file_);
525  std::vector<unsigned int> searched_elements; // stored suspect elements in calculating the intersection
526  std::vector<arma::vec::fixed<3>> q_points; // real coordinates of quadrature points
527  std::vector<double> q_weights; // weights of quadrature points
528  unsigned int quadrature_size=0; // size of quadrature point and weight vector
529  std::vector<double> sum_val(dh_->max_elem_dofs()); // sum of value of one quadrature point
530  unsigned int elem_count; // count of intersect (source) elements of one quadrature point
531  std::vector<double> elem_value(dh_->max_elem_dofs()); // computed value of one (target) element
532  bool contains; // sign if source element contains quadrature point
533 
534  {
535  // set size of vectors to maximal count of quadrature points
536  QGauss quad(3, quadrature_order);
537  q_points.resize(quad.size());
538  q_weights.resize(quad.size());
539  }
540 
541  for (auto cell : dh_->own_range()) {
542  auto ele = cell.elm();
543  std::fill(elem_value.begin(), elem_value.end(), 0.0);
544  switch (cell.dim()) {
545  case 0:
546  quadrature_size = 1;
547  q_points[0] = *ele.node(0);
548  q_weights[0] = 1.0;
549  break;
550  case 1:
551  quadrature_size = value_handler1_.compute_quadrature(q_points, q_weights, ele, quadrature_order);
552  break;
553  case 2:
554  quadrature_size = value_handler2_.compute_quadrature(q_points, q_weights, ele, quadrature_order);
555  break;
556  case 3:
557  quadrature_size = value_handler3_.compute_quadrature(q_points, q_weights, ele, quadrature_order);
558  break;
559  }
560  searched_elements.clear();
561  source_mesh->get_bih_tree().find_bounding_box(ele.bounding_box(), searched_elements);
562 
563  for (unsigned int i=0; i<quadrature_size; ++i) {
564  std::fill(sum_val.begin(), sum_val.end(), 0.0);
565  elem_count = 0;
566  for (std::vector<unsigned int>::iterator it = searched_elements.begin(); it!=searched_elements.end(); it++) {
567  ElementAccessor<3> elm = source_mesh->element_accessor(*it);
568  contains=false;
569  switch (elm->dim()) {
570  case 0:
571  contains = arma::norm(*elm.node(0) - q_points[i], 2) < 4*std::numeric_limits<double>::epsilon();
572  break;
573  case 1:
574  contains = MappingP1<1,3>::contains_point(q_points[i], elm);
575  break;
576  case 2:
577  contains = MappingP1<2,3>::contains_point(q_points[i], elm);
578  break;
579  case 3:
580  contains = MappingP1<3,3>::contains_point(q_points[i], elm);
581  break;
582  default:
583  ASSERT(false).error("Invalid element dimension!");
584  }
585  if ( contains ) {
586  // projection point in element
587  unsigned int index = sum_val.size() * (*it);
588  for (unsigned int j=0; j < sum_val.size(); j++) {
589  sum_val[j] += (*data_vec)[index+j];
590  }
591  ++elem_count;
592  }
593  }
594 
595  if (elem_count > 0) {
596  for (unsigned int j=0; j < sum_val.size(); j++) {
597  elem_value[j] += (sum_val[j] / elem_count) * q_weights[i];
598  }
599  }
600  }
601 
602  LocDofVec loc_dofs;
603  if (this->boundary_domain_) loc_dofs = value_handler1_.get_loc_dof_indices(cell.elm_idx());
604  else loc_dofs = cell.get_loc_dof_indices();
605 
606  ASSERT_LE_DBG(loc_dofs.n_elem, elem_value.size());
607  for (unsigned int i=0; i < elem_value.size(); i++) {
608  ASSERT_LT_DBG( loc_dofs[i], (int)data_vec_.size());
609  data_vec_.set( loc_dofs[i], elem_value[i] * this->unit_conversion_coefficient_ );
610  }
611  }
612 }
613 
614 
615 template <int spacedim, class Value>
617 {
618  std::shared_ptr<Mesh> source_mesh = ReaderCache::get_mesh(reader_file_);
619  std::vector<unsigned int> searched_elements; // stored suspect elements in calculating the intersection
620  std::vector<double> value(dh_->max_elem_dofs());
621  double total_measure;
622  double measure = 0;
623 
624  Mesh *mesh;
625  if (this->boundary_domain_) mesh = dh_->mesh()->get_bc_mesh();
626  else mesh = dh_->mesh();
627  for (auto elm : mesh->elements_range()) {
628  if (elm.dim() == 3) {
629  THROW( ExcInvalidElemeDim() << EI_ElemIdx(elm.idx()) );
630  }
631 
632  double epsilon = 4* numeric_limits<double>::epsilon() * elm.measure();
633 
634  // gets suspect elements
635  if (elm.dim() == 0) {
636  searched_elements.clear();
637  source_mesh->get_bih_tree().find_point(*elm.node(0), searched_elements);
638  } else {
639  BoundingBox bb = elm.bounding_box();
640  searched_elements.clear();
641  source_mesh->get_bih_tree().find_bounding_box(bb, searched_elements);
642  }
643 
644  // set zero values of value object
645  std::fill(value.begin(), value.end(), 0.0);
646  total_measure=0.0;
647 
648  START_TIMER("compute_pressure");
649  ADD_CALLS(searched_elements.size());
650 
651 
652  for (std::vector<unsigned int>::iterator it = searched_elements.begin(); it!=searched_elements.end(); it++)
653  {
654  ElementAccessor<3> ele = source_mesh->element_accessor(*it);
655  if (ele->dim() == 3) {
656  // get intersection (set measure = 0 if intersection doesn't exist)
657  switch (elm.dim()) {
658  case 0: {
659  arma::vec::fixed<3> real_point = *elm.node(0);
660  arma::mat::fixed<3, 4> elm_map = MappingP1<3,3>::element_map(ele);
661  arma::vec::fixed<4> unit_point = MappingP1<3,3>::project_real_to_unit(real_point, elm_map);
662 
663  measure = (std::fabs(arma::sum( unit_point )-1) <= 1e-14
664  && arma::min( unit_point ) >= 0)
665  ? 1.0 : 0.0;
666  break;
667  }
668  case 1: {
670  ComputeIntersection<1,3> CI(elm, ele, source_mesh.get());
671  CI.init();
672  CI.compute(is);
673 
674  IntersectionLocal<1,3> ilc(is);
675  measure = ilc.compute_measure() * elm.measure();
676  break;
677  }
678  case 2: {
680  ComputeIntersection<2,3> CI(elm, ele, source_mesh.get());
681  CI.init();
682  CI.compute(is);
683 
684  IntersectionLocal<2,3> ilc(is);
685  measure = 2 * ilc.compute_measure() * elm.measure();
686  break;
687  }
688  }
689 
690  //adds values to value_ object if intersection exists
691  if (measure > epsilon) {
692  unsigned int index = value.size() * (*it);
693  std::vector<double> &vec = *( data_vec.get() );
694  for (unsigned int i=0; i < value.size(); i++) {
695  value[i] += vec[index+i] * measure;
696  }
697  total_measure += measure;
698  }
699  }
700  }
701 
702  // computes weighted average, store it to data vector
703  if (total_measure > epsilon) {
704  LocDofVec loc_dofs;
705  if (this->boundary_domain_) loc_dofs = value_handler1_.get_loc_dof_indices(elm.idx());
706  else{
707  DHCellAccessor cell = dh_->cell_accessor_from_element(elm.idx());
708  loc_dofs = cell.get_loc_dof_indices();
709  }
710 
711  ASSERT_LE_DBG(loc_dofs.n_elem, value.size());
712  for (unsigned int i=0; i < value.size(); i++) {
713  data_vec_.set(loc_dofs[i], value[i] / total_measure);
714  }
715  } else {
716  WarningOut().fmt("Processed element with idx {} is out of source mesh!\n", elm.idx());
717  }
718  END_TIMER("compute_pressure");
719 
720  }
721 }
722 
723 
724 template <int spacedim, class Value>
726 {
727  // Same algorithm as in output of Node_data. Possibly code reuse.
728  unsigned int dof_size, data_vec_i;
729  std::vector<unsigned int> count_vector(data_vec_.size(), 0);
731  std::vector<LongIdx> global_dof_indices(dh_->max_elem_dofs());
732  std::vector<LongIdx> &source_target_vec = *(source_target_mesh_elm_map_.get());
733 
734  // iterate through cells, assembly MPIVector
735  for (auto cell : dh_->own_range()) {
736  dof_size = cell.get_dof_indices(global_dof_indices);
737  LocDofVec loc_dofs = cell.get_loc_dof_indices();
738  data_vec_i = source_target_vec[cell.elm_idx()] * dof_size;
739  ASSERT_EQ_DBG(dof_size, loc_dofs.n_elem);
740  for (unsigned int i=0; i<dof_size; ++i, ++data_vec_i) {
741  data_vec_.add( loc_dofs[i], (*data_cache)[ data_vec_i ] );
742  ++count_vector[ loc_dofs[i] ];
743  }
744  }
745 
746  // compute averages of values
747  for (unsigned int i=0; i<data_vec_.size(); ++i) {
748  if (count_vector[i]>0) data_vec_.normalize(i, count_vector[i]);
749  }
750 }
751 
752 
753 template <int spacedim, class Value>
755 {
756  // Same algorithm as in output of Node_data. Possibly code reuse.
757  unsigned int data_vec_i, i_elm;
758  std::vector<unsigned int> count_vector(data_vec_.size(), 0);
760 
761  if (this->boundary_domain_) {
762  // iterate through elements, assembly global vector and count number of writes
763  Mesh *mesh = dh_->mesh()->get_bc_mesh();
764  i_elm=0;
765  for (auto ele : mesh->elements_range()) {
766  LocDofVec loc_dofs = value_handler1_.get_loc_dof_indices(ele.idx());
767  data_vec_i = i_elm * dh_->max_elem_dofs();
768  for (unsigned int i=0; i<loc_dofs.n_elem; ++i, ++data_vec_i) {
769  ASSERT_LT_DBG(loc_dofs[i], (LongIdx)data_vec_.size());
770  data_vec_.add( loc_dofs[i], (*data_cache)[data_vec_i] );
771  ++count_vector[ loc_dofs[i] ];
772  }
773  i_elm++;
774  }
775  }
776  else {
777  // iterate through cells, assembly global vector and count number of writes - prepared solution for further development
778  i_elm=0;
779  for (auto cell : dh_->own_range()) {
780  LocDofVec loc_dofs = cell.get_loc_dof_indices();
781  data_vec_i = i_elm * dh_->max_elem_dofs();
782  for (unsigned int i=0; i<loc_dofs.n_elem; ++i, ++data_vec_i) {
783  ASSERT_LT_DBG(loc_dofs[i], (LongIdx)data_vec_.size());
784  data_vec_.add( loc_dofs[i], (*data_cache)[data_vec_i] );
785  ++count_vector[ loc_dofs[i] ];
786  }
787  i_elm++;
788  }
789  }
790 
791  // compute averages of values
792  for (unsigned int i=0; i<data_vec_.size(); ++i) {
793  if (count_vector[i]>0) data_vec_.normalize(i, count_vector[i]);
794  }
795 }
796 
797 
798 template <int spacedim, class Value>
800 {
801  // Same algorithm as in output of Node_data. Possibly code reuse.
802  unsigned int data_vec_i;
803  std::vector<unsigned int> count_vector(data_vec_.size(), 0);
805  std::vector<LongIdx> &source_target_vec = *(source_target_mesh_elm_map_.get());
806 
807  // iterate through elements, assembly global vector and count number of writes
808  if (this->boundary_domain_) {
809  Mesh *mesh = dh_->mesh()->get_bc_mesh();
810  for (auto ele : mesh->elements_range()) {
811  LocDofVec loc_dofs = value_handler1_.get_loc_dof_indices(ele.idx());
812  if (source_target_vec[ele.mesh_idx()] == (int)(Mesh::undef_idx)) { // undefined value in input data mesh
813  if ( std::isnan(default_value_) )
814  THROW( ExcUndefElementValue() << EI_Field(field_name_) );
815  for (unsigned int i=0; i<loc_dofs.n_elem; ++i) {
816  ASSERT_LT_DBG(loc_dofs[i], (LongIdx)data_vec_.size());
818  ++count_vector[ loc_dofs[i] ];
819  }
820  } else {
821  data_vec_i = source_target_vec[ele.mesh_idx()] * dh_->max_elem_dofs();
822  for (unsigned int i=0; i<loc_dofs.n_elem; ++i, ++data_vec_i) {
823  ASSERT_LT_DBG(loc_dofs[i], (LongIdx)data_vec_.size());
824  data_vec_.add( loc_dofs[i], (*data_cache)[data_vec_i] );
825  ++count_vector[ loc_dofs[i] ];
826  }
827  }
828  }
829  }
830  else {
831  // iterate through cells, assembly global vector and count number of writes - prepared solution for further development
832  for (auto cell : dh_->own_range()) {
833  LocDofVec loc_dofs = cell.get_loc_dof_indices();
834  if (source_target_vec[cell.elm_idx()] == (int)(Mesh::undef_idx)) { // undefined value in input data mesh
835  if ( std::isnan(default_value_) )
836  THROW( ExcUndefElementValue() << EI_Field(field_name_) );
837  for (unsigned int i=0; i<loc_dofs.n_elem; ++i) {
838  ASSERT_LT_DBG(loc_dofs[i], (LongIdx)data_vec_.size());
840  ++count_vector[ loc_dofs[i] ];
841  }
842  } else {
843  data_vec_i = source_target_vec[cell.elm_idx()] * dh_->max_elem_dofs();
844  for (unsigned int i=0; i<loc_dofs.n_elem; ++i, ++data_vec_i) {
845  ASSERT_LT_DBG(loc_dofs[i], (LongIdx)data_vec_.size());
846  data_vec_.add( loc_dofs[i], (*data_cache)[data_vec_i] );
847  ++count_vector[ loc_dofs[i] ];
848  }
849  }
850  }
851  }
852 
853  // compute averages of values
854  for (unsigned int i=0; i<data_vec_.size(); ++i) {
855  if (count_vector[i]>0) data_vec_.normalize(i, count_vector[i]);
856  }
857 }
858 
859 
860 template <int spacedim, class Value>
862  ASSERT_EQ(output_data_cache.n_values() * output_data_cache.n_comp(), dh_->distr()->lsize()).error();
863  double loc_values[output_data_cache.n_comp()];
864  unsigned int i;
865 
866  for (auto dh_cell : dh_->own_range()) {
867  LocDofVec loc_dofs = dh_cell.get_loc_dof_indices();
868  for (i=0; i<loc_dofs.n_elem; ++i) loc_values[i] = data_vec_.get( loc_dofs[i] );
869  for ( ; i<output_data_cache.n_comp(); ++i) loc_values[i] = numeric_limits<double>::signaling_NaN();
870  output_data_cache.store_value( dh_cell.local_idx(), loc_values );
871  }
872 
873  output_data_cache.set_dof_handler_hash( dh_->hash() );
874 }
875 
876 
877 
878 template <int spacedim, class Value>
879 inline unsigned int FieldFE<spacedim, Value>::data_size() const {
880  return data_vec_.size();
881 }
882 
883 
884 
885 template <int spacedim, class Value>
888 }
889 
890 
891 
892 template <int spacedim, class Value>
895 }
896 
897 
898 
899 /*template <int spacedim, class Value>
900 Armor::ArmaMat<typename Value::element_type, Value::NRows_, Value::NCols_> FieldFE<spacedim, Value>::handle_fe_shape(unsigned int dim,
901  unsigned int i_dof, unsigned int i_qp, unsigned int comp_index)
902 {
903  Armor::ArmaMat<typename Value::element_type, Value::NCols_, Value::NRows_> v;
904  for (unsigned int c=0; c<Value::NRows_*Value::NCols_; ++c)
905  v(c/spacedim,c%spacedim) = fe_values_[dim].shape_value_component(i_dof, i_qp, comp_index+c);
906  if (Value::NRows_ == Value::NCols_)
907  return v;
908  else
909  return v.t();
910 }*/
911 
912 
913 
914 template <int spacedim, class Value>
916 {}
917 
918 
919 // Instantiations of FieldFE
unsigned int region_chunk_begin(unsigned int region_patch_idx) const
Return begin position of region chunk in FieldValueCache.
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:211
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
Shape function values.
Definition: update_flags.hh:87
virtual ~FieldFE()
Destructor.
Definition: field_fe.cc:915
void init_unit_conversion_coefficient(const Input::Record &rec, const struct FieldAlgoBaseInitData &init_data)
Init value of unit_conversion_coefficient_ from input.
unsigned int size() const
Return size of output data.
Definition: vector_mpi.hh:98
Definition: mixed.hh:25
Bounding box in 3d ambient space.
Definition: bounding_box.hh:54
void set_dof_handler_hash(std::size_t hash)
void interpolate_gauss(ElementDataCache< double >::ComponentDataPtr data_vec)
Interpolate data (use Gaussian distribution) over all elements of target mesh.
Definition: field_fe.cc:521
unsigned int size() const
Definition: armor.hh:728
FEValueHandler< 0, spacedim, Value > value_handler0_
Value handler that allows get value of 0D elements.
Definition: field_fe.hh:247
Armor::Array< double >::ArrayMatSet set(uint i)
Definition: quadrature.hh:98
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:332
std::shared_ptr< std::vector< T > > ComponentDataPtr
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
MixedPtr< FiniteElement > fe_
Definition: field_fe.hh:294
Classes with algorithms for computation of intersections of meshes.
void calculate_equivalent_values(ElementDataCache< double >::ComponentDataPtr data_cache)
Calculate data of equivalent_mesh interpolation on input over all elements of target mesh...
Definition: field_fe.cc:799
Some value(s) is set to NaN.
void initialize(FEValueInitData init_data)
Initialize data members.
std::vector< FEValues< spacedim > > fe_values_
List of FEValues objects of dimensions 0,1,2,3 used for value calculation.
Definition: field_fe.hh:287
unsigned int n_values() const
unsigned int component_idx_
Specify if the field is part of a MultiField and which component it is.
LocDofVec get_loc_dof_indices(unsigned int cell_idx) const
TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
static void get_element_ids(const FilePath &file_path, const Mesh &mesh)
Definition: reader_cache.cc:73
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
std::string field_name_
field name read from input
Definition: field_fe.hh:259
void local_to_ghost_begin()
local_to_ghost_{begin,end} updates the ghost values on neighbouring processors from local values ...
Definition: vector_mpi.hh:147
double compute_measure() const override
Computes the relative measure of intersection object.
int IntIdx
Definition: index_types.hh:25
void zero_entries()
Definition: vector_mpi.hh:82
Abstract linear system class.
Definition: balance.hh:40
unsigned int i_element_
mesh_idx of ElementAccessor appropriate to element
void value_list(const Armor::array &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list)
Returns std::vector of scalar values in several points at once.
FieldFlag::Flags flags_
static const unsigned int undef_idx
Definition: mesh.h:121
FEValueHandler< 2, spacedim, Value > value_handler2_
Value handler that allows get value of 2D elements.
Definition: field_fe.hh:251
const EvalPointData & eval_point_data(unsigned int point_idx) const
Return item of eval_point_data_ specified by its position.
unsigned int range_end
Holds end of component range of evaluation.
VectorMPI data_vec_
Store data of Field.
Definition: field_fe.hh:244
void store_value(unsigned int idx, const T *value)
Fundamental simplicial intersections.
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
static BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map)
Definition: mapping_p1.cc:69
unsigned int dim() const
Return dimension of element appropriate to cell.
FilePath reader_file_
mesh reader file
Definition: field_fe.hh:256
double read_coef(Input::Iterator< Input::Record > unit_it) const
#define INSTANCE_ALL(field)
void init()
Initializes lower dimensional objects. Sets correctly the pointers to Plucker coordinates and product...
Directing class of FieldValueCache.
Definition: mesh.h:77
Iterator< Ret > find(const string &key) const
void initialize(FEValueInitData init_data)
Initialize data members.
Cell accessor allow iterate over DOF handler cells.
Helper struct stores data for initizalize descentants of FieldAlgorithmBase.
Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Returns one value in one given point.
Class FESystem for compound finite elements.
uint n_cols() const
Definition: armor.hh:720
static const Input::Type::Record & get_input_type()
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
void set_mesh(const Mesh *mesh, bool boundary_domain) override
Definition: field_fe.cc:341
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:347
#define ADD_CALLS(n_calls)
Increase number of calls in actual timer.
double default_value_
Default value of element if not set in mesh data file.
Definition: field_fe.hh:271
FEValueHandler< 3, spacedim, Value > value_handler3_
Value handler that allows get value of 3D elements.
Definition: field_fe.hh:253
bool set_time(const TimeStep &time) override
Definition: field_fe.cc:464
static const Input::Type::Selection & get_interp_selection_input_type()
Definition: field_fe.cc:104
Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Returns one value in one given point.
Value::return_type r_value_
void fill_boundary_dofs()
Definition: field_fe.cc:373
NodeAccessor< 3 > node(unsigned int ni) const
Definition: accessors.hh:200
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
Symmetric Gauss-Legendre quadrature formulae on simplices.
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:142
unsigned int data_size() const
Definition: field_fe.cc:879
constexpr bool match(Mask mask) const
Definition: flag_array.hh:163
void init()
Initializes lower dimensional objects. Sets correctly the pointers to Plucker coordinates and product...
double unit_conversion_coefficient_
Coeficient of conversion of user-defined unit.
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:196
unsigned int dim() const
Definition: elements.h:121
static Default optional()
The factory function to make an empty default value which is optional.
Definition: type_record.hh:124
MixedPtr< FiniteElement > mixed_fe
FiniteElement objects of all dimensions.
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:277
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:262
std::shared_ptr< DOFHandlerMultiDim > dh_
DOF handler object.
Definition: field_fe.hh:242
Quadrature init_quad(std::shared_ptr< EvalPoints > eval_points)
Initialize FEValues object of given dimension.
Definition: field_fe.cc:308
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:534
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())
DummyInt isnan(...)
Definition: format.h:306
DataInterpolation interpolation_
Specify type of FE data interpolation.
Definition: field_fe.hh:265
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:101
std::array< FEItem, 4 > fe_item_
Holds specific data of field evaluation over all dimensions.
Definition: field_fe.hh:293
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:524
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
const Ret val(const string &key) const
void resize(unsigned int local_size)
Definition: vector_mpi.cc:50
#define START_TIMER(tag)
Starts a timer with specified tag.
Selection & add_value(const int value, const std::string &key, const std::string &description="", TypeBase::attribute_map attributes=TypeBase::attribute_map())
Adds one new value with name given by key to the Selection.
void interpolate_intersection(ElementDataCache< double >::ComponentDataPtr data_vec)
Interpolate data (use intersection library) over all elements of target mesh.
Definition: field_fe.cc:616
std::shared_ptr< std::vector< IntIdx > > boundary_dofs_
Definition: field_fe.hh:284
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
unsigned int region_chunk_end(unsigned int region_patch_idx) const
Return end position of region chunk in FieldValueCache.
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Definition: field_fe.cc:207
ArrayMatSet set(uint index)
Definition: armor.hh:838
uint n_rows() const
Definition: armor.hh:715
virtual Range< ElementAccessor< 3 > > elements_range() const
Returns range of bulk elements.
Definition: mesh.cc:1241
void local_to_ghost_data_scatter_begin()
Call begin scatter functions (local to ghost) on data vector.
Definition: field_fe.cc:886
double end() const
unsigned int ndofs
number of dofs
static Input::Type::Default get_input_default()
Space< spacedim >::Point Point
void add(unsigned int pos, double val)
Add value to item on given position.
Definition: vector_mpi.hh:123
void local_to_ghost_data_scatter_end()
Call end scatter functions (local to ghost) on data vector.
Definition: field_fe.cc:893
unsigned int i_eval_point_
index of point in EvalPoint object
void set_boundary_dofs_vector(std::shared_ptr< std::vector< IntIdx > > boundary_dofs)
TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
unsigned int compute_quadrature(std::vector< arma::vec::fixed< 3 >> &q_points, std::vector< double > &q_weights, const ElementAccessor< spacedim > &elm, unsigned int order=3)
Compute real coordinates and weights (use QGauss) for given element.
CheckResult
Return type of method that checked data stored in ElementDataCache (NaN values, limits) ...
std::shared_ptr< DOFHandlerMultiDim > dh
DOF handler object.
double read_time(Input::Iterator< Input::Tuple > time_it, double default_time=std::numeric_limits< double >::quiet_NaN()) const
typename arma::Mat< Type >::template fixed< nr, nc > ArmaMat
Definition: armor.hh:502
void local_to_ghost_end()
local_to_ghost_{begin,end} updates the ghost values on neighbouring processors from local values ...
Definition: vector_mpi.hh:151
#define ASSERT_LE_DBG(a, b)
Definition of comparative assert macro (Less or Equal) only for debug mode.
Definition: asserts.hh:308
FEValueHandler< 1, spacedim, Value > value_handler1_
Value handler that allows get value of 1D elements.
Definition: field_fe.hh:249
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:216
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
void set(unsigned int pos, double val)
Set value on given position.
Definition: vector_mpi.hh:111
static const Input::Type::Selection & get_disc_selection_input_type()
Definition: field_fe.cc:92
double get(unsigned int pos) const
Return value on given position.
Definition: vector_mpi.hh:105
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: index_types.hh:24
bool is_constant_in_space_
Flag detects that field is only dependent on time.
void value_list(const Armor::array &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list)
Returns std::vector of scalar values in several points at once.
Value value_
Last value, prevents passing large values (vectors) by value.
Input::Record in_rec_
Accessor to Input::Record.
Definition: field_fe.hh:274
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
Definition: asserts.hh:336
const Selection & close() const
Close the Selection, no more values can be added.
void normalize(unsigned int pos, double divisor)
Normalize value on given position.
Definition: vector_mpi.hh:117
FieldFE(unsigned int n_comp=0)
Definition: field_fe.cc:133
#define ASSERT_DBG(expr)
std::shared_ptr< EvalPoints > eval_points() const
Getter of eval_points object.
MixedPtr< FESystem > mixed_fe_system(MixedPtr< FiniteElement > fe, Args &&...args)
Definition: fe_system.hh:176
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:270
void cache_update(FieldValueCache< typename Value::element_type > &data_cache, ElementCacheMap &cache_map, unsigned int region_patch_idx) override
Definition: field_fe.cc:258
void set_boundary_dofs_vector(std::shared_ptr< std::vector< IntIdx > > boundary_dofs)
TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
#define END_TIMER(tag)
Ends a timer with specified tag.
std::shared_ptr< std::vector< LongIdx > > source_target_mesh_elm_map_
Maps element indices between source (data) and target (computational) mesh if data interpolation is s...
Definition: field_fe.hh:290
VectorMPI & vec()
Definition: field_fe.hh:158
DataInterpolation
Definition: field_fe.hh:58
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:231
virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData &init_data)
Definition: field_fe.cc:318
static std::shared_ptr< Mesh > get_mesh(const FilePath &file_path)
Definition: reader_cache.cc:40
Record type proxy class.
Definition: type_record.hh:182
VectorMPI data_vec
Store data of Field.
void cache_reinit(const ElementCacheMap &cache_map) override
Definition: field_fe.cc:295
Accessor to the data with type Type::Tuple.
Definition: accessors.hh:411
static ElementMap element_map(ElementAccessor< 3 > elm)
Definition: mapping_p1.cc:48
unsigned int n_comp() const
static constexpr Mask equation_input
The field is data parameter of the owning equation. (default on)
Definition: field_flag.hh:33
unsigned int dim() const
Definition: accessors.hh:157
static bool contains_point(arma::vec point, ElementAccessor< 3 > elm)
Test if element contains given point.
Definition: mapping_p1.cc:96
Internal auxiliary class representing intersection object of simplex<dimA> and simplex<dimB>.
static std::shared_ptr< BaseMeshReader > get_reader(const FilePath &file_path)
Definition: reader_cache.cc:36
unsigned int n_comp
number of components
unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:87
void calculate_native_values(ElementDataCache< double >::ComponentDataPtr data_cache)
Calculate native data over all elements of target mesh.
Definition: field_fe.cc:725
FieldFlag::Flags flags_
Field flags.
Definition: field_fe.hh:268
Initialization structure of FEValueHandler class.
unsigned int range_begin
Holds begin of component range of evaluation.
Class for declaration of the input data that are in string format.
Definition: type_base.hh:582
#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.
static std::shared_ptr< std::vector< LongIdx > > get_target_mesh_element_map(const FilePath &file_path, Mesh *computational_mesh)
Definition: reader_cache.cc:79
void make_dof_handler(const Mesh *mesh)
Create DofHandler object.
Definition: field_fe.cc:400
Implementation of range helper class.
void native_data_to_cache(ElementDataCache< double > &output_data_cache)
Definition: field_fe.cc:861
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
Definition: global_defs.h:157
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual)
Definition: asserts.hh:328
static constexpr Mask declare_input
The field can be set from input. The key in input field descriptor is declared. (default on) ...
Definition: field_flag.hh:35
static VectorMPI sequential(unsigned int size)
Definition: vector_mpi.cc:44
Definition: field.hh:63
unsigned int n_comp() const
Internal class representing intersection object.
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:300
void calculate_identic_values(ElementDataCache< double >::ComponentDataPtr data_cache)
Calculate data of identict_mesh interpolation on input data over all elements of target mesh...
Definition: field_fe.cc:754