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