Flow123d  release_3.0.0-910-g56968e3
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"
25 #include "la/vector_mpi.hh"
26 #include "mesh/mesh.h"
27 #include "mesh/point.hh"
28 #include "mesh/bih_tree.hh"
29 #include "mesh/long_idx.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/mapping_p1.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:
55 
56  /**
57  * Possible interpolations of input data.
58  */
60  {
61  identic_msh, //!< identical mesh
62  equivalent_msh, //!< equivalent mesh (default value)
63  gauss_p0, //!< P0 interpolation (with the use of Gaussian distribution)
64  interp_p0 //!< P0 interpolation (with the use of calculation of intersections)
65  };
66 
67  /**
68  * Default constructor, optionally we need number of components @p n_comp in the case of Vector valued fields.
69  */
70  FieldFE(unsigned int n_comp=0);
71 
72  /**
73  * Return Record for initialization of FieldFE that is derived from Abstract
74  */
75  static const Input::Type::Record & get_input_type();
76 
77  /**
78  * Return Input selection for discretization type (determines the section of VTK file).
79  */
81 
82  /**
83  * Return Input selection that allow to set interpolation of input data.
84  */
86 
87  /**
88  * Setter for the finite element data.
89  * @param dh Dof handler.
90  * @param data Vector of dof values, optional (create own vector according to dofhandler).
91  * @param component_index Index of component (for vector_view/tensor_view)
92  * @return Data vector of dof values.
93  */
94  VectorMPI set_fe_data(std::shared_ptr<DOFHandlerMultiDim> dh,
95  unsigned int component_index = 0, VectorMPI dof_values = VectorMPI::sequential(0));
96 
97  /**
98  * Returns one value in one given point. ResultType can be used to avoid some costly calculation if the result is trivial.
99  */
100  virtual typename Value::return_type const &value(const Point &p, const ElementAccessor<spacedim> &elm);
101 
102  /**
103  * Returns std::vector of scalar values in several points at once.
104  */
105  virtual void value_list (const std::vector< Point > &point_list, const ElementAccessor<spacedim> &elm,
107 
108  /**
109  * Initialization from the input interface.
110  */
111  virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData& init_data);
112 
113 
114  /**
115  * Update time and possibly update data from GMSH file.
116  */
117  bool set_time(const TimeStep &time) override;
118 
119 
120  /**
121  * Set target mesh.
122  */
123  void set_mesh(const Mesh *mesh, bool boundary_domain) override;
124 
125 
126  /**
127  * Copy data vector to given output ElementDataCache
128  */
129  void native_data_to_cache(ElementDataCache<double> &output_data_cache);
130 
131 
132  /**
133  * Return size of vector of data stored in Field
134  */
135  unsigned int data_size() const;
136 
137  inline std::shared_ptr<DOFHandlerMultiDim> get_dofhandler() const {
138  return dh_;
139  }
140 
141  inline VectorMPI get_data_vec() const {
142  return data_vec_;
143  }
144 
145 
146  /// Destructor.
147  virtual ~FieldFE();
148 
149 private:
150  /// Create DofHandler object
151  void make_dof_handler(const Mesh *mesh);
152 
153  /// Interpolate data (use Gaussian distribution) over all elements of target mesh.
155 
156  /// Interpolate data (use intersection library) over all elements of target mesh.
158 
159  /// Calculate native data over all elements of target mesh.
161 
162  /**
163  * Fill data to boundary_dofs_ vector.
164  *
165  * TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
166  */
167  void fill_boundary_dofs();
168 
169 
170  /// DOF handler object
171  std::shared_ptr<DOFHandlerMultiDim> dh_;
172  /// Store data of Field
174  /// Array of indexes to data_vec_, used for get/set values
176 
177  /// Value handler that allows get value of 0D elements.
179  /// Value handler that allows get value of 1D elements.
181  /// Value handler that allows get value of 2D elements.
183  /// Value handler that allows get value of 3D elements.
185 
186  /**
187  * Used in DOFHandler::distribute_dofs method. Represents 0D element.
188  *
189  * For correct functionality must be created proper descendant of FiniteElement class.
190  */
192  /// Same as previous, but represents 1D element.
194  /// Same as previous, but represents 2D element.
196  /// Same as previous, but represents 3D element.
198 
199  /// mesh reader file
201 
202  /// field name read from input
203  std::string field_name_;
204 
205  /// Specify section where to find the field data in input mesh file
207 
208  /// Specify type of FE data interpolation
210 
211  /// Field flags.
213 
214  /// Default value of element if not set in mesh data file
216 
217  /// Accessor to Input::Record
219 
220  /// Is set in set_mesh method. Value true means, that we accept only boundary element accessors in the @p value method.
222 
223  /**
224  * Hold dofs of boundary elements.
225  *
226  * TODO: Temporary solution. Fix problem with merge new DOF handler and boundary Mesh. Will be removed in future.
227  */
228  std::shared_ptr< std::vector<LongIdx> > boundary_dofs_;
229 
230  /// Registrar of class to factory
231  static const int registrar;
232 };
233 
234 
235 
236 /**
237  * Method creates FieldFE of existing VectorMPI that represents elementwise field.
238  *
239  * It's necessary to create new VectorMPI of FieldFE, because DOF handler has generally
240  * a different ordering than mesh.
241  * Then is need to call fill_output_data method.
242  *
243  * Temporary solution that will be remove after solving issue 995.
244  */
245 template <int spacedim, class Value>
246 std::shared_ptr<FieldFE<spacedim, Value> > create_field(VectorMPI & vec_seq, Mesh & mesh, unsigned int n_comp)
247 {
248  std::shared_ptr<DOFHandlerMultiDim> dh; // DOF handler object allow create FieldFE
249  FiniteElement<0> *fe0; // Finite element objects (allow to create DOF handler)
250  FiniteElement<1> *fe1;
251  FiniteElement<2> *fe2;
252  FiniteElement<3> *fe3;
253 
254  switch (n_comp) { // prepare FEM objects for DOF handler by number of components
255  case 1: { // scalar
256  fe0 = new FE_P_disc<0>(0);
257  fe1 = new FE_P_disc<1>(0);
258  fe2 = new FE_P_disc<2>(0);
259  fe3 = new FE_P_disc<3>(0);
260  break;
261  }
262  case 3: { // vector
263  std::shared_ptr< FiniteElement<0> > fe0_ptr = std::make_shared< FE_P_disc<0> >(0);
264  std::shared_ptr< FiniteElement<1> > fe1_ptr = std::make_shared< FE_P_disc<1> >(0);
265  std::shared_ptr< FiniteElement<2> > fe2_ptr = std::make_shared< FE_P_disc<2> >(0);
266  std::shared_ptr< FiniteElement<3> > fe3_ptr = std::make_shared< FE_P_disc<3> >(0);
267  fe0 = new FESystem<0>(fe0_ptr, FEType::FEVector, 3);
268  fe1 = new FESystem<1>(fe1_ptr, FEType::FEVector, 3);
269  fe2 = new FESystem<2>(fe2_ptr, FEType::FEVector, 3);
270  fe3 = new FESystem<3>(fe3_ptr, FEType::FEVector, 3);
271  break;
272  }
273  case 9: { // tensor
274  std::shared_ptr< FiniteElement<0> > fe0_ptr = std::make_shared< FE_P_disc<0> >(0);
275  std::shared_ptr< FiniteElement<1> > fe1_ptr = std::make_shared< FE_P_disc<1> >(0);
276  std::shared_ptr< FiniteElement<2> > fe2_ptr = std::make_shared< FE_P_disc<2> >(0);
277  std::shared_ptr< FiniteElement<3> > fe3_ptr = std::make_shared< FE_P_disc<3> >(0);
278  fe0 = new FESystem<0>(fe0_ptr, FEType::FETensor, 9);
279  fe1 = new FESystem<1>(fe1_ptr, FEType::FETensor, 9);
280  fe2 = new FESystem<2>(fe2_ptr, FEType::FETensor, 9);
281  fe3 = new FESystem<3>(fe3_ptr, FEType::FETensor, 9);
282  break;
283  }
284  default:
285  ASSERT(false).error("Should not happen!\n");
286  }
287 
288  // Prepare DOF handler
289  DOFHandlerMultiDim dh_par(mesh);
290  std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( &mesh, fe0, fe1, fe2, fe3);
291  dh_par.distribute_dofs(ds);
292  dh = dh_par.sequential();
293 
294  // Construct FieldFE
295  std::shared_ptr< FieldFE<spacedim, Value> > field_ptr = std::make_shared< FieldFE<spacedim, Value> >();
296  field_ptr->set_fe_data(dh, 0, VectorMPI::sequential(vec_seq.size()) );
297  return field_ptr;
298 }
299 
300 
301 /**
302  * Fill data to VecSeqDouble in order corresponding with element DOFs.
303  *
304  * Set data to data vector of field in correct order according to values of DOF handler indices.
305  *
306  * Temporary solution that will be remove after solving issue 995.
307  */
308 template <int spacedim, class Value>
309 void fill_output_data(VectorMPI & vec_seq, std::shared_ptr<FieldFE<spacedim, Value> > field_ptr)
310 {
311  auto dh = field_ptr->get_dofhandler();
312  unsigned int ndofs = dh->max_elem_dofs();
313  unsigned int idof; // iterate over indices
314  std::vector<LongIdx> indices(ndofs);
315 
316  /*for (auto cell : dh->own_range()) {
317  cell.get_loc_dof_indices(indices);
318  for(idof=0; idof<ndofs; idof++) field_ptr->get_data_vec()[ indices[idof] ] = (*vec_seq.data_ptr())[ ndofs*cell.elm_idx()+idof ];
319  }*/
320 
321  // Fill DOF handler of FieldFE with correct permutation of data corresponding with DOFs.
322  for (auto ele : dh->mesh()->elements_range()) {
323  dh->cell_accessor_from_element(ele.idx()).get_loc_dof_indices(indices);
324  for(idof=0; idof<ndofs; idof++) field_ptr->get_data_vec()[ indices[idof] ] = (*vec_seq.data_ptr())[ ndofs*ele.idx()+idof ];
325  }
326 }
327 
328 
329 #endif /* FIELD_FE_HH_ */
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
std::shared_ptr< DOFHandlerMultiDim > sequential()
Returns sequential version of the current dof handler.
Definition: dofhandler.cc:54
virtual ~FieldFE()
Destructor.
Definition: field_fe.cc:702
unsigned int size() const
Return size of output data.
Definition: vector_mpi.hh:179
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
virtual void value_list(const std::vector< Point > &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list)
Definition: field_fe.cc:200
void interpolate_gauss(ElementDataCache< double >::ComponentDataPtr data_vec)
Interpolate data (use Gaussian distribution) over all elements of target mesh.
Definition: field_fe.cc:430
FEValueHandler< 0, spacedim, Value > value_handler0_
Value handler that allows get value of 0D elements.
Definition: field_fe.hh:178
std::shared_ptr< std::vector< T > > ComponentDataPtr
P0 interpolation (with the use of calculation of intersections)
Definition: field_fe.hh:64
std::string field_name_
field name read from input
Definition: field_fe.hh:203
std::shared_ptr< FieldFE< spacedim, Value > > create_field(VectorMPI &vec_seq, Mesh &mesh, unsigned int n_comp)
Definition: field_fe.hh:246
VectorMPI get_data_vec() const
Definition: field_fe.hh:141
FEValueHandler< 2, spacedim, Value > value_handler2_
Value handler that allows get value of 2D elements.
Definition: field_fe.hh:182
VectorMPI data_vec_
Store data of Field.
Definition: field_fe.hh:173
void distribute_dofs(std::shared_ptr< DiscreteSpace > ds)
Distributes degrees of freedom on the mesh needed for the given discrete space.
Definition: dofhandler.cc:296
FilePath reader_file_
mesh reader file
Definition: field_fe.hh:200
Definition: mesh.h:76
identical mesh
Definition: field_fe.hh:61
Helper struct stores data for initizalize descentants of FieldAlgorithmBase.
FieldAlgorithmBase< spacedim, Value >::Point Point
Definition: field_fe.hh:53
Class FESystem for compound finite elements.
VectorDataPtr data_ptr()
Getter for shared pointer of output data.
Definition: vector_mpi.hh:114
void fill_output_data(VectorMPI &vec_seq, std::shared_ptr< FieldFE< spacedim, Value > > field_ptr)
Definition: field_fe.hh:309
static const int registrar
Registrar of class to factory.
Definition: field_fe.hh:231
void set_mesh(const Mesh *mesh, bool boundary_domain) override
Definition: field_fe.cc:249
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:346
double default_value_
Default value of element if not set in mesh data file.
Definition: field_fe.hh:215
FEValueHandler< 3, spacedim, Value > value_handler3_
Value handler that allows get value of 3D elements.
Definition: field_fe.hh:184
bool set_time(const TimeStep &time) override
Definition: field_fe.cc:380
FiniteElement< 3 > * fe3_
Same as previous, but represents 3D element.
Definition: field_fe.hh:197
static const Input::Type::Selection & get_interp_selection_input_type()
Definition: field_fe.cc:98
std::shared_ptr< std::vector< LongIdx > > boundary_dofs_
Definition: field_fe.hh:228
void fill_boundary_dofs()
Definition: field_fe.cc:289
FieldAlgorithmBase< spacedim, Value > FactoryBaseType
Definition: field_fe.hh:54
unsigned int data_size() const
Definition: field_fe.cc:695
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:221
OutputTime::DiscreteSpace discretization_
Specify section where to find the field data in input mesh file.
Definition: field_fe.hh:206
std::shared_ptr< DOFHandlerMultiDim > dh_
DOF handler object.
Definition: field_fe.hh:171
static const Input::Type::Record & get_input_type()
Implementation.
Definition: field_fe.cc:55
DataInterpolation interpolation_
Specify type of FE data interpolation.
Definition: field_fe.hh:209
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:99
Definitions of basic Lagrangean finite elements with polynomial shape functions.
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Definition: dofhandler.hh:152
void interpolate_intersection(ElementDataCache< double >::ComponentDataPtr data_vec)
Interpolate data (use intersection library) over all elements of target mesh.
Definition: field_fe.cc:527
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Definition: field_fe.cc:176
Space< spacedim >::Point Point
VectorMPI set_fe_data(std::shared_ptr< DOFHandlerMultiDim > dh, unsigned int component_index=0, VectorMPI dof_values=VectorMPI::sequential(0))
Definition: field_fe.cc:136
FEValueHandler< 1, spacedim, Value > value_handler1_
Value handler that allows get value of 1D elements.
Definition: field_fe.hh:180
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
FiniteElement< 1 > * fe1_
Same as previous, but represents 1D element.
Definition: field_fe.hh:193
std::shared_ptr< DOFHandlerMultiDim > get_dofhandler() const
Definition: field_fe.hh:137
static const Input::Type::Selection & get_disc_selection_input_type()
Definition: field_fe.cc:86
FiniteElement< 2 > * fe2_
Same as previous, but represents 2D element.
Definition: field_fe.hh:195
Input::Record in_rec_
Accessor to Input::Record.
Definition: field_fe.hh:218
static VectorMPI sequential(unsigned int size)
Definition: vector_mpi.hh:67
FieldFE(unsigned int n_comp=0)
Definition: field_fe.cc:127
std::vector< LongIdx > dof_indices_
Array of indexes to data_vec_, used for get/set values.
Definition: field_fe.hh:175
FiniteElement< 0 > * fe0_
Definition: field_fe.hh:191
DataInterpolation
Definition: field_fe.hh:59
virtual void init_from_input(const Input::Record &rec, const struct FieldAlgoBaseInitData &init_data)
Definition: field_fe.cc:226
Record type proxy class.
Definition: type_record.hh:182
Abstract class for description of finite elements.
equivalent mesh (default value)
Definition: field_fe.hh:62
void calculate_native_values(ElementDataCache< double >::ComponentDataPtr data_cache)
Calculate native data over all elements of target mesh.
Definition: field_fe.cc:635
FieldFlag::Flags flags_
Field flags.
Definition: field_fe.hh:212
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:316
Implementation of range helper class.
void native_data_to_cache(ElementDataCache< double > &output_data_cache)
Definition: field_fe.cc:675
P0 interpolation (with the use of Gaussian distribution)
Definition: field_fe.hh:63
Definition: field.hh:56
unsigned int n_comp() const