-
Notifications
You must be signed in to change notification settings - Fork 258
/
Copy pathtransient_thermal_element.h
327 lines (279 loc) · 13.9 KB
/
transient_thermal_element.h
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
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
// KRATOS___
// // ) )
// // ___ ___
// // ____ //___) ) // ) )
// // / / // // / /
// ((____/ / ((____ ((___/ / MECHANICS
//
// License: geo_mechanics_application/license.txt
//
// Main authors: Mohamed Nabi
// John van Esch
// Gennady Markelov
//
#pragma once
#include "custom_constitutive/thermal_dispersion_law.h"
#include "custom_constitutive/thermal_filter_law.h"
#include "custom_retention/retention_law_factory.h"
#include "custom_utilities/check_utilities.h"
#include "custom_utilities/dof_utilities.h"
#include "geo_mechanics_application_variables.h"
#include "includes/element.h"
#include "includes/serializer.h"
namespace Kratos
{
template <unsigned int TDim, unsigned int TNumNodes>
class KRATOS_API(GEO_MECHANICS_APPLICATION) TransientThermalElement : public Element
{
public:
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(TransientThermalElement);
explicit TransientThermalElement(IndexType NewId = 0) : Element(NewId) {}
TransientThermalElement(IndexType NewId, GeometryType::Pointer pGeometry)
: Element(NewId, pGeometry)
{
}
TransientThermalElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
: Element(NewId, pGeometry, pProperties)
{
}
Element::Pointer Create(IndexType NewId, const NodesArrayType& rThisNodes, PropertiesType::Pointer pProperties) const override
{
return make_intrusive<TransientThermalElement>(NewId, GetGeometry().Create(rThisNodes), pProperties);
}
Element::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override
{
return make_intrusive<TransientThermalElement>(NewId, pGeom, pProperties);
}
void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo&) const override
{
rElementalDofList = GetDofs();
}
void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo&) const override
{
rResult = Geo::DofUtilities::ExtractEquationIdsFrom(GetDofs());
}
void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,
VectorType& rRightHandSideVector,
const ProcessInfo& rCurrentProcessInfo) override
{
KRATOS_TRY
GeometryType::ShapeFunctionsGradientsType dN_dX_container;
Vector det_J_container;
// ShapreFunctionsIntegrationsPointsGradients does not allow for the line element in 2D/3D
// configuration and will produce errors. To circumvent this, the dN_dX_container is
// separately computed with correct dimensions for the line element.
if (GetGeometry().LocalSpaceDimension() == 1) {
GetGeometry().DeterminantOfJacobian(det_J_container, this->GetIntegrationMethod());
dN_dX_container = GetGeometry().ShapeFunctionsLocalGradients(this->GetIntegrationMethod());
std::transform(dN_dX_container.begin(), dN_dX_container.end(), det_J_container.begin(),
dN_dX_container.begin(), std::divides<>());
} else {
GetGeometry().ShapeFunctionsIntegrationPointsGradients(dN_dX_container, det_J_container,
GetIntegrationMethod());
}
const auto integration_coefficients = CalculateIntegrationCoefficients(det_J_container);
const auto conductivity_matrix = CalculateConductivityMatrix(dN_dX_container, integration_coefficients);
const auto capacity_matrix = CalculateCapacityMatrix(integration_coefficients);
AddContributionsToLhsMatrix(rLeftHandSideMatrix, conductivity_matrix, capacity_matrix,
rCurrentProcessInfo[DT_TEMPERATURE_COEFFICIENT]);
AddContributionsToRhsVector(rRightHandSideVector, conductivity_matrix, capacity_matrix);
KRATOS_CATCH("")
}
GeometryData::IntegrationMethod GetIntegrationMethod() const override
{
switch (this->GetGeometry().GetGeometryOrderType()) {
case GeometryData::Kratos_Cubic_Order:
return this->GetGeometry().GetGeometryType() == GeometryData::KratosGeometryType::Kratos_Triangle2D10
? GeometryData::IntegrationMethod::GI_GAUSS_4
: GeometryData::IntegrationMethod::GI_GAUSS_3;
case GeometryData::Kratos_Quartic_Order:
return GeometryData::IntegrationMethod::GI_GAUSS_5;
default:
return GeometryData::IntegrationMethod::GI_GAUSS_2;
}
}
int Check(const ProcessInfo&) const override
{
KRATOS_TRY
CheckUtilities::CheckDomainSize(GetGeometry().DomainSize(), Id());
CheckHasSolutionStepsDataFor(TEMPERATURE);
CheckHasSolutionStepsDataFor(DT_TEMPERATURE);
CheckHasDofsFor(TEMPERATURE);
CheckProperties();
CheckForNonZeroZCoordinateIn2D();
KRATOS_CATCH("")
return 0;
}
private:
void CheckHasSolutionStepsDataFor(const Variable<double>& rVariable) const
{
for (const auto& node : GetGeometry()) {
KRATOS_ERROR_IF_NOT(node.SolutionStepsDataHas(rVariable))
<< "Missing variable " << rVariable.Name() << " on node " << node.Id() << std::endl;
}
}
void CheckHasDofsFor(const Variable<double>& rVariable) const
{
for (const auto& node : GetGeometry()) {
KRATOS_ERROR_IF_NOT(node.HasDofFor(rVariable))
<< "Missing degree of freedom for " << rVariable.Name() << " on node " << node.Id()
<< std::endl;
}
}
void CheckProperties() const
{
CheckProperty(DENSITY_WATER);
CheckProperty(POROSITY);
CheckProperty(RETENTION_LAW, "SaturatedLaw");
CheckProperty(SATURATED_SATURATION);
CheckProperty(DENSITY_SOLID);
CheckProperty(SPECIFIC_HEAT_CAPACITY_WATER);
CheckProperty(SPECIFIC_HEAT_CAPACITY_SOLID);
CheckProperty(THERMAL_CONDUCTIVITY_WATER);
CheckProperty(THERMAL_CONDUCTIVITY_SOLID_XX);
CheckProperty(THERMAL_CONDUCTIVITY_SOLID_YY);
CheckProperty(THERMAL_CONDUCTIVITY_SOLID_XY);
if constexpr (TDim == 3) {
CheckProperty(THERMAL_CONDUCTIVITY_SOLID_ZZ);
CheckProperty(THERMAL_CONDUCTIVITY_SOLID_YZ);
CheckProperty(THERMAL_CONDUCTIVITY_SOLID_XZ);
}
}
void CheckProperty(const Kratos::Variable<double>& rVariable) const
{
KRATOS_ERROR_IF_NOT(GetProperties().Has(rVariable))
<< rVariable.Name() << " does not exist in the thermal element's properties" << std::endl;
KRATOS_ERROR_IF(GetProperties()[rVariable] < 0.0)
<< rVariable.Name() << " has an invalid value at element " << Id() << std::endl;
}
void CheckProperty(const Kratos::Variable<std::string>& rVariable, const std::string& rName) const
{
KRATOS_ERROR_IF_NOT(GetProperties().Has(rVariable))
<< rVariable.Name() << " does not exist in the thermal element's properties" << std::endl;
KRATOS_ERROR_IF_NOT(GetProperties()[rVariable] == rName)
<< rVariable.Name() << " has a value of (" << GetProperties()[rVariable]
<< ") instead of (" << rName << ") at element " << Id() << std::endl;
}
void CheckForNonZeroZCoordinateIn2D() const
{
if constexpr (TDim == 2) {
const auto& r_geometry = GetGeometry();
auto pos = std::find_if(r_geometry.begin(), r_geometry.end(),
[](const auto& node) { return node.Z() != 0.0; });
KRATOS_ERROR_IF_NOT(pos == r_geometry.end())
<< " Node with non-zero Z coordinate found. Id: " << pos->Id() << std::endl;
}
}
static void AddContributionsToLhsMatrix(MatrixType& rLeftHandSideMatrix,
const BoundedMatrix<double, TNumNodes, TNumNodes>& rConductivityMatrix,
const BoundedMatrix<double, TNumNodes, TNumNodes>& rCapacityMatrix,
double DtTemperatureCoefficient)
{
rLeftHandSideMatrix = rConductivityMatrix;
rLeftHandSideMatrix += (DtTemperatureCoefficient * rCapacityMatrix);
}
void AddContributionsToRhsVector(VectorType& rRightHandSideVector,
const BoundedMatrix<double, TNumNodes, TNumNodes>& rConductivityMatrix,
const BoundedMatrix<double, TNumNodes, TNumNodes>& rCapacityMatrix) const
{
const auto capacity_vector =
array_1d<double, TNumNodes>{-prod(rCapacityMatrix, GetNodalValuesOf(DT_TEMPERATURE))};
rRightHandSideVector = capacity_vector;
const auto conductivity_vector =
array_1d<double, TNumNodes>{-prod(rConductivityMatrix, GetNodalValuesOf(TEMPERATURE))};
rRightHandSideVector += conductivity_vector;
}
Vector CalculateIntegrationCoefficients(const Vector& rDetJContainer) const
{
const auto& r_properties = GetProperties();
const auto& r_integration_points = GetGeometry().IntegrationPoints(GetIntegrationMethod());
auto result = Vector{r_integration_points.size()};
for (unsigned int integration_point_index = 0;
integration_point_index < r_integration_points.size(); ++integration_point_index) {
result[integration_point_index] = r_integration_points[integration_point_index].Weight() *
rDetJContainer[integration_point_index];
if (GetGeometry().LocalSpaceDimension() == 1) {
result[integration_point_index] *= r_properties[CROSS_AREA];
}
}
return result;
}
BoundedMatrix<double, TNumNodes, TNumNodes> CalculateConductivityMatrix(
const GeometryType::ShapeFunctionsGradientsType& rShapeFunctionGradients,
const Vector& rIntegrationCoefficients) const
{
const auto law = CreateThermalLaw();
const auto constitutive_matrix = law->CalculateThermalDispersionMatrix(GetProperties());
auto result = BoundedMatrix<double, TNumNodes, TNumNodes>{ZeroMatrix{TNumNodes, TNumNodes}};
for (unsigned int integration_point_index = 0;
integration_point_index < GetGeometry().IntegrationPointsNumber(GetIntegrationMethod());
++integration_point_index) {
BoundedMatrix<double, TDim, TNumNodes> Temp =
prod(constitutive_matrix, trans(rShapeFunctionGradients[integration_point_index]));
result += prod(rShapeFunctionGradients[integration_point_index], Temp) *
rIntegrationCoefficients[integration_point_index];
}
return result;
}
BoundedMatrix<double, TNumNodes, TNumNodes> CalculateCapacityMatrix(const Vector& rIntegrationCoefficients) const
{
const auto& r_properties = GetProperties();
RetentionLaw::Parameters parameters(r_properties);
auto retention_law = RetentionLawFactory::Clone(r_properties);
const double saturation = retention_law->CalculateSaturation(parameters);
const auto c_water = r_properties[POROSITY] * saturation * r_properties[DENSITY_WATER] *
r_properties[SPECIFIC_HEAT_CAPACITY_WATER];
const auto c_solid = (1.0 - r_properties[POROSITY]) * r_properties[DENSITY_SOLID] *
r_properties[SPECIFIC_HEAT_CAPACITY_SOLID];
auto result = BoundedMatrix<double, TNumNodes, TNumNodes>{ZeroMatrix{TNumNodes, TNumNodes}};
const auto& r_N_container = GetGeometry().ShapeFunctionsValues(GetIntegrationMethod());
for (unsigned int integration_point_index = 0;
integration_point_index < GetGeometry().IntegrationPointsNumber(GetIntegrationMethod());
++integration_point_index) {
const auto N = Vector{row(r_N_container, integration_point_index)};
result += (c_water + c_solid) * outer_prod(N, N) * rIntegrationCoefficients[integration_point_index];
}
return result;
}
array_1d<double, TNumNodes> GetNodalValuesOf(const Variable<double>& rNodalVariable) const
{
auto result = array_1d<double, TNumNodes>{};
const auto& r_geometry = GetGeometry();
std::transform(r_geometry.begin(), r_geometry.end(), result.begin(), [&rNodalVariable](const auto& node) {
return node.FastGetSolutionStepValue(rNodalVariable);
});
return result;
}
std::unique_ptr<GeoThermalLaw> CreateThermalLaw() const
{
const std::size_t number_of_dimensions = GetGeometry().LocalSpaceDimension();
std::unique_ptr<GeoThermalLaw> law;
if (GetProperties().Has(THERMAL_LAW_NAME)) {
const std::string& rThermalLawName = GetProperties()[THERMAL_LAW_NAME];
if (rThermalLawName == "GeoThermalDispersionLaw") {
law = std::make_unique<GeoThermalDispersionLaw>(number_of_dimensions);
} else if (rThermalLawName == "GeoThermalFilterLaw") {
law = std::make_unique<GeoThermalFilterLaw>();
} else {
KRATOS_ERROR << "Undefined THERMAL_LAW_NAME! " << rThermalLawName << std::endl;
}
} else {
law = std::make_unique<GeoThermalDispersionLaw>(number_of_dimensions);
}
return law;
}
[[nodiscard]] DofsVectorType GetDofs() const
{
return Geo::DofUtilities::ExtractDofsFromNodes(GetGeometry(), TEMPERATURE);
}
friend class Serializer;
void save(Serializer& rSerializer) const override
{
KRATOS_SERIALIZE_SAVE_BASE_CLASS(rSerializer, Element)
}
void load(Serializer& rSerializer) override
{
KRATOS_SERIALIZE_LOAD_BASE_CLASS(rSerializer, Element)
}
};
} // namespace Kratos