Flow123d  JS_before_hm-915-gc632be9
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  /**
67  * Default constructor, optionally we need number of components @p n_comp in the case of Vector valued fields.
68  */
69  FieldFE(unsigned int n_comp=0);
70 
71  /**
72  * Return Record for initialization of FieldFE that is derived from Abstract
73  */
74  static const Input::Type::Record & get_input_type();
75 
76  /**
77  * Return Input selection for discretization type (determines the section of VTK file).
78  */
80 
81  /**
82  * Return Input selection that allow to set interpolation of input data.
83  */
85 
86  /**
87  * Setter for the finite element data.
88  * @param dh Dof handler.
89  * @param data Vector of dof values, optional (create own vector according to dofhandler).
90  * @param component_index Index of component (for vector_view/tensor_view)
91  * @return Data vector of dof values.
92  */
93  VectorMPI set_fe_data(std::shared_ptr<DOFHandlerMultiDim> dh,
94  unsigned int component_index = 0, VectorMPI dof_values = VectorMPI::sequential(0));
95 
96  /**
97  * Returns one value in one given point. ResultType can be used to avoid some costly calculation if the result is trivial.
98  */
99  virtual typename Value::return_type const &value(const Point &p, const ElementAccessor<spacedim> &elm);
100 
101  /**
102  * Returns std::vector of scalar values in several points at once.
103  */
104  virtual void value_list (const Armor::array &point_list, const ElementAccessor<spacedim> &elm,
106 
107  /**
108  * Overload @p FieldAlgorithmBase::cache_update
109  */
111  ElementCacheMap &cache_map, unsigned int region_idx) override;
112 
113  /**
114  * Initialization from the input interface.
115  */
116  virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData& init_data);
117 
118 
119  /**
120  * Update time and possibly update data from GMSH file.
121  */
122  bool set_time(const TimeStep &time) override;
123 
124 
125  /**
126  * Set target mesh.
127  */
128  void set_mesh(const Mesh *mesh, bool boundary_domain) override;
129 
130 
131  /**
132  * Copy data vector to given output ElementDataCache
133  */
134  void native_data_to_cache(ElementDataCache<double> &output_data_cache);
135 
136 
137  /**
138  * Return size of vector of data stored in Field
139  */
140  unsigned int data_size() const;
141 
142  inline std::shared_ptr<DOFHandlerMultiDim> get_dofhandler() const {
143  return dh_;
144  }
145 
146  inline VectorMPI get_data_vec() const {
147  return data_vec_;
148  }
149 
150  /// Call begin scatter functions (local to ghost) on data vector
152 
153  /// Call end scatter functions (local to ghost) on data vector
155 
156  /// Destructor.
157  virtual ~FieldFE();
158 
159 private:
160  /// Create DofHandler object
161  void make_dof_handler(const Mesh *mesh);
162 
163  /// Interpolate data (use Gaussian distribution) over all elements of target mesh.
165 
166  /// Interpolate data (use intersection library) over all elements of target mesh.
168 
169  /// Calculate native data over all elements of target mesh.
171 
172  /// Calculate elementwise data over all elements of target mesh.
174 
175  /**
176  * Fill data to boundary_dofs_ vector.
177  *
178  * TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
179  */
180  void fill_boundary_dofs();
181 
182  /// Initialize FEValues object of given dimension.
183  template <unsigned int dim>
184  Quadrature init_quad(std::shared_ptr<EvalPoints> eval_points);
185 
187  unsigned int i_dof, unsigned int i_qp, unsigned int comp_index);
188 
189 
190  /// DOF handler object
191  std::shared_ptr<DOFHandlerMultiDim> dh_;
192  /// Store data of Field
194 
195  /// Value handler that allows get value of 0D elements.
197  /// Value handler that allows get value of 1D elements.
199  /// Value handler that allows get value of 2D elements.
201  /// Value handler that allows get value of 3D elements.
203 
204  /**
205  * Used in DOFHandler::distribute_dofs method. Represents 0D element.
206  *
207  * For correct functionality must be created proper descendant of FiniteElement class.
208  */
210 
211  /// mesh reader file
213 
214  /// field name read from input
215  std::string field_name_;
216 
217  /// Specify section where to find the field data in input mesh file
219 
220  /// Specify type of FE data interpolation
222 
223  /// Field flags.
225 
226  /// Default value of element if not set in mesh data file
228 
229  /// Accessor to Input::Record
231 
232  /// Is set in set_mesh method. Value true means, that we accept only boundary element accessors in the @p value method.
234 
235  /**
236  * Hold dofs of boundary elements.
237  *
238  * TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
239  */
240  std::shared_ptr< std::vector<IntIdx> > boundary_dofs_;
241 
242  /// List of FEValues objects of dimensions 0,1,2,3 used for value calculation
244 
245  /// Registrar of class to factory
246  static const int registrar;
247 };
248 
249 
250 /** Create FieldFE from dhf handler */
251 template <int spacedim, class Value>
252 std::shared_ptr<FieldFE<spacedim, Value> > create_field_fe(std::shared_ptr<DOFHandlerMultiDim> dh,
253  unsigned int comp = 0,
254  VectorMPI *vec = nullptr)
255 {
256  // Construct FieldFE
257  std::shared_ptr< FieldFE<spacedim, Value> > field_ptr = std::make_shared< FieldFE<spacedim, Value> >();
258  if (vec == nullptr)
259  field_ptr->set_fe_data( dh, comp, dh->create_vector() );
260  else
261  field_ptr->set_fe_data( dh, comp, *vec );
262 
263  return field_ptr;
264 }
265 
266 
267 /** Create FieldFE with parallel VectorMPI from finite element */
268 template <int spacedim, class Value>
269 std::shared_ptr<FieldFE<spacedim, Value> > create_field_fe(Mesh & mesh, const MixedPtr<FiniteElement> &fe)
270 {
271  // Prepare DOF handler
272  std::shared_ptr<DOFHandlerMultiDim> dh_par = std::make_shared<DOFHandlerMultiDim>(mesh);
273  std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( &mesh, fe);
274  dh_par->distribute_dofs(ds);
275 
276  return create_field_fe<spacedim,Value>(dh_par);
277 }
278 
279 
280 #endif /* FIELD_FE_HH_ */
virtual ~FieldFE()
Destructor.
Definition: field_fe.cc:807
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:476
FEValueHandler< 0, spacedim, Value > value_handler0_
Value handler that allows get value of 0D elements.
Definition: field_fe.hh:196
std::shared_ptr< std::vector< T > > ComponentDataPtr
P0 interpolation (with the use of calculation of intersections)
Definition: field_fe.hh:63
ArmaVec< double, N > vec
Definition: armor.hh:861
MixedPtr< FiniteElement > fe_
Definition: field_fe.hh:209
std::vector< FEValues< spacedim > > fe_values_
List of FEValues objects of dimensions 0,1,2,3 used for value calculation.
Definition: field_fe.hh:243
std::string field_name_
field name read from input
Definition: field_fe.hh:215
VectorMPI get_data_vec() const
Definition: field_fe.hh:146
FEValueHandler< 2, spacedim, Value > value_handler2_
Value handler that allows get value of 2D elements.
Definition: field_fe.hh:200
VectorMPI data_vec_
Store data of Field.
Definition: field_fe.hh:193
FilePath reader_file_
mesh reader file
Definition: field_fe.hh:212
Directing class of FieldValueCache.
Definition: mesh.h:78
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
Class FESystem for compound finite elements.
static const int registrar
Registrar of class to factory.
Definition: field_fe.hh:246
void set_mesh(const Mesh *mesh, bool boundary_domain) override
Definition: field_fe.cc:309
double default_value_
Default value of element if not set in mesh data file.
Definition: field_fe.hh:227
FEValueHandler< 3, spacedim, Value > value_handler3_
Value handler that allows get value of 3D elements.
Definition: field_fe.hh:202
bool set_time(const TimeStep &time) override
Definition: field_fe.cc:424
static const Input::Type::Selection & get_interp_selection_input_type()
Definition: field_fe.cc:103
void fill_boundary_dofs()
Definition: field_fe.cc:349
FieldAlgorithmBase< spacedim, Value > FactoryBaseType
Definition: field_fe.hh:53
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
unsigned int data_size() const
Definition: field_fe.cc:771
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:233
OutputTime::DiscreteSpace discretization_
Specify section where to find the field data in input mesh file.
Definition: field_fe.hh:218
std::shared_ptr< DOFHandlerMultiDim > dh_
DOF handler object.
Definition: field_fe.hh:191
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:276
Armor::ArmaMat< typename Value::element_type, Value::NRows_, Value::NCols_ > handle_fe_shape(unsigned int dim, unsigned int i_dof, unsigned int i_qp, unsigned int comp_index)
Definition: field_fe.cc:792
DataInterpolation interpolation_
Specify type of FE data interpolation.
Definition: field_fe.hh:221
Definitions of basic Lagrangean finite elements with polynomial shape functions.
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
void interpolate_intersection(ElementDataCache< double >::ComponentDataPtr data_vec)
Interpolate data (use intersection library) over all elements of target mesh.
Definition: field_fe.cc:571
std::shared_ptr< std::vector< IntIdx > > boundary_dofs_
Definition: field_fe.hh:240
void cache_update(FieldValueCache< typename Value::element_type > &data_cache, ElementCacheMap &cache_map, unsigned int region_idx) override
Definition: field_fe.cc:231
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Definition: field_fe.cc:180
void local_to_ghost_data_scatter_begin()
Call begin scatter functions (local to ghost) on data vector.
Definition: field_fe.cc:778
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:785
VectorMPI set_fe_data(std::shared_ptr< DOFHandlerMultiDim > dh, unsigned int component_index=0, VectorMPI dof_values=VectorMPI::sequential(0))
Definition: field_fe.cc:141
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:198
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:142
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:230
FieldFE(unsigned int n_comp=0)
Definition: field_fe.cc:132
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:204
virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData &init_data)
Definition: field_fe.cc:286
Record type proxy class.
Definition: type_record.hh:182
Abstract class for description of finite elements.
equivalent mesh (default value)
Definition: field_fe.hh:61
void calculate_elementwise_values(ElementDataCache< double >::ComponentDataPtr data_cache)
Calculate elementwise data over all elements of target mesh.
Definition: field_fe.cc:711
void calculate_native_values(ElementDataCache< double >::ComponentDataPtr data_cache)
Calculate native data over all elements of target mesh.
Definition: field_fe.cc:682
FieldFlag::Flags flags_
Field flags.
Definition: field_fe.hh:224
std::shared_ptr< FieldFE< spacedim, Value > > create_field_fe(std::shared_ptr< DOFHandlerMultiDim > dh, unsigned int comp=0, VectorMPI *vec=nullptr)
Definition: field_fe.hh:252
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:376
Implementation of range helper class.
void native_data_to_cache(ElementDataCache< double > &output_data_cache)
Definition: field_fe.cc:752
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:60
unsigned int n_comp() const