Flow123d  JS_before_hm-1575-ga41e096
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  static const unsigned int undef_uint = -1;
67 
68  /**
69  * Default constructor, optionally we need number of components @p n_comp in the case of Vector valued fields.
70  */
71  FieldFE(unsigned int n_comp=0);
72 
73  /**
74  * Return Record for initialization of FieldFE that is derived from Abstract
75  */
76  static const Input::Type::Record & get_input_type();
77 
78  /**
79  * Return Input selection for discretization type (determines the section of VTK file).
80  */
82 
83  /**
84  * Return Input selection that allow to set interpolation of input data.
85  */
87 
88  /**
89  * Setter for the finite element data.
90  * @param dh Dof handler.
91  * @param dof_values Vector of dof values, optional (create own vector according to dofhandler).
92  * @param block_index Index of block (in FESystem) or '-1' for FieldFE under simple DOF handler.
93  * @return Data vector of dof values.
94  */
95  VectorMPI set_fe_data(std::shared_ptr<DOFHandlerMultiDim> dh, VectorMPI dof_values = VectorMPI::sequential(0),
96  unsigned int block_index = FieldFE<spacedim, Value>::undef_uint);
97 
98  /**
99  * Returns one value in one given point. ResultType can be used to avoid some costly calculation if the result is trivial.
100  */
101  virtual typename Value::return_type const &value(const Point &p, const ElementAccessor<spacedim> &elm);
102 
103  /**
104  * Returns std::vector of scalar values in several points at once.
105  */
106  virtual void value_list (const Armor::array &point_list, const ElementAccessor<spacedim> &elm,
108 
109  /**
110  * Overload @p FieldAlgorithmBase::cache_update
111  */
113  ElementCacheMap &cache_map, unsigned int region_patch_idx) override;
114 
115  /**
116  * Overload @p FieldAlgorithmBase::cache_reinit
117  *
118  * Reinit fe_values_ data member.
119  */
120  void cache_reinit(const ElementCacheMap &cache_map) override;
121 
122  /**
123  * Initialization from the input interface.
124  */
125  virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData& init_data);
126 
127 
128  /**
129  * Update time and possibly update data from GMSH file.
130  */
131  bool set_time(const TimeStep &time) override;
132 
133 
134  /**
135  * Set target mesh.
136  */
137  void set_mesh(const Mesh *mesh, bool boundary_domain) override;
138 
139 
140  /**
141  * Copy data vector to given output ElementDataCache
142  */
143  void native_data_to_cache(ElementDataCache<double> &output_data_cache);
144 
145 
146  /**
147  * Return size of vector of data stored in Field
148  */
149  unsigned int data_size() const;
150 
151  inline std::shared_ptr<DOFHandlerMultiDim> get_dofhandler() const {
152  return dh_;
153  }
154 
155  inline VectorMPI& vec() {
156  return data_vec_;
157  }
158 
159  /// Call begin scatter functions (local to ghost) on data vector
161 
162  /// Call end scatter functions (local to ghost) on data vector
164 
165  /// Destructor.
166  virtual ~FieldFE();
167 
168 private:
169  /**
170  * Helper class holds specific data of field evaluation.
171  */
172  class FEItem {
173  public:
174  unsigned int comp_index_;
175  unsigned int range_begin_;
176  unsigned int range_end_;
177  };
178 
179  /// Create DofHandler object
180  void make_dof_handler(const Mesh *mesh);
181 
182  /// Interpolate data (use Gaussian distribution) over all elements of target mesh.
184 
185  /// Interpolate data (use intersection library) over all elements of target mesh.
187 
188  /// Calculate native data over all elements of target mesh.
190 
191  /// Calculate data of identict_mesh interpolation on input data over all elements of target mesh.
193 
194  /// Calculate data of equivalent_mesh interpolation on input over all elements of target mesh.
196 
197  /**
198  * Fill data to boundary_dofs_ vector.
199  *
200  * TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
201  */
202  void fill_boundary_dofs();
203 
204  /// Initialize FEValues object of given dimension.
205  template <unsigned int dim>
206  Quadrature init_quad(std::shared_ptr<EvalPoints> eval_points);
207 
209  unsigned int i_dof, unsigned int i_qp)
210  {
212  for (unsigned int c=0; c<Value::NRows_*Value::NCols_; ++c)
213  v(c/spacedim,c%spacedim) = fe_values_[dim].shape_value_component(i_dof, i_qp, c);
214  if (Value::NRows_ == Value::NCols_)
215  return v;
216  else
217  return v.t();
218  }
219 
220  template<unsigned int dim>
221  void fill_fe_system_data(unsigned int block_index) {
222  auto fe_system_ptr = std::dynamic_pointer_cast<FESystem<dim>>( dh_->ds()->fe()[Dim<dim>{}] );
223  ASSERT_DBG(fe_system_ptr != nullptr).error("Wrong type, must be FESystem!\n");
224  this->fe_item_[dim].comp_index_ = fe_system_ptr->function_space()->dof_indices()[block_index].component_offset;
225  this->fe_item_[dim].range_begin_ = fe_system_ptr->fe_dofs(block_index)[0];
226  this->fe_item_[dim].range_end_ = this->fe_item_[dim].range_begin_ + fe_system_ptr->fe()[block_index]->n_dofs();
227  }
228 
229 
230  template<unsigned int dim>
231  void fill_fe_item() {
232  this->fe_item_[dim].comp_index_ = 0;
233  this->fe_item_[dim].range_begin_ = 0;
234  this->fe_item_[dim].range_end_ = dh_->ds()->fe()[Dim<dim>{}]->n_dofs();
235  }
236 
237 
238  /// DOF handler object
239  std::shared_ptr<DOFHandlerMultiDim> dh_;
240  /// Store data of Field
242 
243  /// Value handler that allows get value of 0D elements.
245  /// Value handler that allows get value of 1D elements.
247  /// Value handler that allows get value of 2D elements.
249  /// Value handler that allows get value of 3D elements.
251 
252  /// mesh reader file
254 
255  /// field name read from input
256  std::string field_name_;
257 
258  /// Specify section where to find the field data in input mesh file
260 
261  /// Specify type of FE data interpolation
263 
264  /// Field flags.
266 
267  /// Default value of element if not set in mesh data file
269 
270  /// Accessor to Input::Record
272 
273  /// Is set in set_mesh method. Value true means, that we accept only boundary element accessors in the @p value method.
275 
276  /**
277  * Hold dofs of boundary elements.
278  *
279  * TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
280  */
281  std::shared_ptr< std::vector<IntIdx> > boundary_dofs_;
282 
283  /// List of FEValues objects of dimensions 0,1,2,3 used for value calculation
285 
286  /// Maps element indices between source (data) and target (computational) mesh if data interpolation is set to equivalent_msh
287  std::shared_ptr<std::vector<LongIdx>> source_target_mesh_elm_map_;
288 
289  /// Holds specific data of field evaluation over all dimensions.
290  std::array<FEItem, 4> fe_item_;
292 
293  /// Registrar of class to factory
294  static const int registrar;
295 };
296 
297 
298 /** Create FieldFE from dof handler */
299 template <int spacedim, class Value>
300 std::shared_ptr<FieldFE<spacedim, Value> > create_field_fe(std::shared_ptr<DOFHandlerMultiDim> dh,
301  VectorMPI *vec = nullptr,
302  unsigned int block_index = FieldFE<spacedim, Value>::undef_uint)
303 {
304  // Construct FieldFE
305  std::shared_ptr< FieldFE<spacedim, Value> > field_ptr = std::make_shared< FieldFE<spacedim, Value> >();
306  if (vec == nullptr)
307  field_ptr->set_fe_data( dh, dh->create_vector(), block_index );
308  else
309  field_ptr->set_fe_data( dh, *vec, block_index );
310 
311  return field_ptr;
312 }
313 
314 
315 /** Create FieldFE with parallel VectorMPI from finite element */
316 template <int spacedim, class Value>
317 std::shared_ptr<FieldFE<spacedim, Value> > create_field_fe(Mesh & mesh, const MixedPtr<FiniteElement> &fe)
318 {
319  // Prepare DOF handler
320  std::shared_ptr<DOFHandlerMultiDim> dh_par = std::make_shared<DOFHandlerMultiDim>(mesh);
321  std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( &mesh, fe);
322  dh_par->distribute_dofs(ds);
323 
324  return create_field_fe<spacedim,Value>(dh_par);
325 }
326 
327 
328 #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:208
virtual ~FieldFE()
Destructor.
Definition: field_fe.cc:914
Definition: mixed.hh:25
unsigned int range_begin_
Definition: field_fe.hh:175
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:520
FEValueHandler< 0, spacedim, Value > value_handler0_
Value handler that allows get value of 0D elements.
Definition: field_fe.hh:244
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:291
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:798
std::vector< FEValues< spacedim > > fe_values_
List of FEValues objects of dimensions 0,1,2,3 used for value calculation.
Definition: field_fe.hh:284
void fill_fe_item()
Definition: field_fe.hh:231
std::string field_name_
field name read from input
Definition: field_fe.hh:256
static const unsigned int undef_uint
Definition: field_fe.hh:66
FEValueHandler< 2, spacedim, Value > value_handler2_
Value handler that allows get value of 2D elements.
Definition: field_fe.hh:248
VectorMPI data_vec_
Store data of Field.
Definition: field_fe.hh:241
FilePath reader_file_
mesh reader file
Definition: field_fe.hh:253
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:176
Class FESystem for compound finite elements.
static const int registrar
Registrar of class to factory.
Definition: field_fe.hh:294
void set_mesh(const Mesh *mesh, bool boundary_domain) override
Definition: field_fe.cc:340
double default_value_
Default value of element if not set in mesh data file.
Definition: field_fe.hh:268
FEValueHandler< 3, spacedim, Value > value_handler3_
Value handler that allows get value of 3D elements.
Definition: field_fe.hh:250
bool set_time(const TimeStep &time) override
Definition: field_fe.cc:463
static const Input::Type::Selection & get_interp_selection_input_type()
Definition: field_fe.cc:103
void fill_boundary_dofs()
Definition: field_fe.cc:372
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:141
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:300
unsigned int data_size() const
Definition: field_fe.cc:878
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:274
OutputTime::DiscreteSpace discretization_
Specify section where to find the field data in input mesh file.
Definition: field_fe.hh:259
std::shared_ptr< DOFHandlerMultiDim > dh_
DOF handler object.
Definition: field_fe.hh:239
static const Input::Type::Record & get_input_type()
Implementation.
Definition: field_fe.cc:60
Quadrature init_quad(std::shared_ptr< EvalPoints > eval_points)
Initialize FEValues object of given dimension.
Definition: field_fe.cc:307
DataInterpolation interpolation_
Specify type of FE data interpolation.
Definition: field_fe.hh:262
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:290
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:615
std::shared_ptr< std::vector< IntIdx > > boundary_dofs_
Definition: field_fe.hh:281
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Definition: field_fe.cc:206
void local_to_ghost_data_scatter_begin()
Call begin scatter functions (local to ghost) on data vector.
Definition: field_fe.cc:885
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:892
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:246
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:151
static const Input::Type::Selection & get_disc_selection_input_type()
Definition: field_fe.cc:91
Input::Record in_rec_
Accessor to Input::Record.
Definition: field_fe.hh:271
FieldFE(unsigned int n_comp=0)
Definition: field_fe.cc:132
#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:257
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:287
VectorMPI & vec()
Definition: field_fe.hh:155
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:230
virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData &init_data)
Definition: field_fe.cc:317
Record type proxy class.
Definition: type_record.hh:182
void fill_fe_system_data(unsigned int block_index)
Definition: field_fe.hh:221
Abstract class for description of finite elements.
void cache_reinit(const ElementCacheMap &cache_map) override
Definition: field_fe.cc:294
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:724
FieldFlag::Flags flags_
Field flags.
Definition: field_fe.hh:265
unsigned int comp_index_
Definition: field_fe.hh:174
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:399
Implementation of range helper class.
void native_data_to_cache(ElementDataCache< double > &output_data_cache)
Definition: field_fe.cc:860
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:62
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:753