ExaDG
Loading...
Searching...
No Matches
structure_calculators.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2025 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_OPERATORS_STRUCTURE_CALCULATORS_H_
23#define EXADG_OPERATORS_STRUCTURE_CALCULATORS_H_
24
25// deal.II
26#include <deal.II/lac/la_parallel_vector.h>
27
28// ExaDG
29#include <exadg/matrix_free/integrators.h>
30#include <exadg/structure/spatial_discretization/operators/continuum_mechanics.h>
31#include <exadg/structure/spatial_discretization/operators/elasticity_operator_base.h>
32
33namespace ExaDG
34{
35/*
36 * Calculator for the Jacobian of the vector field u(x) defined as
37 *
38 * J := det(F) with F := I + grad(u),
39 *
40 * where F is the deformation gradient tensor and grad(u) is the gradient of the vector field.
41 *
42 */
43template<int dim, typename Number>
44class DisplacementJacobianCalculator
45{
46public:
47 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
48
49 typedef DisplacementJacobianCalculator<dim, Number> This;
50
51 typedef CellIntegrator<dim, dim, Number> CellIntegratorVector;
52 typedef CellIntegrator<dim, 1, Number> CellIntegratorScalar;
53
54 typedef dealii::VectorizedArray<Number> scalar;
55 typedef dealii::Tensor<2, dim, dealii::VectorizedArray<Number>> tensor;
56
57 DisplacementJacobianCalculator();
58
59 void
60 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
61 unsigned int const dof_index_vector_in,
62 unsigned int const dof_index_scalar_in,
63 unsigned int const quad_index_in);
64
65 /*
66 * Compute the right-hand side of an L2 projection of the Jacobian of the vector field.
67 */
68 void
69 compute_projection_rhs(VectorType & dst, VectorType const & src) const;
70
71private:
72 void
73 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
74 VectorType & dst_scalar_valued,
75 VectorType const & src_vector_valued,
76 std::pair<unsigned int, unsigned int> const & cell_range) const;
77
78 dealii::MatrixFree<dim, Number> const * matrix_free;
79
80 unsigned int dof_index_vector;
81 unsigned int dof_index_scalar;
82 unsigned int quad_index;
83};
84
85/*
86 * Calculator for the maximum principal stress being the largest absolute Eigenvalue of the Cauchy
87 * stress tensor
88 *
89 * sigma_I := max(|lambda_i|) for i = 1,...,dim
90 *
91 * with lambda_i being the Eigenvalues of the Cauchy stress tensor sigma,
92 *
93 * sigma = (1/det(F)) * F * S * F^T
94 *
95 * where F is the deformation gradient tensor and S is the second Piola-Kirchhoff stress tensor.
96 *
97 */
98template<int dim, typename Number>
99class MaxPrincipalStressCalculator
100{
101public:
102 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
103
104 typedef MaxPrincipalStressCalculator<dim, Number> This;
105
106 typedef CellIntegrator<dim, dim, Number> CellIntegratorVector;
107 typedef CellIntegrator<dim, 1, Number> CellIntegratorScalar;
108
109 typedef dealii::VectorizedArray<Number> scalar;
110 typedef dealii::Tensor<2, dim, dealii::VectorizedArray<Number>> tensor;
111 typedef dealii::SymmetricTensor<2, dim, dealii::VectorizedArray<Number>> symmetric_tensor;
112
113 MaxPrincipalStressCalculator();
114
115 void
116 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
117 unsigned int const dof_index_vector_in,
118 unsigned int const dof_index_scalar_in,
119 unsigned int const quad_index_in,
120 Structure::ElasticityOperatorBase<dim, Number> const & elasticity_operator_base_in);
121
122 /*
123 * Compute the right-hand side of an L2 projection of the maximum principal stress.
124 */
125 void
126 compute_projection_rhs(VectorType & dst, VectorType const & src) const;
127
128private:
129 void
130 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
131 VectorType & dst_scalar_valued,
132 VectorType const & src_vector_valued,
133 std::pair<unsigned int, unsigned int> const & cell_range) const;
134
135 dealii::MatrixFree<dim, Number> const * matrix_free;
136
137 unsigned int dof_index_vector;
138 unsigned int dof_index_scalar;
139 unsigned int quad_index;
140
141 // Pointer to the underlying elasticity operator to compute stresses dependent on the material
142 // model.
143 dealii::ObserverPointer<Structure::ElasticityOperatorBase<dim, Number> const>
144 elasticity_operator_base;
145};
146
147} // namespace ExaDG
148
149#endif /* EXADG_OPERATORS_STRUCTURE_CALCULATORS_H_ */
Definition elasticity_operator_base.h:72
Definition driver.cpp:33