ExaDG
Loading...
Searching...
No Matches
st_venant_kirchhoff.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef EXADG_STRUCTURE_MATERIAL_LIBRARY_ST_VENANT_KIRCHHOFF_H_
23#define EXADG_STRUCTURE_MATERIAL_LIBRARY_ST_VENANT_KIRCHHOFF_H_
24
25// deal.II
26#include <deal.II/base/function.h>
27#include <deal.II/base/point.h>
28#include <deal.II/matrix_free/matrix_free.h>
29
30// ExaDG
31#include <exadg/matrix_free/integrators.h>
32#include <exadg/operators/variable_coefficients.h>
33#include <exadg/structure/material/material.h>
34
35namespace ExaDG
36{
37namespace Structure
38{
39template<int dim>
40struct StVenantKirchhoffData : public MaterialData
41{
42 StVenantKirchhoffData(
43 MaterialType const & type,
44 double const & youngs_modulus,
45 double const & poissons_ratio,
46 Type2D const & type_two_dim,
47 std::shared_ptr<dealii::Function<dim>> const youngs_modulus_function = nullptr)
48 : MaterialData(type),
49 youngs_modulus(youngs_modulus),
50 youngs_modulus_function(youngs_modulus_function),
51 poissons_ratio(poissons_ratio),
52 type_two_dim(type_two_dim)
53 {
54 }
55
56 double youngs_modulus;
57 std::shared_ptr<dealii::Function<dim>> youngs_modulus_function;
58
59 double poissons_ratio;
60 Type2D type_two_dim;
61};
62
63template<int dim, typename Number>
64class StVenantKirchhoff : public Material<dim, Number>
65{
66public:
67 typedef typename Material<dim, Number>::VectorType VectorType;
68 typedef typename Material<dim, Number>::Range Range;
69 typedef typename Material<dim, Number>::IntegratorCell IntegratorCell;
70
71 typedef typename Material<dim, Number>::scalar scalar;
72 typedef typename Material<dim, Number>::tensor tensor;
73 typedef typename Material<dim, Number>::symmetric_tensor symmetric_tensor;
74
75 StVenantKirchhoff(dealii::MatrixFree<dim, Number> const & matrix_free,
76 unsigned int const dof_index,
77 unsigned int const quad_index,
78 StVenantKirchhoffData<dim> const & data,
79 bool const large_deformation);
80
81 symmetric_tensor
82 second_piola_kirchhoff_stress(tensor const & gradient_displacement,
83 unsigned int const cell,
84 unsigned int const q) const final;
85
86 symmetric_tensor
87 second_piola_kirchhoff_stress_displacement_derivative(tensor const & gradient_increment,
88 tensor const & deformation_gradient,
89 unsigned int const cell,
90 unsigned int const q) const final;
91
92private:
93 /*
94 * Factor out coefficients for faster computation. Note that these factors do not contain the
95 * (potentially variable) Young's modulus.
96 */
97 Number
98 get_f0_factor(Number const & poissons_ratio, Type2D const type_two_dim) const;
99
100 Number
101 get_f1_factor(Number const & poissons_ratio, Type2D const type_two_dim) const;
102
103 Number
104 get_f2_factor(Number const & poissons_ratio, Type2D const type_two_dim) const;
105
106 /*
107 * The second Piola-Kirchhoff stress tensor S is given as S = lambda * I * tr(E) + 2 mu E, with E
108 * being the Green-Lagrange strain tensor and Lamee parameters lambda and mu. This leads to
109 * Sii = f0 * Eii + f1 * sum_{j = 1, ..., dim; i!=j} Eij, for i = 1, ..., dim, and
110 * Sij = f2 * (Eij + Eji), for i, j = 1, ..., dim and i != j.
111 * The latter symmetrizes the off-diagonal entries in the strain argument to reduce computations.
112 */
113 symmetric_tensor
114 second_piola_kirchhoff_stress_symmetrize(tensor const & strain,
115 unsigned int const cell,
116 unsigned int const q) const;
117
118 /*
119 * Store factors involving (potentially variable) Young's modulus.
120 */
121 void
122 cell_loop_set_coefficients(dealii::MatrixFree<dim, Number> const & matrix_free,
123 VectorType &,
124 VectorType const & src,
125 Range const & cell_range) const;
126
127 unsigned int dof_index;
128 unsigned int quad_index;
129
130 StVenantKirchhoffData<dim> const & data;
131
132 bool large_deformation;
133
134 Number const f0_factor;
135 Number const f1_factor;
136 Number const f2_factor;
137
138 mutable dealii::VectorizedArray<Number> f0;
139 mutable dealii::VectorizedArray<Number> f1;
140 mutable dealii::VectorizedArray<Number> f2;
141
142 // cache coefficients for spatially varying material parameters
143 bool youngs_modulus_is_variable;
147};
148
149} // namespace Structure
150} // namespace ExaDG
151
152#endif /* EXADG_STRUCTURE_MATERIAL_LIBRARY_ST_VENANT_KIRCHHOFF_H_ */
Definition material.h:38
Definition variable_coefficients.h:40
Definition driver.cpp:33
Definition st_venant_kirchhoff.h:41