ExaDG
Loading...
Searching...
No Matches
elasticity_operator_base.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_ELASTICITY_OPERATOR_BASE_H_
23#define EXADG_STRUCTURE_SPATIAL_DISCRETIZATION_OPERATORS_ELASTICITY_OPERATOR_BASE_H_
24
25// ExaDG
26#include <exadg/operators/operator_base.h>
27#include <exadg/structure/material/material_handler.h>
28#include <exadg/structure/user_interface/boundary_descriptor.h>
29#include <exadg/structure/user_interface/material_descriptor.h>
30
31namespace ExaDG
32{
33namespace Structure
34{
35template<int dim>
36struct OperatorData : public OperatorBaseData
37{
38 OperatorData()
39 : OperatorBaseData(),
40 large_deformation(false),
41 pull_back_traction(false),
42 unsteady(false),
43 density(1.0),
44 quad_index_gauss_lobatto(0)
45 {
46 }
47
48 std::shared_ptr<BoundaryDescriptor<dim> const> bc;
49 std::shared_ptr<MaterialDescriptor const> material_descriptor;
50
51 // Boolean parameter differentiating between linear elasticity and finite strain theory
52 bool large_deformation;
53
54 // This parameter is only relevant for nonlinear operator
55 // with large deformations. When set to true, the traction t
56 // is pulled back to the reference configuration, t_0 = da/dA t.
57 bool pull_back_traction;
58
59 // activates mass operator in operator evaluation for unsteady problems
60 bool unsteady;
61
62 // density
63 double density;
64
65 // for DirichletCached boundary conditions, another quadrature rule
66 // is needed to set the constrained DoFs.
67 unsigned int quad_index_gauss_lobatto;
68};
69
70template<int dim, typename Number>
71class ElasticityOperatorBase : public OperatorBase<dim, Number, dim>
72{
73public:
74 typedef Number value_type;
75
76protected:
77 typedef OperatorBase<dim, Number, dim> Base;
78 typedef typename Base::IntegratorCell IntegratorCell;
79 typedef typename Base::VectorType VectorType;
80 typedef typename Base::IntegratorFace IntegratorFace;
81
82public:
83 ElasticityOperatorBase();
84
85 virtual ~ElasticityOperatorBase()
86 {
87 }
88
90 get_integrator_flags(bool const unsteady) const;
91
92 static MappingFlags
93 get_mapping_flags();
94
95 void
96 initialize(dealii::MatrixFree<dim, Number> const & matrix_free,
97 dealii::AffineConstraints<Number> const & affine_constraints,
98 OperatorData<dim> const & data,
99 bool const assemble_matrix);
100
101 OperatorData<dim> const &
102 get_data() const;
103
104 /*
105 * Provide near null space basis vectors, that is, the rigid body modes, used e.g. in AMG setup.
106 */
107 void
108 get_constant_modes(std::vector<std::vector<bool>> & constant_modes,
109 std::vector<std::vector<double>> & constant_modes_values) const override
110 {
111 (void)constant_modes;
112
113 dealii::DoFHandler<dim> const & dof_handler =
114 this->matrix_free->get_dof_handler(this->get_dof_index());
115
116 if(dof_handler.has_level_dofs())
117 {
118 // Extract coarse level constant modes.
119 constant_modes_values = dealii::DoFTools::extract_level_rigid_body_modes(
120 this->matrix_free->get_mg_level(),
121 *this->matrix_free->get_mapping_info().mapping,
122 dof_handler,
123 dealii::ComponentMask(dim, true));
124 }
125 else
126 {
127 // Extract finest level constant modes.
128 constant_modes_values =
129 dealii::DoFTools::extract_rigid_body_modes(*this->matrix_free->get_mapping_info().mapping,
130 dof_handler,
131 dealii::ComponentMask(dim, true));
132 }
133 }
134
135 void
136 set_scaling_factor_mass_operator(double const scaling_factor) const;
137
138 double
139 get_scaling_factor_mass_operator() const;
140
141 void
142 set_inhomogeneous_constrained_values(VectorType & dst) const final;
143
144protected:
145 void
146 reinit_cell_derived(IntegratorCell & integrator, unsigned int const cell) const override;
147
148 virtual void
149 initialize_derived(){};
150
151 OperatorData<dim> operator_data;
152
153 mutable MaterialHandler<dim, Number> material_handler;
154
155 mutable double scaling_factor_mass;
156};
157
158} // namespace Structure
159} // namespace ExaDG
160
161#endif /* EXADG_STRUCTURE_SPATIAL_DISCRETIZATION_OPERATORS_ELASTICITY_OPERATOR_BASE_H_ */
Definition material_handler.h:39
Definition driver.cpp:33
Definition integrator_flags.h:31
Definition mapping_flags.h:31
Definition elasticity_operator_base.h:37