Flow123d  master-7bf36fe
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  /// Declaration of exception.
69  TYPEDEF_ERR_INFO( EI_Field, std::string);
70  TYPEDEF_ERR_INFO( EI_File, std::string);
71  TYPEDEF_ERR_INFO( EI_ElemIdx, unsigned int);
72  TYPEDEF_ERR_INFO( EI_Region, std::string);
73  DECLARE_EXCEPTION( ExcInvalidElemeDim, << "Dimension of element in target mesh must be 0, 1 or 2! elm.idx() = " << EI_ElemIdx::val << ".\n" );
74  DECLARE_INPUT_EXCEPTION( ExcUndefElementValue,
75  << "FieldFE " << EI_Field::qval << " on region " << EI_Region::qval << " have invalid value .\n"
76  << "Provided by file " << EI_File::qval << " at element ID " << EI_ElemIdx::val << ".\n"
77  << "Please specify in default_value key.\n");
78 
79 
80 
81  static const unsigned int undef_uint = -1;
82 
83  /**
84  * Default constructor, optionally we need number of components @p n_comp in the case of Vector valued fields.
85  */
86  FieldFE(unsigned int n_comp=0);
87 
88  /**
89  * Return Record for initialization of FieldFE that is derived from Abstract
90  */
91  static const Input::Type::Record & get_input_type();
92 
93  /**
94  * Return Input selection for discretization type (determines the section of VTK file).
95  */
97 
98  /**
99  * Return Input selection that allow to set interpolation of input data.
100  */
102 
103  /**
104  * Factory class (descendant of @p Field<...>::FactoryBase) that is necessary
105  * for setting pressure values are piezometric head values.
106  */
108  public:
109  /// Constructor.
110  NativeFactory(unsigned int index, std::shared_ptr<DOFHandlerMultiDim> conc_dof_handler, VectorMPI dof_vector = VectorMPI::sequential(0))
111  : index_(index),
112  conc_dof_handler_(conc_dof_handler),
113  dof_vector_(dof_vector)
114  {}
115 
116  typename Field<spacedim,Value>::FieldBasePtr create_field(Input::Record rec, const FieldCommon &field) override;
117 
118  unsigned int index_;
119  std::shared_ptr<DOFHandlerMultiDim> conc_dof_handler_;
121  };
122 
123 
124  /**
125  * Setter for the finite element data.
126  * @param dh Dof handler.
127  * @param dof_values Vector of dof values, optional (create own vector according to dofhandler).
128  * @param block_index Index of block (in FESystem) or '-1' for FieldFE under simple DOF handler.
129  * @return Data vector of dof values.
130  */
131  VectorMPI set_fe_data(std::shared_ptr<DOFHandlerMultiDim> dh, VectorMPI dof_values = VectorMPI::sequential(0),
132  unsigned int block_index = FieldFE<spacedim, Value>::undef_uint);
133 
134  /**
135  * Overload @p FieldAlgorithmBase::cache_update
136  */
138  ElementCacheMap &cache_map, unsigned int region_patch_idx) override;
139 
140  /**
141  * Overload @p FieldAlgorithmBase::cache_reinit
142  *
143  * Reinit fe_values_ data member.
144  */
145  void cache_reinit(const ElementCacheMap &cache_map) override;
146 
147  /**
148  * Initialization from the input interface.
149  */
150  virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData& init_data);
151 
152 
153  /**
154  * Update time and possibly update data from GMSH file.
155  */
156  bool set_time(const TimeStep &time) override;
157 
158 
159  /**
160  * Set target mesh.
161  */
162  void set_mesh(const Mesh *mesh, bool boundary_domain) override;
163 
164 
165  /**
166  * Copy data vector to given output ElementDataCache
167  */
168  void native_data_to_cache(ElementDataCache<double> &output_data_cache);
169 
170 
171  /**
172  * Return size of vector of data stored in Field
173  */
174  unsigned int data_size() const;
175 
176  inline std::shared_ptr<DOFHandlerMultiDim> get_dofhandler() const {
177  return dh_;
178  }
179 
180  inline VectorMPI& vec() {
181  return data_vec_;
182  }
183 
184  inline const VectorMPI& vec() const {
185  return data_vec_;
186  }
187 
188  /// Call begin scatter functions (local to ghost) on data vector
190 
191  /// Call end scatter functions (local to ghost) on data vector
193 
194  /// Destructor.
195  virtual ~FieldFE();
196 
197 private:
198  /**
199  * Helper class holds specific data of field evaluation.
200  */
201  class FEItem {
202  public:
203  unsigned int comp_index_;
204  unsigned int range_begin_;
205  unsigned int range_end_;
206  };
207 
208  /**
209  * Helper class holds data of invalid values of all regions.
210  *
211  * If region contains invalid element value (typically 'not a number') is_invalid_ flag is set to true
212  * and other information (element id an value) are stored. Check of invalid values is performed
213  * during processing data of reader cache and possible exception is thrown only if FieldFE is defined
214  * on appropriate region.
215  *
216  * Data of all regions are stored in vector of RegionValueErr instances.
217  */
219  public:
220  /// Default constructor, sets valid region
222 
223  /// Constructor, sets invalid region, element and value specification
224  RegionValueErr(const std::string &region_name, unsigned int elm_id, double value) {
225  is_invalid_ = true;
226  region_name_ = region_name;
227  elm_id_ = elm_id;
228  value_ = value;
229  }
230 
232  std::string region_name_;
233  unsigned int elm_id_;
234  double value_;
235  };
236 
237  /// Create DofHandler object
238  void make_dof_handler(const MeshBase *mesh);
239 
240  /// Interpolate data (use Gaussian distribution) over all elements of target mesh.
241  void interpolate_gauss();
242 
243  /// Interpolate data (use intersection library) over all elements of target mesh.
245 
246 // /// Calculate native data over all elements of target mesh.
247 // void calculate_native_values(ElementDataCache<double>::CacheData data_cache);
248 //
249 // /// Calculate data of equivalent_mesh interpolation on input over all elements of target mesh.
250 // void calculate_equivalent_values(ElementDataCache<double>::CacheData data_cache);
251 
252  /// Calculate data of equivalent_mesh interpolation or native data on input over all elements of target mesh.
254 
255  /// Initialize FEValues object of given dimension.
256  template <unsigned int dim>
257  Quadrature init_quad(std::shared_ptr<EvalPoints> eval_points);
258 
260  unsigned int i_dof, unsigned int i_qp)
261  {
263  for (unsigned int c=0; c<Value::NRows_*Value::NCols_; ++c)
264  v(c/spacedim,c%spacedim) = fe_values_[dim].shape_value_component(i_dof, i_qp, c);
265  if (Value::NRows_ == Value::NCols_)
266  return v;
267  else
268  return v.t();
269  }
270 
271  template<unsigned int dim>
272  void fill_fe_system_data(unsigned int block_index) {
273  auto fe_system_ptr = std::dynamic_pointer_cast<FESystem<dim>>( dh_->ds()->fe()[Dim<dim>{}] );
274  ASSERT(fe_system_ptr != nullptr).error("Wrong type, must be FESystem!\n");
275  this->fe_item_[dim].comp_index_ = fe_system_ptr->function_space()->dof_indices()[block_index].component_offset;
276  this->fe_item_[dim].range_begin_ = fe_system_ptr->fe_dofs(block_index)[0];
277  this->fe_item_[dim].range_end_ = this->fe_item_[dim].range_begin_ + fe_system_ptr->fe()[block_index]->n_dofs();
278  }
279 
280 
281  template<unsigned int dim>
282  void fill_fe_item() {
283  this->fe_item_[dim].comp_index_ = 0;
284  this->fe_item_[dim].range_begin_ = 0;
285  this->fe_item_[dim].range_end_ = dh_->ds()->fe()[Dim<dim>{}]->n_dofs();
286  }
287 
288  /**
289  * Method computes value of given input cache element.
290  *
291  * If computed value is invalid (e.g. NaN value) sets the data specifying error value.
292  * @param i_cache_el Index of element of input ElementDataCache
293  * @param elm_idx Idx of element of computational mesh
294  * @param region_name Region of computational mesh
295  * @param actual_compute_region_error Data object holding data of region with invalid value.
296  */
297  double get_scaled_value(int i_cache_el, unsigned int elm_idx, const std::string &region_name, RegionValueErr &actual_compute_region_error);
298 
299 
300  /// DOF handler object
301  std::shared_ptr<DOFHandlerMultiDim> dh_;
302  /// Store data of Field
304 
305  /// Value handler that allows get value of 0D elements.
307  /// Value handler that allows get value of 1D elements.
309  /// Value handler that allows get value of 2D elements.
311  /// Value handler that allows get value of 3D elements.
313 
314  /// mesh reader file
316 
317  /// field name read from input
318  std::string field_name_;
319 
320  /// Specify section where to find the field data in input mesh file
322 
323  /// Specify type of FE data interpolation
325 
326  /// Field flags.
328 
329  /// Default value of element if not set in mesh data file
331 
332  /// Accessor to Input::Record
334 
335  /// Is set in set_mesh method. Value true means, that we accept only boundary element accessors in the @p value method.
337 
338  /// List of FEValues objects of dimensions 0,1,2,3 used for value calculation
340 
341  /// Maps element indices from computational mesh to the source (data).
342  std::shared_ptr<EquivalentMeshMap> source_target_mesh_elm_map_;
343 
344  /// Holds specific data of field evaluation over all dimensions.
345  std::array<FEItem, 4> fe_item_;
347 
348  /// Set holds data of valid / invalid element values on all regions
350 
351  /// Input ElementDataCache is stored in set_time and used in all evaluation and interpolation methods.
353 
354  /// Registrar of class to factory
355  static const int registrar;
356 };
357 
358 
359 /** Create FieldFE from dof handler */
360 template <int spacedim, class Value>
361 std::shared_ptr<FieldFE<spacedim, Value> > create_field_fe(std::shared_ptr<DOFHandlerMultiDim> dh,
362  VectorMPI *vec = nullptr,
363  unsigned int block_index = FieldFE<spacedim, Value>::undef_uint)
364 {
365  // Construct FieldFE
366  std::shared_ptr< FieldFE<spacedim, Value> > field_ptr = std::make_shared< FieldFE<spacedim, Value> >();
367  if (vec == nullptr)
368  field_ptr->set_fe_data( dh, dh->create_vector(), block_index );
369  else
370  field_ptr->set_fe_data( dh, *vec, block_index );
371 
372  return field_ptr;
373 }
374 
375 
376 /** Create FieldFE with parallel VectorMPI from finite element */
377 template <int spacedim, class Value>
378 std::shared_ptr<FieldFE<spacedim, Value> > create_field_fe(Mesh & mesh, const MixedPtr<FiniteElement> &fe)
379 {
380  // Prepare DOF handler
381  std::shared_ptr<DOFHandlerMultiDim> dh_par = std::make_shared<DOFHandlerMultiDim>(mesh);
382  std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( &mesh, fe);
383  dh_par->distribute_dofs(ds);
384 
385  return create_field_fe<spacedim,Value>(dh_par);
386 }
387 
388 
389 #endif /* FIELD_FE_HH_ */
FieldFE::data_size
unsigned int data_size() const
Definition: field_fe.cc:815
FieldFE::fill_fe_item
void fill_fe_item()
Definition: field_fe.hh:282
FieldFE::field_name_
std::string field_name_
field name read from input
Definition: field_fe.hh:318
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:342
FieldFE::flags_
FieldFlag::Flags flags_
Field flags.
Definition: field_fe.hh:327
element_data_cache.hh
FieldFE::FieldFE
FieldFE(unsigned int n_comp=0)
Definition: field_fe.cc:133
FieldFE::RegionValueErr::elm_id_
unsigned int elm_id_
Definition: field_fe.hh:233
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:119
FieldFE::TYPEDEF_ERR_INFO
TYPEDEF_ERR_INFO(EI_Field, std::string)
Declaration of exception.
factory.hh
bih_tree.hh
vector_mpi.hh
FieldFE::FEItem::range_begin_
unsigned int range_begin_
Definition: field_fe.hh:204
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:829
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:306
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:237
msh_basereader.hh
FieldFE::RegionValueErr::RegionValueErr
RegionValueErr()
Default constructor, sets valid region.
Definition: field_fe.hh:221
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:259
ElementDataCache< double >
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:355
FieldFE::RegionValueErr
Definition: field_fe.hh:218
point.hh
ElementCacheMap
Directing class of FieldValueCache.
Definition: field_value_cache.hh:152
value
static constexpr bool value
Definition: json.hpp:87
FieldFE::input_data_cache_
ElementDataCache< double >::CacheData input_data_cache_
Input ElementDataCache is stored in set_time and used in all evaluation and interpolation methods.
Definition: field_fe.hh:352
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
FieldFE::in_rec_
Input::Record in_rec_
Accessor to Input::Record.
Definition: field_fe.hh:333
system.hh
FieldFE::region_value_err_
std::vector< RegionValueErr > region_value_err_
Set holds data of valid / invalid element values on all regions.
Definition: field_fe.hh:349
FieldFE::NativeFactory::dof_vector_
VectorMPI dof_vector_
Definition: field_fe.hh:120
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:822
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:336
FieldAlgorithmBase::Point
Space< spacedim >::Point Point
Definition: field_algo_base.hh:108
index_types.hh
Field::FactoryBase
Definition: field.hh:122
fe_system.hh
Class FESystem for compound finite elements.
FieldFE::FEItem::range_end_
unsigned int range_end_
Definition: field_fe.hh:205
FieldFE::NativeFactory
Definition: field_fe.hh:107
FieldFE::data_vec_
VectorMPI data_vec_
Store data of Field.
Definition: field_fe.hh:303
FieldFE::interpolate_gauss
void interpolate_gauss()
Interpolate data (use Gaussian distribution) over all elements of target mesh.
Definition: field_fe.cc:475
Dim
Definition: mixed.hh:25
FieldFE::cache_reinit
void cache_reinit(const ElementCacheMap &cache_map) override
Definition: field_fe.cc:278
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:301
FieldFE::get_dofhandler
std::shared_ptr< DOFHandlerMultiDim > get_dofhandler() const
Definition: field_fe.hh:176
FieldFE::fe_
MixedPtr< FiniteElement > fe_
Definition: field_fe.hh:346
FieldFE::dh_
std::shared_ptr< DOFHandlerMultiDim > dh_
DOF handler object.
Definition: field_fe.hh:301
FieldFE::NativeFactory::index_
unsigned int index_
Definition: field_fe.hh:118
FieldFE::calculate_element_values
void calculate_element_values()
Calculate data of equivalent_mesh interpolation or native data on input over all elements of target m...
Definition: field_fe.cc:674
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::vec
const VectorMPI & vec() const
Definition: field_fe.hh:184
FieldFE::set_mesh
void set_mesh(const Mesh *mesh, bool boundary_domain) override
Definition: field_fe.cc:324
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:330
FieldFE::value_handler3_
FEValueHandler< 3, spacedim, Value > value_handler3_
Value handler that allows get value of 3D elements.
Definition: field_fe.hh:312
FieldAlgoBaseInitData
Helper struct stores data for initizalize descentants of FieldAlgorithmBase.
Definition: field_algo_base.hh:81
FieldFE::RegionValueErr::is_invalid_
bool is_invalid_
Definition: field_fe.hh:231
FEValueHandler< 0, spacedim, Value >
Definition: fe_value_handler.hh:95
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:180
FieldCommon
Common abstract parent of all Field<...> classes.
Definition: field_common.hh:77
FieldFE::interpolation_
DataInterpolation interpolation_
Specify type of FE data interpolation.
Definition: field_fe.hh:324
FieldFE::fe_item_
std::array< FEItem, 4 > fe_item_
Holds specific data of field evaluation over all dimensions.
Definition: field_fe.hh:345
FieldFE::discretization_
OutputTime::DiscreteSpace discretization_
Specify section where to find the field data in input mesh file.
Definition: field_fe.hh:321
FieldFE::RegionValueErr::RegionValueErr
RegionValueErr(const std::string &region_name, unsigned int elm_id, double value)
Constructor, sets invalid region, element and value specification.
Definition: field_fe.hh:224
mesh.h
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
Field::FieldBasePtr
std::shared_ptr< FieldBaseType > FieldBasePtr
Definition: field.hh:96
FieldFE::native_data_to_cache
void native_data_to_cache(ElementDataCache< double > &output_data_cache)
Definition: field_fe.cc:796
FieldFE
Definition: field.hh:60
Mesh
Definition: mesh.h:362
FieldFE::FEItem
Definition: field_fe.hh:201
FieldFE::value_handler2_
FEValueHandler< 2, spacedim, Value > value_handler2_
Value handler that allows get value of 2D elements.
Definition: field_fe.hh:310
FieldFE::FactoryBaseType
FieldAlgorithmBase< spacedim, Value > FactoryBaseType
Definition: field_fe.hh:54
FieldAlgorithmBase
Definition: field_algo_base.hh:105
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
FieldFE::interpolate_intersection
void interpolate_intersection()
Interpolate data (use intersection library) over all elements of target mesh.
Definition: field_fe.cc:571
FieldFE::DECLARE_INPUT_EXCEPTION
DECLARE_INPUT_EXCEPTION(ExcUndefElementValue,<< "FieldFE "<< EI_Field::qval<< " on region "<< EI_Region::qval<< " have invalid value .\n"<< "Provided by file "<< EI_File::qval<< " at element ID "<< EI_ElemIdx::val<< ".\n"<< "Please specify in default_value key.\n")
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:315
MeshBase
Base class for Mesh and BCMesh.
Definition: mesh.h:96
MixedPtr< FiniteElement >
FieldFE::FEItem::comp_index_
unsigned int comp_index_
Definition: field_fe.hh:203
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:361
FieldFE::fill_fe_system_data
void fill_fe_system_data(unsigned int block_index)
Definition: field_fe.hh:272
Armor::Array
Definition: armor.hh:597
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:339
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::set_time
bool set_time(const TimeStep &time) override
Definition: field_fe.cc:424
OutputTime::DiscreteSpace
DiscreteSpace
Definition: output_time.hh:108
FieldFE::undef_uint
static const unsigned int undef_uint
Definition: field_fe.hh:81
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:872
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:361
Armor::ArmaMat
typename arma::Mat< Type >::template fixed< nr, nc > ArmaMat
Definition: armor.hh:502
FieldFE::get_scaled_value
double get_scaled_value(int i_cache_el, unsigned int elm_idx, const std::string &region_name, RegionValueErr &actual_compute_region_error)
Definition: field_fe.cc:836
FieldFE::RegionValueErr::value_
double value_
Definition: field_fe.hh:234
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:110
FieldFE::value_handler1_
FEValueHandler< 1, spacedim, Value > value_handler1_
Value handler that allows get value of 1D elements.
Definition: field_fe.hh:308
FieldFE::init_quad
Quadrature init_quad(std::shared_ptr< EvalPoints > eval_points)
Initialize FEValues object of given dimension.
Definition: field_fe.cc:291
FieldFE::RegionValueErr::region_name_
std::string region_name_
Definition: field_fe.hh:232
range_wrapper.hh
Implementation of range helper class.