Flow123d  JS_before_hm-1626-gde32303
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 "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 "input/factory.hh"
39 
40 #include <memory>
41 
42 
43 
44 /**
45  * Class representing fields given by finite element approximation.
46  *
47  */
48 template <int spacedim, class Value>
49 class FieldFE : public FieldAlgorithmBase<spacedim, Value>
50 {
51 public:
54 
55  /**
56  * Possible interpolations of input data.
57  */
59  {
60  identic_msh, //!< identical mesh
61  equivalent_msh, //!< equivalent mesh (default value)
62  gauss_p0, //!< P0 interpolation (with the use of Gaussian distribution)
63  interp_p0 //!< P0 interpolation (with the use of calculation of intersections)
64  };
65 
66  TYPEDEF_ERR_INFO( EI_ElemIdx, unsigned int);
67  DECLARE_EXCEPTION( ExcInvalidElemeDim, << "Dimension of element in target mesh must be 0, 1 or 2! elm.idx() = " << EI_ElemIdx::val << ".\n" );
68 
69  static const unsigned int undef_uint = -1;
70 
71  /**
72  * Default constructor, optionally we need number of components @p n_comp in the case of Vector valued fields.
73  */
74  FieldFE(unsigned int n_comp=0);
75 
76  /**
77  * Return Record for initialization of FieldFE that is derived from Abstract
78  */
79  static const Input::Type::Record & get_input_type();
80 
81  /**
82  * Return Input selection for discretization type (determines the section of VTK file).
83  */
85 
86  /**
87  * Return Input selection that allow to set interpolation of input data.
88  */
90 
91  /**
92  * Setter for the finite element data.
93  * @param dh Dof handler.
94  * @param dof_values Vector of dof values, optional (create own vector according to dofhandler).
95  * @param block_index Index of block (in FESystem) or '-1' for FieldFE under simple DOF handler.
96  * @return Data vector of dof values.
97  */
98  VectorMPI set_fe_data(std::shared_ptr<DOFHandlerMultiDim> dh, VectorMPI dof_values = VectorMPI::sequential(0),
99  unsigned int block_index = FieldFE<spacedim, Value>::undef_uint);
100 
101  /**
102  * Returns one value in one given point. ResultType can be used to avoid some costly calculation if the result is trivial.
103  */
104  virtual typename Value::return_type const &value(const Point &p, const ElementAccessor<spacedim> &elm);
105 
106  /**
107  * Returns std::vector of scalar values in several points at once.
108  */
109  virtual void value_list (const Armor::array &point_list, const ElementAccessor<spacedim> &elm,
111 
112  /**
113  * Overload @p FieldAlgorithmBase::cache_update
114  */
116  ElementCacheMap &cache_map, unsigned int region_patch_idx) override;
117 
118  /**
119  * Overload @p FieldAlgorithmBase::cache_reinit
120  *
121  * Reinit fe_values_ data member.
122  */
123  void cache_reinit(const ElementCacheMap &cache_map) override;
124 
125  /**
126  * Initialization from the input interface.
127  */
128  virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData& init_data);
129 
130 
131  /**
132  * Update time and possibly update data from GMSH file.
133  */
134  bool set_time(const TimeStep &time) override;
135 
136 
137  /**
138  * Set target mesh.
139  */
140  void set_mesh(const Mesh *mesh, bool boundary_domain) override;
141 
142 
143  /**
144  * Copy data vector to given output ElementDataCache
145  */
146  void native_data_to_cache(ElementDataCache<double> &output_data_cache);
147 
148 
149  /**
150  * Return size of vector of data stored in Field
151  */
152  unsigned int data_size() const;
153 
154  inline std::shared_ptr<DOFHandlerMultiDim> get_dofhandler() const {
155  return dh_;
156  }
157 
158  inline VectorMPI& vec() {
159  return data_vec_;
160  }
161 
162  /// Call begin scatter functions (local to ghost) on data vector
164 
165  /// Call end scatter functions (local to ghost) on data vector
167 
168  /// Destructor.
169  virtual ~FieldFE();
170 
171 private:
172  /**
173  * Helper class holds specific data of field evaluation.
174  */
175  class FEItem {
176  public:
177  unsigned int comp_index_;
178  unsigned int range_begin_;
179  unsigned int range_end_;
180  };
181 
182  /// Create DofHandler object
183  void make_dof_handler(const Mesh *mesh);
184 
185  /// Interpolate data (use Gaussian distribution) over all elements of target mesh.
187 
188  /// Interpolate data (use intersection library) over all elements of target mesh.
190 
191  /// Calculate native data over all elements of target mesh.
193 
194  /// Calculate data of identict_mesh interpolation on input data over all elements of target mesh.
196 
197  /// Calculate data of equivalent_mesh interpolation on input over all elements of target mesh.
199 
200  /**
201  * Fill data to boundary_dofs_ vector.
202  *
203  * TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
204  */
205  void fill_boundary_dofs();
206 
207  /// Initialize FEValues object of given dimension.
208  template <unsigned int dim>
209  Quadrature init_quad(std::shared_ptr<EvalPoints> eval_points);
210 
212  unsigned int i_dof, unsigned int i_qp)
213  {
215  for (unsigned int c=0; c<Value::NRows_*Value::NCols_; ++c)
216  v(c/spacedim,c%spacedim) = fe_values_[dim].shape_value_component(i_dof, i_qp, c);
217  if (Value::NRows_ == Value::NCols_)
218  return v;
219  else
220  return v.t();
221  }
222 
223  template<unsigned int dim>
224  void fill_fe_system_data(unsigned int block_index) {
225  auto fe_system_ptr = std::dynamic_pointer_cast<FESystem<dim>>( dh_->ds()->fe()[Dim<dim>{}] );
226  ASSERT_DBG(fe_system_ptr != nullptr).error("Wrong type, must be FESystem!\n");
227  this->fe_item_[dim].comp_index_ = fe_system_ptr->function_space()->dof_indices()[block_index].component_offset;
228  this->fe_item_[dim].range_begin_ = fe_system_ptr->fe_dofs(block_index)[0];
229  this->fe_item_[dim].range_end_ = this->fe_item_[dim].range_begin_ + fe_system_ptr->fe()[block_index]->n_dofs();
230  }
231 
232 
233  template<unsigned int dim>
234  void fill_fe_item() {
235  this->fe_item_[dim].comp_index_ = 0;
236  this->fe_item_[dim].range_begin_ = 0;
237  this->fe_item_[dim].range_end_ = dh_->ds()->fe()[Dim<dim>{}]->n_dofs();
238  }
239 
240 
241  /// DOF handler object
242  std::shared_ptr<DOFHandlerMultiDim> dh_;
243  /// Store data of Field
245 
246  /// Value handler that allows get value of 0D elements.
248  /// Value handler that allows get value of 1D elements.
250  /// Value handler that allows get value of 2D elements.
252  /// Value handler that allows get value of 3D elements.
254 
255  /// mesh reader file
257 
258  /// field name read from input
259  std::string field_name_;
260 
261  /// Specify section where to find the field data in input mesh file
263 
264  /// Specify type of FE data interpolation
266 
267  /// Field flags.
269 
270  /// Default value of element if not set in mesh data file
272 
273  /// Accessor to Input::Record
275 
276  /// Is set in set_mesh method. Value true means, that we accept only boundary element accessors in the @p value method.
278 
279  /**
280  * Hold dofs of boundary elements.
281  *
282  * TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
283  */
284  std::shared_ptr< std::vector<IntIdx> > boundary_dofs_;
285 
286  /// List of FEValues objects of dimensions 0,1,2,3 used for value calculation
288 
289  /// Maps element indices between source (data) and target (computational) mesh if data interpolation is set to equivalent_msh
290  std::shared_ptr<std::vector<LongIdx>> source_target_mesh_elm_map_;
291 
292  /// Holds specific data of field evaluation over all dimensions.
293  std::array<FEItem, 4> fe_item_;
295 
296  /// Registrar of class to factory
297  static const int registrar;
298 };
299 
300 
301 /** Create FieldFE from dof handler */
302 template <int spacedim, class Value>
303 std::shared_ptr<FieldFE<spacedim, Value> > create_field_fe(std::shared_ptr<DOFHandlerMultiDim> dh,
304  VectorMPI *vec = nullptr,
305  unsigned int block_index = FieldFE<spacedim, Value>::undef_uint)
306 {
307  // Construct FieldFE
308  std::shared_ptr< FieldFE<spacedim, Value> > field_ptr = std::make_shared< FieldFE<spacedim, Value> >();
309  if (vec == nullptr)
310  field_ptr->set_fe_data( dh, dh->create_vector(), block_index );
311  else
312  field_ptr->set_fe_data( dh, *vec, block_index );
313 
314  return field_ptr;
315 }
316 
317 
318 /** Create FieldFE with parallel VectorMPI from finite element */
319 template <int spacedim, class Value>
320 std::shared_ptr<FieldFE<spacedim, Value> > create_field_fe(Mesh & mesh, const MixedPtr<FiniteElement> &fe)
321 {
322  // Prepare DOF handler
323  std::shared_ptr<DOFHandlerMultiDim> dh_par = std::make_shared<DOFHandlerMultiDim>(mesh);
324  std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( &mesh, fe);
325  dh_par->distribute_dofs(ds);
326 
327  return create_field_fe<spacedim,Value>(dh_par);
328 }
329 
330 
331 #endif /* FIELD_FE_HH_ */
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
virtual ~FieldFE()
Destructor.
Definition: field_fe.cc:915
Definition: mixed.hh:25
unsigned int range_begin_
Definition: field_fe.hh:178
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
void interpolate_gauss(ElementDataCache< double >::ComponentDataPtr data_vec)
Interpolate data (use Gaussian distribution) over all elements of target mesh.
Definition: field_fe.cc:521
FEValueHandler< 0, spacedim, Value > value_handler0_
Value handler that allows get value of 0D elements.
Definition: field_fe.hh:247
std::shared_ptr< std::vector< T > > ComponentDataPtr
P0 interpolation (with the use of calculation of intersections)
Definition: field_fe.hh:63
MixedPtr< FiniteElement > fe_
Definition: field_fe.hh:294
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
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
void fill_fe_item()
Definition: field_fe.hh:234
std::string field_name_
field name read from input
Definition: field_fe.hh:259
static const unsigned int undef_uint
Definition: field_fe.hh:69
FEValueHandler< 2, spacedim, Value > value_handler2_
Value handler that allows get value of 2D elements.
Definition: field_fe.hh:251
VectorMPI data_vec_
Store data of Field.
Definition: field_fe.hh:244
FilePath reader_file_
mesh reader file
Definition: field_fe.hh:256
Directing class of FieldValueCache.
Definition: mesh.h:77
identical mesh
Definition: field_fe.hh:60
Helper struct stores data for initizalize descentants of FieldAlgorithmBase.
FieldAlgorithmBase< spacedim, Value >::Point Point
Definition: field_fe.hh:52
unsigned int range_end_
Definition: field_fe.hh:179
Class FESystem for compound finite elements.
static const int registrar
Registrar of class to factory.
Definition: field_fe.hh:297
void set_mesh(const Mesh *mesh, bool boundary_domain) override
Definition: field_fe.cc:341
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
void fill_boundary_dofs()
Definition: field_fe.cc:373
FieldAlgorithmBase< spacedim, Value > FactoryBaseType
Definition: field_fe.hh:53
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
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
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:303
unsigned int data_size() const
Definition: field_fe.cc:879
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
TYPEDEF_ERR_INFO(EI_ElemIdx, unsigned int)
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
static const Input::Type::Record & get_input_type()
Implementation.
Definition: field_fe.cc:61
Quadrature init_quad(std::shared_ptr< EvalPoints > eval_points)
Initialize FEValues object of given dimension.
Definition: field_fe.cc:308
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.
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
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
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Definition: field_fe.cc:207
void local_to_ghost_data_scatter_begin()
Call begin scatter functions (local to ghost) on data vector.
Definition: field_fe.cc:886
Space< spacedim >::Point Point
void local_to_ghost_data_scatter_end()
Call end scatter functions (local to ghost) on data vector.
Definition: field_fe.cc:893
typename arma::Mat< Type >::template fixed< nr, nc > ArmaMat
Definition: armor.hh:502
FEValueHandler< 1, spacedim, Value > value_handler1_
Value handler that allows get value of 1D elements.
Definition: field_fe.hh:249
DECLARE_EXCEPTION(ExcInvalidElemeDim,<< "Dimension of element in target mesh must be 0, 1 or 2! elm.idx() = "<< EI_ElemIdx::val<< ".\n")
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
std::shared_ptr< DOFHandlerMultiDim > get_dofhandler() const
Definition: field_fe.hh:154
static const Input::Type::Selection & get_disc_selection_input_type()
Definition: field_fe.cc:92
Input::Record in_rec_
Accessor to Input::Record.
Definition: field_fe.hh:274
FieldFE(unsigned int n_comp=0)
Definition: field_fe.cc:133
#define ASSERT_DBG(expr)
void cache_update(FieldValueCache< typename Value::element_type > &data_cache, ElementCacheMap &cache_map, unsigned int region_patch_idx) override
Definition: field_fe.cc:258
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
Record type proxy class.
Definition: type_record.hh:182
void fill_fe_system_data(unsigned int block_index)
Definition: field_fe.hh:224
Abstract class for description of finite elements.
void cache_reinit(const ElementCacheMap &cache_map) override
Definition: field_fe.cc:295
equivalent mesh (default value)
Definition: field_fe.hh:61
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
unsigned int comp_index_
Definition: field_fe.hh:177
Representation of one time step..
Template for classes storing finite set of named values.
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
static VectorMPI sequential(unsigned int size)
Definition: vector_mpi.cc:44
P0 interpolation (with the use of Gaussian distribution)
Definition: field_fe.hh:62
Definition: field.hh:63
unsigned int n_comp() const
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