ExaDG
Loading...
Searching...
No Matches
nonlinear_operator.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_SPATIAL_DISCRETIZATION_OPERATORS_NONLINEAR_OPERATOR_H_
23#define EXADG_STRUCTURE_SPATIAL_DISCRETIZATION_OPERATORS_NONLINEAR_OPERATOR_H_
24
25// ExaDG
26#include <exadg/structure/spatial_discretization/operators/elasticity_operator_base.h>
27
28namespace ExaDG
29{
30namespace Structure
31{
32template<int dim, typename Number>
33class NonLinearOperator : public ElasticityOperatorBase<dim, Number>
34{
35private:
36 typedef ElasticityOperatorBase<dim, Number> Base;
37
38 typedef typename Base::VectorType VectorType;
39 typedef typename Base::IntegratorCell IntegratorCell;
40 typedef typename Base::IntegratorFace IntegratorFace;
41
43
44 typedef std::pair<unsigned int, unsigned int> Range;
45
46 typedef dealii::VectorizedArray<Number> scalar;
47 typedef dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> vector;
48 typedef dealii::Tensor<2, dim, dealii::VectorizedArray<Number>> tensor;
49
50public:
54 void
55 initialize(dealii::MatrixFree<dim, Number> const & matrix_free,
56 dealii::AffineConstraints<Number> const & affine_constraints,
57 OperatorData<dim> const & data) override;
58
68 void
69 evaluate_nonlinear(VectorType & dst, VectorType const & src) const;
70
74 bool
75 valid_deformation(VectorType const & displacement) const;
76
80 void
81 set_solution_linearization(VectorType const & vector) const;
82
86 VectorType const &
88
89private:
90 /*
91 * Non-linear operator.
92 */
93 void
94 reinit_cell_nonlinear(IntegratorCell & integrator, unsigned int const cell) const;
95
96 void
97 cell_loop_nonlinear(dealii::MatrixFree<dim, Number> const & matrix_free,
98 VectorType & dst,
99 VectorType const & src,
100 Range const & range) const;
101
102 void
103 face_loop_nonlinear(dealii::MatrixFree<dim, Number> const & matrix_free,
104 VectorType & dst,
105 VectorType const & src,
106 Range const & range) const;
107
108 void
109 boundary_face_loop_nonlinear(dealii::MatrixFree<dim, Number> const & matrix_free,
110 VectorType & dst,
111 VectorType const & src,
112 Range const & range) const;
113
114 /*
115 * A cell loop that checks whether the Jacobian determinant is positive for all q-points.
116 * Prior to calling this function, inhomogeneous Dirichlet degrees of freedom need to be
117 * set correctly.
118 */
119 void
120 cell_loop_valid_deformation(dealii::MatrixFree<dim, Number> const & matrix_free,
121 Number & dst,
122 VectorType const & src,
123 Range const & range) const;
124
125 /*
126 * Calculates the integral
127 *
128 * (v_h, factor * d_h)_Omega + (Grad(v_h), P_h)_Omega
129 *
130 * with 1st Piola-Kirchhoff stress tensor P_h
131 *
132 * P_h = F_h * S_h ,
133 *
134 * 2nd Piola-Kirchhoff stress tensor S_h
135 *
136 * S_h = function(E_h) ,
137 *
138 * Green-Lagrange strain tensor E_h
139 *
140 * E_h = 1/2 (F_h^T * F_h - 1) ,
141 *
142 * material deformation gradient F_h
143 *
144 * F_h = 1 + Grad(d_h) ,
145 *
146 * where
147 *
148 * d_h denotes the displacement vector.
149 */
150 void
151 do_cell_integral_nonlinear(IntegratorCell & integrator) const;
152
153 /*
154 * Computes Neumann BC integral
155 *
156 * - (v_h, t_0)_{Gamma_N}
157 *
158 * with traction
159 *
160 * t_0 = da/dA t .
161 *
162 * If the traction is specified as force per surface area of the underformed
163 * body, the specified traction t is interpreted as t_0 = t, and no pull-back
164 * is necessary.
165 */
166 void
167 do_boundary_integral_continuous(IntegratorFace & integrator_m,
168 dealii::types::boundary_id const & boundary_id) const final;
169
170 /*
171 * Linearized operator.
172 */
173 void
174 reinit_cell_derived(IntegratorCell & integrator, unsigned int const cell) const final;
175
176 /*
177 * Calculates the integral
178 *
179 * (v_h, factor * delta d_h)_Omega + (Grad(v_h), delta P_h)_Omega
180 *
181 * with the directional derivative of the 1st Piola-Kirchhoff stress tensor P_h
182 *
183 * delta P_h = d(P)/d(d)|_{d_lin} * delta d_h ,
184 *
185 * with the point of linearization
186 *
187 * d_lin ,
188 *
189 * and displacement increment
190 *
191 * delta d_h .
192 *
193 * Computing the linearization yields
194 *
195 * delta P_h = + Grad(delta_d) * S(d_lin)
196 * + F(d_lin) * (C_lin : (F^T(d_lin) * Grad(delta d))) .
197 *
198 * Note that a dependency of the Neumann BC on the displacements d through
199 * the area ratio da/dA = function(d) is neglected in the linearization.
200 */
201 void
202 do_cell_integral(IntegratorCell & integrator) const override;
203
204 mutable std::shared_ptr<IntegratorCell> integrator_lin;
205 mutable VectorType displacement_lin;
206};
207
208} // namespace Structure
209} // namespace ExaDG
210
211#endif /* EXADG_STRUCTURE_SPATIAL_DISCRETIZATION_OPERATORS_NONLINEAR_OPERATOR_H_ */
Definition nonlinear_operator.h:34
VectorType const & get_solution_linearization() const
Definition nonlinear_operator.cpp:99
bool valid_deformation(VectorType const &displacement) const
Definition nonlinear_operator.cpp:62
void initialize(dealii::MatrixFree< dim, Number > const &matrix_free, dealii::AffineConstraints< Number > const &affine_constraints, OperatorData< dim > const &data) override
Definition nonlinear_operator.cpp:32
void evaluate_nonlinear(VectorType &dst, VectorType const &src) const
Definition nonlinear_operator.cpp:49
void set_solution_linearization(VectorType const &vector) const
Definition nonlinear_operator.cpp:84
Definition driver.cpp:33
Definition elasticity_operator_base.h:37