Flow123d  master-3768d5dec
field_fe.hh
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.hh
15  * @brief
16  */
17 
18 #ifndef FIELD_FE_HH_
19 #define FIELD_FE_HH_
20 
21 #include "petscmat.h"
22 #include "system/system.hh"
23 #include "system/index_types.hh"
26 #include "fields/field.hh"
27 #include "la/vector_mpi.hh"
28 #include "mesh/mesh.h"
29 #include "mesh/point.hh"
30 #include "mesh/bih_tree.hh"
31 #include "mesh/range_wrapper.hh"
32 #include "io/element_data_cache.hh"
33 #include "io/msh_basereader.hh"
34 #include "fem/fe_p.hh"
35 #include "fem/fe_system.hh"
36 #include "fem/dofhandler.hh"
37 #include "fem/finite_element.hh"
38 #include "fem/dh_cell_accessor.hh"
39 #include "input/factory.hh"
40 
41 #include <memory>
42 
43 
44 
45 /**
46  * Class representing fields given by finite element approximation.
47  *
48  */
49 template <int spacedim, class Value>
50 class FieldFE : public FieldAlgorithmBase<spacedim, Value>
51 {
52 public:
56 
57  /**
58  * Possible interpolations of input data.
59  */
61  {
62  identic_msh, //!< identical mesh
63  equivalent_msh, //!< equivalent mesh (default value)
64  gauss_p0, //!< P0 interpolation (with the use of Gaussian distribution)
65  interp_p0 //!< P0 interpolation (with the use of calculation of intersections)
66  };
67 
68  TYPEDEF_ERR_INFO( EI_ElemIdx, unsigned int);
69  DECLARE_EXCEPTION( ExcInvalidElemeDim, << "Dimension of element in target mesh must be 0, 1 or 2! elm.idx() = " << EI_ElemIdx::val << ".\n" );
70 
71  static const unsigned int undef_uint = -1;
72 
73  /**
74  * Default constructor, optionally we need number of components @p n_comp in the case of Vector valued fields.
75  */
76  FieldFE(unsigned int n_comp=0);
77 
78  /**
79  * Return Record for initialization of FieldFE that is derived from Abstract
80  */
81  static const Input::Type::Record & get_input_type();
82 
83  /**
84  * Return Input selection for discretization type (determines the section of VTK file).
85  */
87 
88  /**
89  * Return Input selection that allow to set interpolation of input data.
90  */
92 
93  /**
94  * Factory class (descendant of @p Field<...>::FactoryBase) that is necessary
95  * for setting pressure values are piezometric head values.
96  */
98  public:
99  /// Constructor.
100  NativeFactory(unsigned int index, std::shared_ptr<DOFHandlerMultiDim> conc_dof_handler, VectorMPI dof_vector = VectorMPI::sequential(0))
101  : index_(index),
102  conc_dof_handler_(conc_dof_handler),
103  dof_vector_(dof_vector)
104  {}
105 
106  typename Field<spacedim,Value>::FieldBasePtr create_field(Input::Record rec, const FieldCommon &field) override;
107 
108  unsigned int index_;
109  std::shared_ptr<DOFHandlerMultiDim> conc_dof_handler_;
111  };
112 
113 
114  /**
115  * Setter for the finite element data.
116  * @param dh Dof handler.
117  * @param dof_values Vector of dof values, optional (create own vector according to dofhandler).
118  * @param block_index Index of block (in FESystem) or '-1' for FieldFE under simple DOF handler.
119  * @return Data vector of dof values.
120  */
121  VectorMPI set_fe_data(std::shared_ptr<DOFHandlerMultiDim> dh, VectorMPI dof_values = VectorMPI::sequential(0),
122  unsigned int block_index = FieldFE<spacedim, Value>::undef_uint);
123 
124  /**
125  * Returns one value in one given point. ResultType can be used to avoid some costly calculation if the result is trivial.
126  */
127  virtual typename Value::return_type const &value(const Point &p, const ElementAccessor<spacedim> &elm);
128 
129  /**
130  * Returns std::vector of scalar values in several points at once.
131  */
132  virtual void value_list (const Armor::array &point_list, const ElementAccessor<spacedim> &elm,
134 
135  /**
136  * Overload @p FieldAlgorithmBase::cache_update
137  */
139  ElementCacheMap &cache_map, unsigned int region_patch_idx) override;
140 
141  /**
142  * Overload @p FieldAlgorithmBase::cache_reinit
143  *
144  * Reinit fe_values_ data member.
145  */
146  void cache_reinit(const ElementCacheMap &cache_map) override;
147 
148  /**
149  * Initialization from the input interface.
150  */
151  virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData& init_data);
152 
153 
154  /**
155  * Update time and possibly update data from GMSH file.
156  */
157  bool set_time(const TimeStep &time) override;
158 
159 
160  /**
161  * Set target mesh.
162  */
163  void set_mesh(const Mesh *mesh, bool boundary_domain) override;
164 
165 
166  /**
167  * Copy data vector to given output ElementDataCache
168  */
169  void native_data_to_cache(ElementDataCache<double> &output_data_cache);
170 
171 
172  /**
173  * Return size of vector of data stored in Field
174  */
175  unsigned int data_size() const;
176 
177  inline std::shared_ptr<DOFHandlerMultiDim> get_dofhandler() const {
178  return dh_;
179  }
180 
181  inline VectorMPI& vec() {
182  return data_vec_;
183  }
184 
185  /// Call begin scatter functions (local to ghost) on data vector
187 
188  /// Call end scatter functions (local to ghost) on data vector
190 
191  /// Destructor.
192  virtual ~FieldFE();
193 
194 private:
195  /**
196  * Helper class holds specific data of field evaluation.
197  */
198  class FEItem {
199  public:
200  unsigned int comp_index_;
201  unsigned int range_begin_;
202  unsigned int range_end_;
203  };
204 
205  /// Create DofHandler object
206  void make_dof_handler(const MeshBase *mesh);
207 
208  /// Interpolate data (use Gaussian distribution) over all elements of target mesh.
210 
211  /// Interpolate data (use intersection library) over all elements of target mesh.
213 
214  /// Calculate native data over all elements of target mesh.
216 
217  /// Calculate data of identict_mesh interpolation on input data over all elements of target mesh.
219 
220  /// Calculate data of equivalent_mesh interpolation on input over all elements of target mesh.
222 
223  /**
224  * Fill data to boundary_dofs_ vector.
225  *
226  * TODO: Temporary solution. REMOVE this method and fix all places where is boundary_dofs_ vector used.
227  */
228  void fill_boundary_dofs();
229 
230  /// Initialize FEValues object of given dimension.
231  template <unsigned int dim>
232  Quadrature init_quad(std::shared_ptr<EvalPoints> eval_points);
233 
235  unsigned int i_dof, unsigned int i_qp)
236  {
238  for (unsigned int c=0; c<Value::NRows_*Value::NCols_; ++c)
239  v(c/spacedim,c%spacedim) = fe_values_[dim].shape_value_component(i_dof, i_qp, c);
240  if (Value::NRows_ == Value::NCols_)
241  return v;
242  else
243  return v.t();
244  }
245 
246  template<unsigned int dim>
247  void fill_fe_system_data(unsigned int block_index) {
248  auto fe_system_ptr = std::dynamic_pointer_cast<FESystem<dim>>( dh_->ds()->fe()[Dim<dim>{}] );
249  ASSERT(fe_system_ptr != nullptr).error("Wrong type, must be FESystem!\n");
250  this->fe_item_[dim].comp_index_ = fe_system_ptr->function_space()->dof_indices()[block_index].component_offset;
251  this->fe_item_[dim].range_begin_ = fe_system_ptr->fe_dofs(block_index)[0];
252  this->fe_item_[dim].range_end_ = this->fe_item_[dim].range_begin_ + fe_system_ptr->fe()[block_index]->n_dofs();
253  }
254 
255 
256  template<unsigned int dim>
257  void fill_fe_item() {
258  this->fe_item_[dim].comp_index_ = 0;
259  this->fe_item_[dim].range_begin_ = 0;
260  this->fe_item_[dim].range_end_ = dh_->ds()->fe()[Dim<dim>{}]->n_dofs();
261  }
262 
263 
264  /// DOF handler object
265  std::shared_ptr<DOFHandlerMultiDim> dh_;
266  /// Store data of Field
268 
269  /// Value handler that allows get value of 0D elements.
271  /// Value handler that allows get value of 1D elements.
273  /// Value handler that allows get value of 2D elements.
275  /// Value handler that allows get value of 3D elements.
277 
278  /// mesh reader file
280 
281  /// field name read from input
282  std::string field_name_;
283 
284  /// Specify section where to find the field data in input mesh file
286 
287  /// Specify type of FE data interpolation
289 
290  /// Field flags.
292 
293  /// Default value of element if not set in mesh data file
295 
296  /// Accessor to Input::Record
298 
299  /// Is set in set_mesh method. Value true means, that we accept only boundary element accessors in the @p value method.
301 
302  /**
303  * Hold dofs of boundary elements.
304  *
305  * TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
306  */
307  std::shared_ptr< std::vector<IntIdx> > boundary_dofs_;
308 
309  /// List of FEValues objects of dimensions 0,1,2,3 used for value calculation
311 
312  /// Maps element indices from computational mesh to the source (data).
313  std::shared_ptr<EquivalentMeshMap> source_target_mesh_elm_map_;
314 
315  /// Holds specific data of field evaluation over all dimensions.
316  std::array<FEItem, 4> fe_item_;
318 
319  /// Registrar of class to factory
320  static const int registrar;
321 };
322 
323 
324 /** Create FieldFE from dof handler */
325 template <int spacedim, class Value>
326 std::shared_ptr<FieldFE<spacedim, Value> > create_field_fe(std::shared_ptr<DOFHandlerMultiDim> dh,
327  VectorMPI *vec = nullptr,
328  unsigned int block_index = FieldFE<spacedim, Value>::undef_uint)
329 {
330  // Construct FieldFE
331  std::shared_ptr< FieldFE<spacedim, Value> > field_ptr = std::make_shared< FieldFE<spacedim, Value> >();
332  if (vec == nullptr)
333  field_ptr->set_fe_data( dh, dh->create_vector(), block_index );
334  else
335  field_ptr->set_fe_data( dh, *vec, block_index );
336 
337  return field_ptr;
338 }
339 
340 
341 /** Create FieldFE with parallel VectorMPI from finite element */
342 template <int spacedim, class Value>
343 std::shared_ptr<FieldFE<spacedim, Value> > create_field_fe(Mesh & mesh, const MixedPtr<FiniteElement> &fe)
344 {
345  // Prepare DOF handler
346  std::shared_ptr<DOFHandlerMultiDim> dh_par = std::make_shared<DOFHandlerMultiDim>(mesh);
347  std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( &mesh, fe);
348  dh_par->distribute_dofs(ds);
349 
350  return create_field_fe<spacedim,Value>(dh_par);
351 }
352 
353 
354 #endif /* FIELD_FE_HH_ */
FieldFE::data_size
unsigned int data_size() const
Definition: field_fe.cc:868
FieldFE::fill_fe_item
void fill_fe_item()
Definition: field_fe.hh:257
FieldFE::field_name_
std::string field_name_
field name read from input
Definition: field_fe.hh:282
FieldFE::boundary_dofs_
std::shared_ptr< std::vector< IntIdx > > boundary_dofs_
Definition: field_fe.hh:307
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:313
FieldFE::flags_
FieldFlag::Flags flags_
Field flags.
Definition: field_fe.hh:291
element_data_cache.hh
FieldFE::FieldFE
FieldFE(unsigned int n_comp=0)
Definition: field_fe.cc:133
Armor::vec
ArmaVec< double, N > vec
Definition: armor.hh:885
FieldFE::NativeFactory::conc_dof_handler_
std::shared_ptr< DOFHandlerMultiDim > conc_dof_handler_
Definition: field_fe.hh:109
factory.hh
bih_tree.hh
vector_mpi.hh
FieldFE::FEItem::range_begin_
unsigned int range_begin_
Definition: field_fe.hh:201
FlagArray< FieldFlag >
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:882
field_algo_base.hh
FieldFE::value_handler0_
FEValueHandler< 0, spacedim, Value > value_handler0_
Value handler that allows get value of 0D elements.
Definition: field_fe.hh:270
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
msh_basereader.hh
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:234
ElementDataCache
Definition: element_data_cache.hh:44
ASSERT
#define ASSERT(expr)
Definition: asserts.hh:351
FieldFE::DataInterpolation
DataInterpolation
Definition: field_fe.hh:60
FieldFE::registrar
static const int registrar
Registrar of class to factory.
Definition: field_fe.hh:320
point.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:755
ElementCacheMap
Directing class of FieldValueCache.
Definition: field_value_cache.hh:151
FilePath
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
FieldFE::equivalent_msh
@ equivalent_msh
equivalent mesh (default value)
Definition: field_fe.hh:63
std::vector
Definition: doxy_dummy_defs.hh:7
ElementAccessor
Definition: dh_cell_accessor.hh:32
FieldFE::in_rec_
Input::Record in_rec_
Accessor to Input::Record.
Definition: field_fe.hh:297
system.hh
FieldFE::NativeFactory::dof_vector_
VectorMPI dof_vector_
Definition: field_fe.hh:110
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:875
dofhandler.hh
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
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:300
FieldAlgorithmBase::Point
Space< spacedim >::Point Point
Definition: field_algo_base.hh:115
index_types.hh
Field::FactoryBase
Definition: field.hh:121
fe_system.hh
Class FESystem for compound finite elements.
FieldFE::FEItem::range_end_
unsigned int range_end_
Definition: field_fe.hh:202
FieldFE::NativeFactory
Definition: field_fe.hh:97
FieldFE::data_vec_
VectorMPI data_vec_
Store data of Field.
Definition: field_fe.hh:267
Dim
Definition: mixed.hh:25
FieldFE::cache_reinit
void cache_reinit(const ElementCacheMap &cache_map) override
Definition: field_fe.cc:325
fe_p.hh
Definitions of basic Lagrangean finite elements with polynomial shape functions.
FieldFE::init_from_input
virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData &init_data)
Definition: field_fe.cc:348
FieldFE::get_dofhandler
std::shared_ptr< DOFHandlerMultiDim > get_dofhandler() const
Definition: field_fe.hh:177
FieldFE::fe_
MixedPtr< FiniteElement > fe_
Definition: field_fe.hh:317
FieldFE::dh_
std::shared_ptr< DOFHandlerMultiDim > dh_
DOF handler object.
Definition: field_fe.hh:265
FieldFE::NativeFactory::index_
unsigned int index_
Definition: field_fe.hh:108
FieldFE::value
virtual const Value::return_type & value(const Point &p, const ElementAccessor< spacedim > &elm)
Definition: field_fe.cc:237
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
FieldFE::set_mesh
void set_mesh(const Mesh *mesh, bool boundary_domain) override
Definition: field_fe.cc:371
FieldFE::interp_p0
@ interp_p0
P0 interpolation (with the use of calculation of intersections)
Definition: field_fe.hh:65
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:294
FieldFE::value_handler3_
FEValueHandler< 3, spacedim, Value > value_handler3_
Value handler that allows get value of 3D elements.
Definition: field_fe.hh:276
FieldAlgoBaseInitData
Helper struct stores data for initizalize descentants of FieldAlgorithmBase.
Definition: field_algo_base.hh:81
FEValueHandler< 0, spacedim, Value >
Definition: fe_value_handler.hh:121
TimeStep
Representation of one time step..
Definition: time_governor.hh:123
finite_element.hh
Abstract class for description of finite elements.
FieldFE::vec
VectorMPI & vec()
Definition: field_fe.hh:181
FieldCommon
Common abstract parent of all Field<...> classes.
Definition: field_common.hh:76
FieldFE::interpolation_
DataInterpolation interpolation_
Specify type of FE data interpolation.
Definition: field_fe.hh:288
FieldFE::fe_item_
std::array< FEItem, 4 > fe_item_
Holds specific data of field evaluation over all dimensions.
Definition: field_fe.hh:316
FieldFE::discretization_
OutputTime::DiscreteSpace discretization_
Specify section where to find the field data in input mesh file.
Definition: field_fe.hh:285
mesh.h
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:559
Input::Type::Selection
Template for classes storing finite set of named values.
Definition: type_selection.hh:65
FieldFE::gauss_p0
@ gauss_p0
P0 interpolation (with the use of Gaussian distribution)
Definition: field_fe.hh:64
Input::Type::Record
Record type proxy class.
Definition: type_record.hh:182
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:784
Field::FieldBasePtr
std::shared_ptr< FieldBaseType > FieldBasePtr
Definition: field.hh:95
FieldFE::native_data_to_cache
void native_data_to_cache(ElementDataCache< double > &output_data_cache)
Definition: field_fe.cc:849
FieldFE
Definition: field.hh:59
Mesh
Definition: mesh.h:361
FieldFE::FEItem
Definition: field_fe.hh:198
FieldFE::value_handler2_
FEValueHandler< 2, spacedim, Value > value_handler2_
Value handler that allows get value of 2D elements.
Definition: field_fe.hh:274
FieldFE::FactoryBaseType
FieldAlgorithmBase< spacedim, Value > FactoryBaseType
Definition: field_fe.hh:54
FieldAlgorithmBase
Definition: field_algo_base.hh:112
FEValueHandler< 1, spacedim, Value >
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
ElementDataCache::CacheData
std::shared_ptr< std::vector< T > > CacheData
Definition: element_data_cache.hh:52
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:279
MeshBase
Base class for Mesh and BCMesh.
Definition: mesh.h:96
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
FieldFE::FEItem::comp_index_
unsigned int comp_index_
Definition: field_fe.hh:200
FieldFE::get_input_type
static const Input::Type::Record & get_input_type()
Implementation.
Definition: field_fe.cc:61
create_field_fe
std::shared_ptr< FieldFE< spacedim, Value > > create_field_fe(std::shared_ptr< DOFHandlerMultiDim > dh, VectorMPI *vec=nullptr, unsigned int block_index=FieldFE< spacedim, Value >::undef_uint)
Definition: field_fe.hh:326
FieldFE::fill_fe_system_data
void fill_fe_system_data(unsigned int block_index)
Definition: field_fe.hh:247
Armor::Array< double >
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:310
FieldFE::DECLARE_EXCEPTION
DECLARE_EXCEPTION(ExcInvalidElemeDim,<< "Dimension of element in target mesh must be 0, 1 or 2! elm.idx() = "<< EI_ElemIdx::val<< ".\n")
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
FieldFE::identic_msh
@ identic_msh
identical mesh
Definition: field_fe.hh:62
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:653
FieldFE::set_time
bool set_time(const TimeStep &time) override
Definition: field_fe.cc:497
OutputTime::DiscreteSpace
DiscreteSpace
Definition: output_time.hh:108
FieldFE::undef_uint
static const unsigned int undef_uint
Definition: field_fe.hh:71
FieldFE::fill_boundary_dofs
void fill_boundary_dofs()
Definition: field_fe.cc:407
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
FieldFE::~FieldFE
virtual ~FieldFE()
Destructor.
Definition: field_fe.cc:904
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:434
FieldFE::TYPEDEF_ERR_INFO
TYPEDEF_ERR_INFO(EI_ElemIdx, unsigned int)
Armor::ArmaMat
typename arma::Mat< Type >::template fixed< nr, nc > ArmaMat
Definition: armor.hh:502
FieldFE::FieldFactoryBaseType
Field< spacedim, Value >::FactoryBase FieldFactoryBaseType
Definition: field_fe.hh:55
field.hh
FieldFE::NativeFactory::NativeFactory
NativeFactory(unsigned int index, std::shared_ptr< DOFHandlerMultiDim > conc_dof_handler, VectorMPI dof_vector=VectorMPI::sequential(0))
Constructor.
Definition: field_fe.hh:100
FieldFE::value_handler1_
FEValueHandler< 1, spacedim, Value > value_handler1_
Value handler that allows get value of 1D elements.
Definition: field_fe.hh:272
FieldFE::init_quad
Quadrature init_quad(std::shared_ptr< EvalPoints > eval_points)
Initialize FEValues object of given dimension.
Definition: field_fe.cc:338
range_wrapper.hh
Implementation of range helper class.
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:810