Flow123d  master-a4da60c
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"
25 #include "fields/field.hh"
26 #include "la/vector_mpi.hh"
27 #include "mesh/mesh.h"
28 #include "mesh/point.hh"
29 #include "mesh/bih_tree.hh"
30 #include "mesh/range_wrapper.hh"
31 #include "io/element_data_cache.hh"
32 #include "io/msh_basereader.hh"
33 #include "fem/fe_p.hh"
34 #include "fem/fe_system.hh"
35 #include "fem/dofhandler.hh"
36 #include "fem/finite_element.hh"
37 #include "fem/dh_cell_accessor.hh"
38 #include "fem/mapping_p1.hh"
39 #include "input/factory.hh"
41 
42 #include <memory>
43 
44 
45 
46 /**
47  * Class representing fields given by finite element approximation.
48  *
49  */
50 template <int spacedim, class Value>
51 class FieldFE : public FieldAlgorithmBase<spacedim, Value>
52 {
53 public:
57 
58  /**
59  * Possible interpolations of input data.
60  */
62  {
63  identic_msh, //!< identical mesh
64  equivalent_msh, //!< equivalent mesh (default value)
65  gauss_p0, //!< P0 interpolation (with the use of Gaussian distribution)
66  interp_p0 //!< P0 interpolation (with the use of calculation of intersections)
67  };
68 
69  /// Declaration of exception.
70  TYPEDEF_ERR_INFO( EI_Field, std::string);
71  TYPEDEF_ERR_INFO( EI_File, std::string);
72  TYPEDEF_ERR_INFO( EI_ElemIdx, unsigned int);
73  TYPEDEF_ERR_INFO( EI_Region, std::string);
74  DECLARE_EXCEPTION( ExcInvalidElemeDim, << "Dimension of element in target mesh must be 0, 1 or 2! elm.idx() = " << EI_ElemIdx::val << ".\n" );
75  DECLARE_INPUT_EXCEPTION( ExcUndefElementValue,
76  << "FieldFE " << EI_Field::qval << " on region " << EI_Region::qval << " have invalid value .\n"
77  << "Provided by file " << EI_File::qval << " at element ID " << EI_ElemIdx::val << ".\n"
78  << "Please specify in default_value key.\n");
79 
80 
81 
82  static const unsigned int undef_uint = -1;
83 
84  /**
85  * Default constructor, optionally we need number of components @p n_comp in the case of Vector valued fields.
86  */
87  FieldFE(unsigned int n_comp=0);
88 
89  /**
90  * Return Record for initialization of FieldFE that is derived from Abstract
91  */
92  static const Input::Type::Record & get_input_type();
93 
94  /**
95  * Return Input selection for discretization type (determines the section of VTK file).
96  */
98 
99  /**
100  * Return Input selection that allow to set interpolation of input data.
101  */
103 
104  /**
105  * Factory class (descendant of @p Field<...>::FactoryBase) that is necessary
106  * for setting pressure values are piezometric head values.
107  */
109  public:
110  /// Constructor.
111  NativeFactory(unsigned int index, std::shared_ptr<DOFHandlerMultiDim> conc_dof_handler, VectorMPI dof_vector = VectorMPI::sequential(0))
112  : index_(index),
113  conc_dof_handler_(conc_dof_handler),
114  dof_vector_(dof_vector)
115  {}
116 
117  typename Field<spacedim,Value>::FieldBasePtr create_field(Input::Record rec, const FieldCommon &field) override;
118 
119  unsigned int index_;
120  std::shared_ptr<DOFHandlerMultiDim> conc_dof_handler_;
122  };
123 
124 
125  /**
126  * Setter for the finite element data.
127  * @param dh Dof handler.
128  * @param dof_values Vector of dof values, optional (create own vector according to dofhandler).
129  * @param block_index Index of block (in FESystem) or '-1' for FieldFE under simple DOF handler.
130  * @return Data vector of dof values.
131  */
132  VectorMPI set_fe_data(std::shared_ptr<DOFHandlerMultiDim> dh, VectorMPI dof_values = VectorMPI::sequential(0),
133  unsigned int block_index = FieldFE<spacedim, Value>::undef_uint);
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(AssemblyInternals &asm_internals) override;
147 
148 
149  /**
150  * Initialization from the input interface.
151  */
152  virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData& init_data);
153 
154 
155  /**
156  * Update time and possibly update data from GMSH file.
157  */
158  bool set_time(const TimeStep &time) override;
159 
160 
161  /**
162  * Set target mesh.
163  */
164  void set_mesh(const Mesh *mesh) override;
165 
166 
167  /**
168  * Copy data vector to given output ElementDataCache
169  */
170  void native_data_to_cache(ElementDataCache<double> &output_data_cache);
171 
172 
173  /**
174  * Return size of vector of data stored in Field
175  */
176  unsigned int data_size() const;
177 
178  inline std::shared_ptr<DOFHandlerMultiDim> get_dofhandler() const {
179  return dh_;
180  }
181 
182  inline VectorMPI& vec() {
183  return data_vec_;
184  }
185 
186  inline const VectorMPI& vec() const {
187  return data_vec_;
188  }
189 
190  /// Call begin scatter functions (local to ghost) on data vector
192 
193  /// Call end scatter functions (local to ghost) on data vector
195 
196  /// Destructor.
197  virtual ~FieldFE();
198 
199 private:
200  /**
201  * Helper class holds specific data of field evaluation.
202  */
203  class FEItem {
204  public:
205  unsigned int comp_index_;
206  unsigned int range_begin_;
207  unsigned int range_end_;
208  };
209 
210  /**
211  * Helper class holds data of invalid values of all regions.
212  *
213  * If region contains invalid element value (typically 'not a number') is_invalid_ flag is set to true
214  * and other information (element id an value) are stored. Check of invalid values is performed
215  * during processing data of reader cache and possible exception is thrown only if FieldFE is defined
216  * on appropriate region.
217  *
218  * Data of all regions are stored in vector of RegionValueErr instances.
219  */
221  public:
222  /// Default constructor, sets valid region
224 
225  /// Constructor, sets invalid region, element and value specification
226  RegionValueErr(const std::string &region_name, unsigned int elm_id, double value) {
227  is_invalid_ = true;
228  region_name_ = region_name;
229  elm_id_ = elm_id;
230  value_ = value;
231  }
232 
234  std::string region_name_;
235  unsigned int elm_id_;
236  double value_;
237  };
238 
239  /// Create DofHandler object
240  void make_dof_handler(const MeshBase *mesh);
241 
242  /// Interpolate data (use Gaussian distribution) over all elements of target mesh.
243  void interpolate_gauss();
244 
245  /// Interpolate data (use intersection library) over all elements of target mesh.
247 
248 // /// Calculate native data over all elements of target mesh.
249 // void calculate_native_values(ElementDataCache<double>::CacheData data_cache);
250 //
251 // /// Calculate data of equivalent_mesh interpolation on input over all elements of target mesh.
252 // void calculate_equivalent_values(ElementDataCache<double>::CacheData data_cache);
253 
254  /// Calculate data of equivalent_mesh interpolation or native data on input over all elements of target mesh.
256 
257  /// Initialize FEValues object of given dimension.
258  template <unsigned int dim>
259  Quadrature init_quad(std::shared_ptr<EvalPoints> eval_points);
260 
262  unsigned int i_dof, unsigned int i_qp)
263  {
265  for (unsigned int c=0; c<Value::NRows_*Value::NCols_; ++c)
266  v(c/spacedim,c%spacedim) = fe_values_[dim].shape_value_component(i_dof, i_qp, c);
267  if (Value::NRows_ == Value::NCols_)
268  return v;
269  else
270  return v.t();
271  }
272 
273  template<unsigned int dim>
274  void fill_fe_system_data(unsigned int block_index) {
275  auto fe_system_ptr = std::dynamic_pointer_cast<FESystem<dim>>( dh_->ds()->fe()[Dim<dim>{}] );
276  ASSERT(fe_system_ptr != nullptr).error("Wrong type, must be FESystem!\n");
277  this->fe_item_[dim].comp_index_ = fe_system_ptr->function_space()->dof_indices()[block_index].component_offset;
278  this->fe_item_[dim].range_begin_ = fe_system_ptr->fe_dofs(block_index)[0];
279  this->fe_item_[dim].range_end_ = this->fe_item_[dim].range_begin_ + fe_system_ptr->fe()[block_index]->n_dofs();
280  }
281 
282 
283  template<unsigned int dim>
284  void fill_fe_item() {
285  this->fe_item_[dim].comp_index_ = 0;
286  this->fe_item_[dim].range_begin_ = 0;
287  this->fe_item_[dim].range_end_ = dh_->ds()->fe()[Dim<dim>{}]->n_dofs();
288  }
289 
290  /**
291  * Method computes value of given input cache element.
292  *
293  * If computed value is invalid (e.g. NaN value) sets the data specifying error value.
294  * @param i_cache_el Index of element of input ElementDataCache
295  * @param elm_idx Idx of element of computational mesh
296  * @param region_name Region of computational mesh
297  * @param actual_compute_region_error Data object holding data of region with invalid value.
298  */
299  double get_scaled_value(int i_cache_el, unsigned int elm_idx, const std::string &region_name, RegionValueErr &actual_compute_region_error);
300 
301  /**
302  * Helper method. Compute real coordinates and weights (use QGauss) of given element.
303  *
304  * Method is needs in Gauss interpolation.
305  */
306  template<int elemdim>
307  unsigned int compute_fe_quadrature(std::vector<arma::vec::fixed<3>> & q_points, std::vector<double> & q_weights,
308  const ElementAccessor<spacedim> &elm, unsigned int order=3)
309  {
310  static_assert(elemdim <= spacedim, "Dimension of element must be less equal than spacedim.");
311  static const double weight_coefs[] = { 1., 1., 2., 6. };
312 
313  QGauss qgauss(elemdim, order);
315 
316  for(unsigned i=0; i<qgauss.size(); ++i) {
317  q_weights[i] = qgauss.weight(i)*weight_coefs[elemdim];
319  }
320 
321  return qgauss.size();
322  }
323 
324 
325  /// DOF handler object
326  std::shared_ptr<DOFHandlerMultiDim> dh_;
327  /// Store data of Field
329 
330  /// mesh reader file
332 
333  /// field name read from input
334  std::string field_name_;
335 
336  /// Specify section where to find the field data in input mesh file
338 
339  /// Specify type of FE data interpolation
341 
342  /// Field flags.
344 
345  /// Default value of element if not set in mesh data file
347 
348  /// Accessor to Input::Record
350 
351  /// Is set in set_mesh method. Value true means, that we accept only boundary element accessors in the @p value method.
353 
354  /// List of FEValues objects of dimensions 0,1,2,3 used for value calculation
356 
357  /// Maps element indices from computational mesh to the source (data).
358  std::shared_ptr<EquivalentMeshMap> source_target_mesh_elm_map_;
359 
360  /// Holds specific data of field evaluation over all dimensions.
361  std::array<FEItem, 4> fe_item_;
363 
364  /// Set holds data of valid / invalid element values on all regions
366 
367  /// Input ElementDataCache is stored in set_time and used in all evaluation and interpolation methods.
369 
370  /// Registrar of class to factory
371  static const int registrar;
372 };
373 
374 
375 /** Create FieldFE from dof handler */
376 template <int spacedim, class Value>
377 std::shared_ptr<FieldFE<spacedim, Value> > create_field_fe(std::shared_ptr<DOFHandlerMultiDim> dh,
378  VectorMPI *vec = nullptr,
379  unsigned int block_index = FieldFE<spacedim, Value>::undef_uint)
380 {
381  // Construct FieldFE
382  std::shared_ptr< FieldFE<spacedim, Value> > field_ptr = std::make_shared< FieldFE<spacedim, Value> >();
383  if (vec == nullptr)
384  field_ptr->set_fe_data( dh, dh->create_vector(), block_index );
385  else
386  field_ptr->set_fe_data( dh, *vec, block_index );
387 
388  return field_ptr;
389 }
390 
391 
392 /** Create FieldFE with parallel VectorMPI from finite element */
393 template <int spacedim, class Value>
394 std::shared_ptr<FieldFE<spacedim, Value> > create_field_fe(Mesh & mesh, const MixedPtr<FiniteElement> &fe)
395 {
396  // Prepare DOF handler
397  std::shared_ptr<DOFHandlerMultiDim> dh_par = std::make_shared<DOFHandlerMultiDim>(mesh);
398  std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( &mesh, fe);
399  dh_par->distribute_dofs(ds);
400 
401  return create_field_fe<spacedim,Value>(dh_par);
402 }
403 
404 
405 #endif /* FIELD_FE_HH_ */
#define ASSERT(expr)
Definition: asserts.hh:351
Directing class of FieldValueCache.
unsigned int n_comp() const
Space< spacedim >::Point Point
Common abstract parent of all Field<...> classes.
Definition: field_common.hh:78
unsigned int range_end_
Definition: field_fe.hh:207
unsigned int range_begin_
Definition: field_fe.hh:206
unsigned int comp_index_
Definition: field_fe.hh:205
Field< spacedim, Value >::FieldBasePtr create_field(Input::Record rec, const FieldCommon &field) override
Definition: field_fe.cc:144
std::shared_ptr< DOFHandlerMultiDim > conc_dof_handler_
Definition: field_fe.hh:120
NativeFactory(unsigned int index, std::shared_ptr< DOFHandlerMultiDim > conc_dof_handler, VectorMPI dof_vector=VectorMPI::sequential(0))
Constructor.
Definition: field_fe.hh:111
unsigned int index_
Definition: field_fe.hh:119
RegionValueErr(const std::string &region_name, unsigned int elm_id, double value)
Constructor, sets invalid region, element and value specification.
Definition: field_fe.hh:226
std::string region_name_
Definition: field_fe.hh:234
unsigned int elm_id_
Definition: field_fe.hh:235
RegionValueErr()
Default constructor, sets valid region.
Definition: field_fe.hh:223
FieldAlgorithmBase< spacedim, Value >::Point Point
Definition: field_fe.hh:54
DECLARE_EXCEPTION(ExcInvalidElemeDim,<< "Dimension of element in target mesh must be 0, 1 or 2! elm.idx() = "<< EI_ElemIdx::val<< ".\n")
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:261
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:797
MixedPtr< FiniteElement > fe_
Definition: field_fe.hh:362
FieldAlgorithmBase< spacedim, Value > FactoryBaseType
Definition: field_fe.hh:55
VectorMPI data_vec_
Store data of Field.
Definition: field_fe.hh:328
static const Input::Type::Selection & get_interp_selection_input_type()
Definition: field_fe.cc:105
unsigned int data_size() const
Definition: field_fe.cc:776
std::vector< RegionValueErr > region_value_err_
Set holds data of valid / invalid element values on all regions.
Definition: field_fe.hh:365
void local_to_ghost_data_scatter_begin()
Call begin scatter functions (local to ghost) on data vector.
Definition: field_fe.cc:783
TYPEDEF_ERR_INFO(EI_File, std::string)
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:635
virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData &init_data)
Definition: field_fe.cc:285
virtual ~FieldFE()
Destructor.
Definition: field_fe.cc:833
Field< spacedim, Value >::FactoryBase FieldFactoryBaseType
Definition: field_fe.hh:56
DataInterpolation
Definition: field_fe.hh:62
@ interp_p0
P0 interpolation (with the use of calculation of intersections)
Definition: field_fe.hh:66
@ equivalent_msh
equivalent mesh (default value)
Definition: field_fe.hh:64
@ gauss_p0
P0 interpolation (with the use of Gaussian distribution)
Definition: field_fe.hh:65
@ identic_msh
identical mesh
Definition: field_fe.hh:63
void local_to_ghost_data_scatter_end()
Call end scatter functions (local to ghost) on data vector.
Definition: field_fe.cc:790
static const Input::Type::Record & get_input_type()
Implementation.
Definition: field_fe.cc:60
unsigned int compute_fe_quadrature(std::vector< arma::vec::fixed< 3 >> &q_points, std::vector< double > &q_weights, const ElementAccessor< spacedim > &elm, unsigned int order=3)
Definition: field_fe.hh:307
std::vector< FEValues< spacedim > > fe_values_
List of FEValues objects of dimensions 0,1,2,3 used for value calculation.
Definition: field_fe.hh:355
void interpolate_gauss()
Interpolate data (use Gaussian distribution) over all elements of target mesh.
Definition: field_fe.cc:436
TYPEDEF_ERR_INFO(EI_Field, std::string)
Declaration of exception.
void make_dof_handler(const MeshBase *mesh)
Create DofHandler object.
Definition: field_fe.cc:345
void cache_reinit(AssemblyInternals &asm_internals) override
Definition: field_fe.cc:262
Quadrature init_quad(std::shared_ptr< EvalPoints > eval_points)
Initialize FEValues object of given dimension.
Definition: field_fe.cc:275
OutputTime::DiscreteSpace discretization_
Specify section where to find the field data in input mesh file.
Definition: field_fe.hh:337
std::shared_ptr< DOFHandlerMultiDim > get_dofhandler() const
Definition: field_fe.hh:178
TYPEDEF_ERR_INFO(EI_Region, std::string)
FilePath reader_file_
mesh reader file
Definition: field_fe.hh:331
VectorMPI & vec()
Definition: field_fe.hh:182
std::array< FEItem, 4 > fe_item_
Holds specific data of field evaluation over all dimensions.
Definition: field_fe.hh:361
TYPEDEF_ERR_INFO(EI_ElemIdx, unsigned int)
bool set_time(const TimeStep &time) override
Definition: field_fe.cc:385
static const int registrar
Registrar of class to factory.
Definition: field_fe.hh:371
void interpolate_intersection()
Interpolate data (use intersection library) over all elements of target mesh.
Definition: field_fe.cc:532
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:175
double default_value_
Default value of element if not set in mesh data file.
Definition: field_fe.hh:346
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")
const VectorMPI & vec() const
Definition: field_fe.hh:186
void cache_update(FieldValueCache< typename Value::element_type > &data_cache, ElementCacheMap &cache_map, unsigned int region_patch_idx) override
Definition: field_fe.cc:215
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:352
DataInterpolation interpolation_
Specify type of FE data interpolation.
Definition: field_fe.hh:340
Input::Record in_rec_
Accessor to Input::Record.
Definition: field_fe.hh:349
FieldFlag::Flags flags_
Field flags.
Definition: field_fe.hh:343
void native_data_to_cache(ElementDataCache< double > &output_data_cache)
Definition: field_fe.cc:757
std::string field_name_
field name read from input
Definition: field_fe.hh:334
FieldFE(unsigned int n_comp=0)
Definition: field_fe.cc:134
std::shared_ptr< DOFHandlerMultiDim > dh_
DOF handler object.
Definition: field_fe.hh:326
void fill_fe_item()
Definition: field_fe.hh:284
std::shared_ptr< EquivalentMeshMap > source_target_mesh_elm_map_
Maps element indices from computational mesh to the source (data).
Definition: field_fe.hh:358
void fill_fe_system_data(unsigned int block_index)
Definition: field_fe.hh:274
static const unsigned int undef_uint
Definition: field_fe.hh:82
static const Input::Type::Selection & get_disc_selection_input_type()
Definition: field_fe.cc:93
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:368
void set_mesh(const Mesh *mesh) override
Definition: field_fe.cc:309
std::shared_ptr< FieldBaseType > FieldBasePtr
Definition: field.hh:97
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
Record type proxy class.
Definition: type_record.hh:182
Template for classes storing finite set of named values.
static RealPoint project_unit_to_real(const BaryPoint &point, const ElementMap &map)
Definition: mapping_p1.cc:85
static ElementMap element_map(ElementAccessor< 3 > elm)
Definition: mapping_p1.cc:48
Base class for Mesh and BCMesh.
Definition: mesh.h:96
Definition: mesh.h:362
Symmetric Gauss-Legendre quadrature formulae on simplices.
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:86
Armor::ArmaVec< double, point_dim > point(unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:151
double weight(unsigned int i) const
Returns the ith weight.
Definition: quadrature.hh:103
Representation of one time step..
static VectorMPI sequential(unsigned int size)
Definition: vector_mpi.cc:44
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
Definitions of basic Lagrangean finite elements with polynomial shape functions.
Class FESystem for compound finite elements.
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:377
Abstract class for description of finite elements.
static constexpr bool value
Definition: json.hpp:87
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell.
typename arma::Mat< Type >::template fixed< nr, nc > ArmaMat
Definition: armor.hh:502
ArmaMat< double, N, M > mat
Definition: armor.hh:936
ArmaVec< double, N > vec
Definition: armor.hh:933
Implementation of range helper class.
Holds common data shared between GenericAssemblz and Assembly<dim> classes.
Definition: mixed.hh:25
Helper struct stores data for initizalize descentants of FieldAlgorithmBase.