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 STRUCTURE_MATERIAL_LIBRARY_STVENANTKIRCHHOFF
23#define STRUCTURE_MATERIAL_LIBRARY_STVENANTKIRCHHOFF
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>
41{
42 StVenantKirchhoffData(MaterialType const & type,
43 double const & E,
44 double const & nu,
45 Type2D const & type_two_dim,
46 std::shared_ptr<dealii::Function<dim>> const E_function = nullptr)
47 : MaterialData(type), E(E), E_function(E_function), nu(nu), type_two_dim(type_two_dim)
48 {
49 }
50
51 double E;
52 std::shared_ptr<dealii::Function<dim>> E_function;
53
54 double nu;
55 Type2D type_two_dim;
56};
57
58template<int dim, typename Number>
59class StVenantKirchhoff : public Material<dim, Number>
60{
61public:
62 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
63 typedef std::pair<unsigned int, unsigned int> Range;
64 typedef CellIntegrator<dim, dim, Number> IntegratorCell;
65
66 StVenantKirchhoff(dealii::MatrixFree<dim, Number> const & matrix_free,
67 unsigned int const dof_index,
68 unsigned int const quad_index,
69 StVenantKirchhoffData<dim> const & data,
70 bool const large_deformation);
71
72 dealii::Tensor<2, dim, dealii::VectorizedArray<Number>>
73 second_piola_kirchhoff_stress(
74 dealii::Tensor<2, dim, dealii::VectorizedArray<Number>> const & gradient_displacement,
75 unsigned int const cell,
76 unsigned int const q) const final;
77
78 dealii::Tensor<2, dim, dealii::VectorizedArray<Number>>
79 second_piola_kirchhoff_stress_displacement_derivative(
80 dealii::Tensor<2, dim, dealii::VectorizedArray<Number>> const & gradient_increment,
81 dealii::Tensor<2, dim, dealii::VectorizedArray<Number>> const & deformation_gradient,
82 unsigned int const cell,
83 unsigned int const q) const final;
84
85private:
86 /*
87 * Factor out coefficients for faster computation. Note that these factors do not contain the
88 * (potentially variable) Young's modulus.
89 */
90 Number
91 get_f0_factor() const;
92
93 Number
94 get_f1_factor() const;
95
96 Number
97 get_f2_factor() const;
98
99 /*
100 * The second Piola-Kirchhoff stress tensor S is given as S = lambda * I * tr(E) + 2 mu E, with E
101 * being the Green-Lagrange strain tensor and Lamee parameters lambda and mu. This leads to
102 * Sii = f0 * Eii + f1 * sum_{j = 1, ..., dim; i!=j} Eij, for i = 1, ..., dim, and
103 * Sij = f2 * (Eij + Eji), for i, j = 1, ..., dim and i != j.
104 * The latter symmetrizes the off-diagonal entries in the strain argument to reduce computations.
105 */
106 dealii::Tensor<2, dim, dealii::VectorizedArray<Number>>
107 second_piola_kirchhoff_stress_symmetrize(
108 dealii::Tensor<2, dim, dealii::VectorizedArray<Number>> const & strain,
109 unsigned int const cell,
110 unsigned int const q) const;
111
112 /*
113 * Store factors involving (potentially variable) Young's modulus.
114 */
115 void
116 cell_loop_set_coefficients(dealii::MatrixFree<dim, Number> const & matrix_free,
117 VectorType &,
118 VectorType const & src,
119 Range const & cell_range) const;
120
121 unsigned int dof_index;
122 unsigned int quad_index;
123
124 StVenantKirchhoffData<dim> const & data;
125
126 bool large_deformation;
127
128 mutable dealii::VectorizedArray<Number> f0;
129 mutable dealii::VectorizedArray<Number> f1;
130 mutable dealii::VectorizedArray<Number> f2;
131
132 // cache coefficients for spatially varying material parameters
133 bool E_is_variable;
137};
138} // namespace Structure
139} // namespace ExaDG
140
141#endif
Definition material.h:38
Definition st_venant_kirchhoff.h:60
Definition variable_coefficients.h:40
Definition driver.cpp:33
Definition material_data.h:32
Definition st_venant_kirchhoff.h:41