Flow123d  release_3.0.0-1120-g2970c33
elasticity.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 elasticity.hh
15  * @brief FEM for linear elasticity.
16  * @author Jan Stebel
17  */
18 
19 #ifndef ELASTICITY_HH_
20 #define ELASTICITY_HH_
21 
22 #include "fields/bc_field.hh"
23 #include "fields/field.hh"
24 #include "fields/multi_field.hh"
25 #include "la/linsys.hh"
26 #include "la/vector_mpi.hh"
27 #include "flow/mh_dofhandler.hh"
29 #include "fem/mapping_p1.hh"
30 #include "coupling/equation.hh"
31 #include "fem/fe_values_views.hh"
32 
33 class Distribution;
34 class OutputTime;
35 class DOFHandlerMultiDim;
36 template<unsigned int dim, unsigned int spacedim> class FEValuesBase;
37 template<unsigned int dim> class FiniteElement;
38 template<unsigned int dim, unsigned int spacedim> class Mapping;
39 template<unsigned int dim> class Quadrature;
40 class Elasticity;
41 
42 
43 namespace Mechanics {
44 
45 /**
46  * Auxiliary container class for Finite element and related objects of all dimensions.
47  * Its purpose is to provide templated access to these objects, applicable in
48  * the assembling methods.
49  */
50 class FEObjects {
51 public:
52 
53  FEObjects(Mesh *mesh_, unsigned int fe_order);
54  ~FEObjects();
55 
56  template<unsigned int dim>
57  inline FiniteElement<dim> *fe();
58 
59  template<unsigned int dim>
60  inline Quadrature<dim> *q();
61 
62  template<unsigned int dim>
63  inline MappingP1<dim,3> *mapping();
64 
65  inline std::shared_ptr<DOFHandlerMultiDim> dh();
66 
67 // const FEValuesViews::Vector<dim,3> vec;
68 
69 private:
70 
71  /// Finite elements for the solution of the advection-diffusion equation.
76 
77  /// Quadratures used in assembling methods.
82 
83  /// Auxiliary mappings of reference elements.
87 
88  std::shared_ptr<DiscreteSpace> ds_;
89 
90  /// Object for distribution of dofs.
91  std::shared_ptr<DOFHandlerMultiDim> dh_;
92 };
93 
94 
95 } // namespace Mechanics
96 
97 
98 
99 
100 
101 class Elasticity : public EquationBase
102 {
103 public:
104 
105  class EqData : public FieldSet {
106  public:
107 
108  enum Bc_types {
110  bc_type_traction
111  };
112 
113  EqData();
114 
115  static constexpr const char * name() { return "Mechanics_LinearElasticity"; }
116 
117  static string default_output_field() { return "\"displacement\""; }
118 
119  static const Input::Type::Selection & get_bc_type_selection();
120 
121  static IT::Selection get_output_selection();
122 
123 
130  Field<3, FieldValue<3>::Scalar> fracture_sigma; ///< Transition parameter for diffusive transfer on fractures.
131 
132  /// Pointer to DarcyFlow field cross_section
136 
138 
140 
141  };
142 
143 
144 
145  /**
146  * @brief Constructor.
147  * @param init_mesh computational mesh
148  * @param in_rec input record
149  */
150  Elasticity(Mesh &init_mesh, const Input::Record in_rec);
151  /**
152 
153  * @brief Declare input record type for the equation TransportDG.
154  */
155  static const Input::Type::Record & get_input_type();
156 
157  /**
158  * @brief Initialize solution in the zero time.
159  */
160  void zero_time_step() override;
161 
162  bool evaluate_time_constraint(double &time_constraint)
163  { return false; }
164 
165  /**
166  * @brief Computes the solution in one time instant.
167  */
168  void update_solution() override;
169 
170  /**
171  * @brief Postprocesses the solution and writes to output file.
172  */
173  void output_data();
174 
175  /**
176  * @brief Destructor.
177  */
178  ~Elasticity();
179 
180  void initialize() override;
181 
182  void calculate_cumulative_balance();
183 
184  const Vec &get_solution()
185  { return ls->get_solution(); }
186 
187  inline EqData &data() { return data_; }
188 
189 
191 
192 
193 
194 
195 private:
196  /// Registrar of class to factory
197  static const int registrar;
198 
199  void output_vector_gather();
200 
201  void preallocate();
202 
203  /**
204  * @brief Assembles the stiffness matrix.
205  *
206  * This routine just calls assemble_volume_integrals(), assemble_fluxes_boundary(),
207  * assemble_fluxes_element_element() and assemble_fluxes_element_side() for each
208  * space dimension.
209  */
210  void assemble_stiffness_matrix();
211 
212  /**
213  * @brief Assembles the volume integrals into the stiffness matrix.
214  */
215  template<unsigned int dim>
216  void assemble_volume_integrals();
217 
218  /**
219  * @brief Assembles the right hand side due to volume sources.
220  *
221  * This method just calls set_sources() for each space dimension.
222  */
223  void set_sources();
224 
225  /**
226  * @brief Assembles the right hand side vector due to volume sources.
227  */
228  template<unsigned int dim>
229  void set_sources();
230 
231  /**
232  * @brief Assembles the fluxes on the boundary.
233  */
234  template<unsigned int dim>
235  void assemble_fluxes_boundary();
236 
237  /**
238  * @brief Assembles the fluxes between elements of different dimensions.
239  */
240  template<unsigned int dim>
241  void assemble_fluxes_element_side();
242 
243 
244  /**
245  * @brief Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions.
246  *
247  * The routine just calls templated method set_boundary_condition() for each space dimension.
248  */
249  void set_boundary_conditions();
250 
251  /**
252  * @brief Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions
253  * for a given space dimension.
254  */
255  template<unsigned int dim>
256  void set_boundary_conditions();
257 
258 
259 
260 
261  /// @name Physical parameters
262  // @{
263 
264  /// Field data for model parameters.
266 
267  // @}
268 
269 
270  /// @name Parameters of the numerical method
271  // @{
272 
273  /// Finite element objects
275 
276  // @}
277 
278 
279 
280  /// @name Solution of algebraic system
281  // @{
282 
283  /// Vector of right hand side.
284  Vec rhs;
285 
286  /// The stiffness matrix.
288 
289  /// Linear algebra system for the transport equation.
291 
292  // @}
293 
294 
295  /// @name Output to file
296  // @{
297 
298  /// Vector of solution data.
300 
301  std::shared_ptr<OutputTime> output_stream_;
302 
303  /// Record with input specification.
305 
306 
307  // @}
308 
309 
310 
311 
312  /// @name Other
313  // @{
314 
315  /// Indicates whether matrices have been preallocated.
317 
318  /// Indicator of change in advection vector field.
320 
321  // @}
322 };
323 
324 
325 
326 
327 
328 
329 #endif /* ELASTICITY_HH_ */
BCField< 3, FieldValue< 3 >::Enum > bc_type
Definition: elasticity.hh:124
std::shared_ptr< DiscreteSpace > ds_
Definition: elasticity.hh:88
BCField< 3, FieldValue< 3 >::VectorFixed > bc_displacement
Definition: elasticity.hh:125
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
Field< 3, FieldValue< 3 >::Scalar > young_modulus
Definition: elasticity.hh:128
Abstract base class for equation clasess.
static constexpr const char * name()
Definition: elasticity.hh:115
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:71
FiniteElement< 3 > * fe3_
Definition: elasticity.hh:75
Field< 3, FieldValue< 3 >::Scalar > region_id
Definition: elasticity.hh:134
const Vec & get_solution()
Definition: elasticity.hh:184
FEObjects(Mesh *mesh_, unsigned int fe_order)
Definition: elasticity.cc:77
Wrappers for linear systems based on MPIAIJ and MATIS format.
Class template representing a field with values dependent on: point, element, and region...
Definition: field.hh:83
VectorMPI output_vec
Vector of solution data.
Definition: elasticity.hh:299
FiniteElement< dim > * fe()
LinSys * ls
Linear algebra system for the transport equation.
Definition: elasticity.hh:290
Definition: mesh.h:76
Field< 3, FieldValue< 3 >::Scalar > cross_section
Pointer to DarcyFlow field cross_section.
Definition: elasticity.hh:133
MappingP1< 2, 3 > * map2_
Definition: elasticity.hh:85
bool flux_changed
Indicator of change in advection vector field.
Definition: elasticity.hh:319
Input::Record input_rec
Record with input specification.
Definition: elasticity.hh:304
Quadrature< dim > * q()
Vec rhs
Vector of right hand side.
Definition: elasticity.hh:284
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:35
bool allocation_done
Indicates whether matrices have been preallocated.
Definition: elasticity.hh:316
MappingP1< 3, 3 > * map3_
Definition: elasticity.hh:86
Quadrature< 0 > * q0_
Quadratures used in assembling methods.
Definition: elasticity.hh:78
MappingP1< 1, 3 > * map1_
Auxiliary mappings of reference elements.
Definition: elasticity.hh:84
BCField< 3, FieldValue< 3 >::VectorFixed > bc_traction
Definition: elasticity.hh:126
Abstract class for the mapping between reference and actual cell.
Definition: fe_values.hh:38
Quadrature< 2 > * q2_
Definition: elasticity.hh:80
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
FiniteElement< 0 > * fe0_
Finite elements for the solution of the advection-diffusion equation.
Definition: elasticity.hh:72
std::shared_ptr< DOFHandlerMultiDim > dh()
Definition: elasticity.cc:141
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Definition: dofhandler.hh:152
FiniteElement< 1 > * fe1_
Definition: elasticity.hh:73
The class for outputting data during time.
Definition: output_time.hh:50
Mat stiffness_matrix
The stiffness matrix.
Definition: elasticity.hh:287
EqData & data()
Definition: elasticity.hh:187
Field< 3, FieldValue< 3 >::VectorFixed > load
Definition: elasticity.hh:127
Field< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures.
Definition: elasticity.hh:130
Elasticity FactoryBaseType
Definition: elasticity.hh:190
Field< 3, FieldValue< 3 >::Scalar > poisson_ratio
Definition: elasticity.hh:129
EquationOutput output_fields
Definition: elasticity.hh:139
static const int registrar
Registrar of class to factory.
Definition: elasticity.hh:197
EqData data_
Field data for model parameters.
Definition: elasticity.hh:265
Quadrature< 1 > * q1_
Definition: elasticity.hh:79
static string default_output_field()
Definition: elasticity.hh:117
std::shared_ptr< DOFHandlerMultiDim > dh_
Object for distribution of dofs.
Definition: elasticity.hh:91
Record type proxy class.
Definition: type_record.hh:182
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Base class for FEValues and FESideValues.
Definition: fe_values.hh:37
MappingP1< dim, 3 > * mapping()
Quadrature< 3 > * q3_
Definition: elasticity.hh:81
FiniteElement< 2 > * fe2_
Definition: elasticity.hh:74
bool evaluate_time_constraint(double &time_constraint)
Definition: elasticity.hh:162
Field< 3, FieldValue< 3 >::VectorFixed > output_field
Definition: elasticity.hh:137
Template for classes storing finite set of named values.
Field< 3, FieldValue< 3 >::Scalar > subdomain
Definition: elasticity.hh:135
std::shared_ptr< OutputTime > output_stream_
Definition: elasticity.hh:301
Mechanics::FEObjects * feo
Finite element objects.
Definition: elasticity.hh:274