-
Notifications
You must be signed in to change notification settings - Fork 284
Expand file tree
/
Copy pathU_Pw_base_element.h
More file actions
186 lines (131 loc) · 7.98 KB
/
U_Pw_base_element.h
File metadata and controls
186 lines (131 loc) · 7.98 KB
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
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
// KRATOS___
// // ) )
// // ___ ___
// // ____ //___) ) // ) )
// // / / // // / /
// ((____/ / ((____ ((___/ / MECHANICS
//
// License: geo_mechanics_application/license.txt
//
// Main authors: Ignasi de Pouplana,
// Vahid Galavi
//
#pragma once
// Project includes
#include "containers/array_1d.h"
#include "custom_retention/retention_law.h"
#include "custom_retention/retention_law_factory.h"
#include "geometries/geometry.h"
#include "includes/constitutive_law.h"
#include "includes/define.h"
#include "includes/element.h"
#include "utilities/math_utils.h"
// Application includes
#include "custom_utilities/element_utilities.hpp"
#include "integration_coefficients_calculator.hpp"
#include "stress_state_policy.h"
namespace Kratos
{
class KRATOS_API(GEO_MECHANICS_APPLICATION) UPwBaseElement : public Element
{
public:
using SizeType = std::size_t;
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(UPwBaseElement);
using Element::Element;
/// Constructor using an array of nodes
UPwBaseElement(IndexType NewId,
const NodesArrayType& ThisNodes,
std::unique_ptr<StressStatePolicy> pStressStatePolicy,
std::unique_ptr<IntegrationCoefficientModifier> pCoefficientModifier = nullptr);
UPwBaseElement(IndexType NewId,
GeometryType::Pointer pGeometry,
std::unique_ptr<StressStatePolicy> pStressStatePolicy,
std::unique_ptr<IntegrationCoefficientModifier> pCoefficientModifier = nullptr);
UPwBaseElement(IndexType NewId,
GeometryType::Pointer pGeometry,
PropertiesType::Pointer pProperties,
std::unique_ptr<StressStatePolicy> pStressStatePolicy,
std::unique_ptr<IntegrationCoefficientModifier> pCoefficientModifier = nullptr);
~UPwBaseElement() override = default;
UPwBaseElement(const UPwBaseElement&) = delete;
UPwBaseElement& operator=(const UPwBaseElement&) = delete;
UPwBaseElement(UPwBaseElement&&) noexcept = delete;
UPwBaseElement& operator=(UPwBaseElement&&) noexcept = delete;
int Check(const ProcessInfo& rCurrentProcessInfo) const override;
void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo&) const override;
void ResetConstitutiveLaw() override;
GeometryData::IntegrationMethod GetIntegrationMethod() const override;
void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,
VectorType& rRightHandSideVector,
const ProcessInfo& rCurrentProcessInfo) override;
void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix, const ProcessInfo& rCurrentProcessInfo) override;
void CalculateRightHandSide(VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo&) const override;
void CalculateMassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo) override;
void CalculateDampingMatrix(MatrixType& rDampingMatrix, const ProcessInfo& rCurrentProcessInfo) override;
void GetValuesVector(Vector& rValues, int Step = 0) const override;
void GetFirstDerivativesVector(Vector& rValues, int Step = 0) const override;
void GetSecondDerivativesVector(Vector& rValues, int Step = 0) const override;
void CalculateOnIntegrationPoints(const Variable<ConstitutiveLaw::Pointer>& rVariable,
std::vector<ConstitutiveLaw::Pointer>& rValues,
const ProcessInfo& rCurrentProcessInfo) override;
void CalculateOnIntegrationPoints(const Variable<array_1d<double, 3>>& rVariable,
std::vector<array_1d<double, 3>>& rValues,
const ProcessInfo& rCurrentProcessInfo) override;
void CalculateOnIntegrationPoints(const Variable<Matrix>& rVariable,
std::vector<Matrix>& rValues,
const ProcessInfo& rCurrentProcessInfo) override;
void CalculateOnIntegrationPoints(const Variable<double>& rVariable,
std::vector<double>& rValues,
const ProcessInfo& rCurrentProcessInfo) override;
void CalculateOnIntegrationPoints(const Variable<Vector>& rVariable,
std::vector<Vector>& rValues,
const ProcessInfo& rCurrentProcessInfo) override;
using Element::CalculateOnIntegrationPoints;
void SetValuesOnIntegrationPoints(const Variable<double>& rVariable,
const std::vector<double>& rValues,
const ProcessInfo& rCurrentProcessInfo) override;
void SetValuesOnIntegrationPoints(const Variable<Vector>& rVariable,
const std::vector<Vector>& rValues,
const ProcessInfo& rCurrentProcessInfo) override;
void SetValuesOnIntegrationPoints(const Variable<Matrix>& rVariable,
const std::vector<Matrix>& rValues,
const ProcessInfo& rCurrentProcessInfo) override;
using Element::SetValuesOnIntegrationPoints;
std::string Info() const override;
void PrintInfo(std::ostream& rOStream) const override;
protected:
/// Member Variables
GeometryData::IntegrationMethod mThisIntegrationMethod;
std::vector<ConstitutiveLaw::Pointer> mConstitutiveLawVector;
std::vector<RetentionLaw::Pointer> mRetentionLawVector;
std::vector<Vector> mStressVector;
std::vector<Vector> mStateVariablesFinalized;
virtual void CalculateMaterialStiffnessMatrix(MatrixType& rStiffnessMatrix, const ProcessInfo& CurrentProcessInfo);
virtual void CalculateAll(MatrixType& rLeftHandSideMatrix,
VectorType& rRightHandSideVector,
const ProcessInfo& CurrentProcessInfo,
bool CalculateStiffnessMatrixFlag,
bool CalculateResidualVectorFlag);
std::vector<double> CalculateIntegrationCoefficients(const GeometryType::IntegrationPointsArrayType& rIntegrationPoints,
const Vector& rDetJs) const;
void CalculateDerivativesOnInitialConfiguration(double& rDetJ,
Matrix& rJ0,
Matrix& rInvJ0,
Matrix& rDNu_DX0,
unsigned int IntegrationPointIndex) const;
void CalculateJacobianOnCurrentConfiguration(double& detJ, Matrix& rJ, Matrix& rInvJ, unsigned int GPoint) const;
virtual std::size_t GetNumberOfDOF() const;
StressStatePolicy& GetStressStatePolicy() const;
std::unique_ptr<IntegrationCoefficientModifier> CloneIntegrationCoefficientModifier() const;
private:
[[nodiscard]] virtual DofsVectorType GetDofs() const;
friend class Serializer;
void save(Serializer& rSerializer) const override;
void load(Serializer& rSerializer) override;
std::unique_ptr<StressStatePolicy> mpStressStatePolicy;
IntegrationCoefficientsCalculator mIntegrationCoefficientsCalculator;
};
// Class UPwBaseElement
} // namespace Kratos