-
Notifications
You must be signed in to change notification settings - Fork 9
/
MechanicalPhysics.hh
170 lines (149 loc) · 4.74 KB
/
MechanicalPhysics.hh
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
/* Copyright (c) 2022 - 2024, the adamantine authors.
*
* This file is subject to the Modified BSD License and may not be distributed
* without copyright and license information. Please refer to the file LICENSE
* for the text and further information on this license.
*/
#ifndef MECHANICAL_PHYSICS_HH
#define MECHANICAL_PHYSICS_HH
#include <Geometry.hh>
#include <MechanicalOperator.hh>
#include <deal.II/base/quadrature_lib.h>
#include <deal.II/base/symmetric_tensor.h>
#include <deal.II/fe/fe_system.h>
#include <deal.II/hp/fe_collection.h>
namespace adamantine
{
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
class MechanicalPhysics
{
public:
/**
* Constructor.
*/
MechanicalPhysics(MPI_Comm const &communicator, unsigned int const fe_degree,
Geometry<dim> &geometry,
MaterialProperty<dim, p_order, MaterialStates,
MemorySpaceType> &material_properties,
std::vector<double> const &initial_temperatures);
/**
* Setup the DoFHandler, the AffineConstraints, and the
* MechanicalOperator.
*/
void
setup_dofs(std::vector<std::shared_ptr<BodyForce<dim>>> const &body_forces =
std::vector<std::shared_ptr<BodyForce<dim>>>());
/**
* Same as above when solving a thermo-mechanical problem.
*/
void setup_dofs(
dealii::DoFHandler<dim> const &thermal_dof_handler,
dealii::LA::distributed::Vector<double, dealii::MemorySpace::Host> const
&temperature,
std::vector<bool> const &has_melted,
std::vector<std::shared_ptr<BodyForce<dim>>> const &body_forces =
std::vector<std::shared_ptr<BodyForce<dim>>>());
/**
* Solve the mechanical problem and return the displacement.
*/
dealii::LA::distributed::Vector<double, dealii::MemorySpace::Host> solve();
/**
* Return the DoFHandler.
*/
dealii::DoFHandler<dim> &get_dof_handler();
/**
* Return the AffineConstraints<double>.
*/
dealii::AffineConstraints<double> &get_affine_constraints();
/**
* Return the stress tensor associated to each quadrature point.
*/
std::vector<std::vector<dealii::SymmetricTensor<2, dim>>> &
get_stress_tensor();
private:
// Compute the stress using linear combination of isotropic and kinematic
// hardening in the book Plasticity Modeling & Computation from Ronaldo I.
// Borja.
void compute_stress(
dealii::LA::distributed::Vector<double, dealii::MemorySpace::Host> const
&displacement);
/**
* Associated Geometry.
*/
Geometry<dim> &_geometry;
/**
* Associated MaterialProperty.
*/
MaterialProperty<dim, p_order, MaterialStates, MemorySpaceType>
&_material_properties;
/**
* Associated FECollection.
*/
dealii::hp::FECollection<dim> _fe_collection;
/**
* Associated DoFHandler.
*/
dealii::DoFHandler<dim> _dof_handler;
/**
* Associated AffineConstraints.
*/
dealii::AffineConstraints<double> _affine_constraints;
/**
* Associated QCollection.
*/
dealii::hp::QCollection<dim> _q_collection;
/**
* Pointer to the MechanicalOperator
*/
std::unique_ptr<
MechanicalOperator<dim, p_order, MaterialStates, MemorySpaceType>>
_mechanical_operator;
/**
* Whether to include a gravitional body force in the calculation.
*/
bool _include_gravity;
/**
* Save displacement from the previous time step.
*/
dealii::LA::distributed::Vector<double, dealii::MemorySpace::Host>
_old_displacement;
/**
* Plastic internal variable related to the strain
*/
std::vector<std::vector<double>> _plastic_internal_variable;
/**
* Stress tensor at each (cell, quadrature point).
*/
std::vector<std::vector<dealii::SymmetricTensor<2, dim>>> _stress;
/**
* Back stress tensor at each (cell, quadrature point).
*/
std::vector<std::vector<dealii::SymmetricTensor<2, dim>>> _back_stress;
};
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
inline dealii::DoFHandler<dim> &
MechanicalPhysics<dim, p_order, MaterialStates,
MemorySpaceType>::get_dof_handler()
{
return _dof_handler;
}
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
inline dealii::AffineConstraints<double> &
MechanicalPhysics<dim, p_order, MaterialStates,
MemorySpaceType>::get_affine_constraints()
{
return _affine_constraints;
}
template <int dim, int p_order, typename MaterialStates,
typename MemorySpaceType>
inline std::vector<std::vector<dealii::SymmetricTensor<2, dim>>> &
MechanicalPhysics<dim, p_order, MaterialStates,
MemorySpaceType>::get_stress_tensor()
{
return _stress;
}
} // namespace adamantine
#endif