/
SmallDeformationFEM.h
467 lines (391 loc) · 17.2 KB
/
SmallDeformationFEM.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
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
/**
* \copyright
* Copyright (c) 2012-2018, OpenGeoSys Community (http://www.opengeosys.org)
* Distributed under a Modified BSD License.
* See accompanying file LICENSE.txt or
* http://www.opengeosys.org/project/license
*
*/
#pragma once
#include <memory>
#include <vector>
#include "MaterialLib/SolidModels/LinearElasticIsotropic.h"
#include "MathLib/LinAlg/Eigen/EigenMapTools.h"
#include "NumLib/Extrapolation/ExtrapolatableElement.h"
#include "NumLib/Fem/FiniteElement/TemplateIsoparametric.h"
#include "NumLib/Fem/ShapeMatrixPolicy.h"
#include "ProcessLib/Deformation/BMatrixPolicy.h"
#include "ProcessLib/Deformation/GMatrixPolicy.h"
#include "ProcessLib/Deformation/LinearBMatrix.h"
#include "ProcessLib/LocalAssemblerInterface.h"
#include "ProcessLib/LocalAssemblerTraits.h"
#include "ProcessLib/Parameter/Parameter.h"
#include "ProcessLib/Utils/InitShapeMatrices.h"
#include "LocalAssemblerInterface.h"
#include "SmallDeformationProcessData.h"
namespace ProcessLib
{
namespace SmallDeformation
{
template <typename BMatricesType, typename ShapeMatricesType,
int DisplacementDim>
struct IntegrationPointData final
{
explicit IntegrationPointData(
MaterialLib::Solids::MechanicsBase<DisplacementDim>& solid_material)
: solid_material(solid_material),
material_state_variables(
solid_material.createMaterialStateVariables())
{
}
typename BMatricesType::KelvinVectorType sigma, sigma_prev;
typename BMatricesType::KelvinVectorType eps, eps_prev;
double free_energy_density = 0;
MaterialLib::Solids::MechanicsBase<DisplacementDim>& solid_material;
std::unique_ptr<typename MaterialLib::Solids::MechanicsBase<
DisplacementDim>::MaterialStateVariables>
material_state_variables;
double integration_weight;
typename ShapeMatricesType::NodalRowVectorType N;
typename ShapeMatricesType::GlobalDimNodalMatrixType dNdx;
void pushBackState()
{
eps_prev = eps;
sigma_prev = sigma;
material_state_variables->pushBackState();
}
EIGEN_MAKE_ALIGNED_OPERATOR_NEW;
};
/// Used by for extrapolation of the integration point values. It is ordered
/// (and stored) by integration points.
template <typename ShapeMatrixType>
struct SecondaryData
{
std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N;
};
template <typename ShapeFunction, typename IntegrationMethod,
int DisplacementDim>
class SmallDeformationLocalAssembler
: public SmallDeformationLocalAssemblerInterface<DisplacementDim>
{
public:
using ShapeMatricesType =
ShapeMatrixPolicyType<ShapeFunction, DisplacementDim>;
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType;
using NodalVectorType = typename ShapeMatricesType::NodalVectorType;
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices;
using BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>;
using BMatrixType = typename BMatricesType::BMatrixType;
using StiffnessMatrixType = typename BMatricesType::StiffnessMatrixType;
using NodalForceVectorType = typename BMatricesType::NodalForceVectorType;
using NodalDisplacementVectorType =
typename BMatricesType::NodalForceVectorType;
using GMatricesType = GMatrixPolicyType<ShapeFunction, DisplacementDim>;
using GradientVectorType = typename GMatricesType::GradientVectorType;
using GradientMatrixType = typename GMatricesType::GradientMatrixType;
SmallDeformationLocalAssembler(SmallDeformationLocalAssembler const&) =
delete;
SmallDeformationLocalAssembler(SmallDeformationLocalAssembler&&) = delete;
SmallDeformationLocalAssembler(
MeshLib::Element const& e,
std::size_t const /*local_matrix_size*/,
bool const is_axially_symmetric,
unsigned const integration_order,
SmallDeformationProcessData<DisplacementDim>& process_data)
: _process_data(process_data),
_integration_method(integration_order),
_element(e),
_is_axially_symmetric(is_axially_symmetric)
{
unsigned const n_integration_points =
_integration_method.getNumberOfPoints();
_ip_data.reserve(n_integration_points);
_secondary_data.N.resize(n_integration_points);
auto const shape_matrices =
initShapeMatrices<ShapeFunction, ShapeMatricesType,
IntegrationMethod, DisplacementDim>(
e, is_axially_symmetric, _integration_method);
for (unsigned ip = 0; ip < n_integration_points; ip++)
{
_ip_data.emplace_back(*_process_data.material);
auto& ip_data = _ip_data[ip];
auto const& sm = shape_matrices[ip];
_ip_data[ip].integration_weight =
_integration_method.getWeightedPoint(ip).getWeight() *
sm.integralMeasure * sm.detJ;
ip_data.N = sm.N;
ip_data.dNdx = sm.dNdx;
static const int kelvin_vector_size =
MathLib::KelvinVector::KelvinVectorDimensions<
DisplacementDim>::value;
// Initialize current time step values
ip_data.sigma.setZero(kelvin_vector_size);
ip_data.eps.setZero(kelvin_vector_size);
// Previous time step values are not initialized and are set later.
ip_data.sigma_prev.resize(kelvin_vector_size);
ip_data.eps_prev.resize(kelvin_vector_size);
_secondary_data.N[ip] = shape_matrices[ip].N;
}
}
void setIPDataInitialConditions(std::string const& name,
double const* values) override
{
if (name == "sigma_ip")
{
setSigma(values);
}
}
void assemble(double const /*t*/, std::vector<double> const& /*local_x*/,
std::vector<double>& /*local_M_data*/,
std::vector<double>& /*local_K_data*/,
std::vector<double>& /*local_b_data*/) override
{
OGS_FATAL(
"SmallDeformationLocalAssembler: assembly without jacobian is not "
"implemented.");
}
void assembleWithJacobian(double const t,
std::vector<double> const& local_x,
std::vector<double> const& /*local_xdot*/,
const double /*dxdot_dx*/, const double /*dx_dx*/,
std::vector<double>& /*local_M_data*/,
std::vector<double>& /*local_K_data*/,
std::vector<double>& local_b_data,
std::vector<double>& local_Jac_data) override
{
auto const local_matrix_size = local_x.size();
auto local_Jac = MathLib::createZeroedMatrix<StiffnessMatrixType>(
local_Jac_data, local_matrix_size, local_matrix_size);
auto local_b = MathLib::createZeroedVector<NodalDisplacementVectorType>(
local_b_data, local_matrix_size);
unsigned const n_integration_points =
_integration_method.getNumberOfPoints();
SpatialPosition x_position;
x_position.setElementID(_element.getID());
for (unsigned ip = 0; ip < n_integration_points; ip++)
{
x_position.setIntegrationPoint(ip);
auto const& w = _ip_data[ip].integration_weight;
auto const& N = _ip_data[ip].N;
auto const& dNdx = _ip_data[ip].dNdx;
typename ShapeMatricesType::template MatrixType<DisplacementDim,
displacement_size>
N_u_op = ShapeMatricesType::template MatrixType<
DisplacementDim,
displacement_size>::Zero(DisplacementDim,
displacement_size);
for (int i = 0; i < DisplacementDim; ++i)
N_u_op
.template block<1, displacement_size / DisplacementDim>(
i, i * displacement_size / DisplacementDim)
.noalias() = N;
auto const x_coord =
interpolateXCoordinate<ShapeFunction, ShapeMatricesType>(
_element, N);
auto const B = LinearBMatrix::computeBMatrix<
DisplacementDim, ShapeFunction::NPOINTS,
typename BMatricesType::BMatrixType>(dNdx, N, x_coord,
_is_axially_symmetric);
auto const& eps_prev = _ip_data[ip].eps_prev;
auto const& sigma_prev = _ip_data[ip].sigma_prev;
auto& eps = _ip_data[ip].eps;
auto& sigma = _ip_data[ip].sigma;
auto& state = _ip_data[ip].material_state_variables;
eps.noalias() =
B *
Eigen::Map<typename BMatricesType::NodalForceVectorType const>(
local_x.data(), ShapeFunction::NPOINTS * DisplacementDim);
auto&& solution = _ip_data[ip].solid_material.integrateStress(
t, x_position, _process_data.dt, eps_prev, eps, sigma_prev,
*state);
if (!solution)
OGS_FATAL("Computation of local constitutive relation failed.");
MathLib::KelvinVector::KelvinMatrixType<DisplacementDim> C;
std::tie(sigma, state, C) = std::move(*solution);
auto const rho = _process_data.solid_density(t, x_position)[0];
auto const& b = _process_data.specific_body_force;
local_b.noalias() -=
(B.transpose() * sigma - N_u_op.transpose() * rho * b) * w;
local_Jac.noalias() += B.transpose() * C * B * w;
}
}
void preTimestepConcrete(std::vector<double> const& /*local_x*/,
double const /*t*/,
double const /*delta_t*/) override
{
unsigned const n_integration_points =
_integration_method.getNumberOfPoints();
for (unsigned ip = 0; ip < n_integration_points; ip++)
{
_ip_data[ip].pushBackState();
}
}
void postTimestepConcrete(std::vector<double> const& /*local_x*/) override
{
unsigned const n_integration_points =
_integration_method.getNumberOfPoints();
SpatialPosition x_position;
x_position.setElementID(_element.getID());
for (unsigned ip = 0; ip < n_integration_points; ip++)
{
x_position.setIntegrationPoint(ip);
auto& ip_data = _ip_data[ip];
// Update free energy density needed for material forces.
ip_data.free_energy_density =
ip_data.solid_material.computeFreeEnergyDensity(
_process_data.t, x_position, _process_data.dt, ip_data.eps,
ip_data.sigma, *ip_data.material_state_variables);
}
}
std::vector<double> const& getMaterialForces(
std::vector<double> const& local_x,
std::vector<double>& nodal_values) override
{
return ProcessLib::SmallDeformation::getMaterialForces<
DisplacementDim, ShapeFunction, ShapeMatricesType,
typename BMatricesType::NodalForceVectorType,
NodalDisplacementVectorType, GradientVectorType,
GradientMatrixType>(local_x, nodal_values, _integration_method,
_ip_data, _element, _is_axially_symmetric);
}
Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
const unsigned integration_point) const override
{
auto const& N = _secondary_data.N[integration_point];
// assumes N is stored contiguously in memory
return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
}
std::vector<double> const& getIntPtFreeEnergyDensity(
const double /*t*/,
GlobalVector const& /*current_solution*/,
NumLib::LocalToGlobalIndexMap const& /*dof_table*/,
std::vector<double>& cache) const override
{
cache.clear();
cache.reserve(_ip_data.size());
for (auto const& ip_data : _ip_data)
{
cache.push_back(ip_data.free_energy_density);
}
return cache;
}
void setSigma(double const* values)
{
auto const kelvin_vector_size =
MathLib::KelvinVector::KelvinVectorDimensions<
DisplacementDim>::value;
auto const n_integration_points = _ip_data.size();
std::vector<double> ip_sigma_values;
auto sigma_values =
Eigen::Map<Eigen::Matrix<double, kelvin_vector_size, Eigen::Dynamic,
Eigen::ColMajor> const>(
values, kelvin_vector_size, n_integration_points);
for (unsigned ip = 0; ip < n_integration_points; ++ip)
{
_ip_data[ip].sigma =
MathLib::KelvinVector::symmetricTensorToKelvinVector(
sigma_values.col(ip));
}
}
std::vector<double> getSigma() const override
{
using KelvinVectorType = typename BMatricesType::KelvinVectorType;
auto const kelvin_vector_size =
MathLib::KelvinVector::KelvinVectorDimensions<
DisplacementDim>::value;
auto const n_integration_points = _ip_data.size();
std::vector<double> ip_sigma_values;
//ip_sigma_values.resize(kelvin_vector_size * n_integration_points);
auto cache_mat = MathLib::createZeroedMatrix<Eigen::Matrix<
double, Eigen::Dynamic, kelvin_vector_size, Eigen::RowMajor>>(
ip_sigma_values, n_integration_points, kelvin_vector_size);
// TODO make a general implementation for converting KelvinVectors
// back to symmetric rank-2 tensors.
for (unsigned ip = 0; ip < n_integration_points; ++ip)
{
auto const& sigma = _ip_data[ip].sigma;
for (typename KelvinVectorType::Index component = 0;
component < kelvin_vector_size && component < 3;
++component)
{ // xx, yy, zz components
cache_mat(ip, component) = sigma[component];
}
for (typename KelvinVectorType::Index component = 3;
component < kelvin_vector_size;
++component)
{ // mixed xy, yz, xz components
cache_mat(ip, component) = sigma[component] / std::sqrt(2);
}
}
return ip_sigma_values;
}
std::vector<double> const& getIntPtSigma(
const double /*t*/,
GlobalVector const& /*current_solution*/,
NumLib::LocalToGlobalIndexMap const& /*dof_table*/,
std::vector<double>& cache) const override
{
static const int kelvin_vector_size =
MathLib::KelvinVector::KelvinVectorDimensions<
DisplacementDim>::value;
auto const num_intpts = _ip_data.size();
cache.clear();
auto cache_mat = MathLib::createZeroedMatrix<Eigen::Matrix<
double, kelvin_vector_size, Eigen::Dynamic, Eigen::RowMajor>>(
cache, kelvin_vector_size, num_intpts);
for (unsigned ip = 0; ip < num_intpts; ++ip)
{
auto const& sigma = _ip_data[ip].sigma;
cache_mat.col(ip) =
MathLib::KelvinVector::kelvinVectorToSymmetricTensor(sigma);
}
return cache;
}
virtual std::vector<double> const& getIntPtEpsilon(
const double /*t*/,
GlobalVector const& /*current_solution*/,
NumLib::LocalToGlobalIndexMap const& /*dof_table*/,
std::vector<double>& cache) const override
{
auto const kelvin_vector_size =
MathLib::KelvinVector::KelvinVectorDimensions<
DisplacementDim>::value;
auto const num_intpts = _ip_data.size();
cache.clear();
auto cache_mat = MathLib::createZeroedMatrix<Eigen::Matrix<
double, kelvin_vector_size, Eigen::Dynamic, Eigen::RowMajor>>(
cache, kelvin_vector_size, num_intpts);
for (unsigned ip = 0; ip < num_intpts; ++ip)
{
auto const& eps = _ip_data[ip].eps;
cache_mat.col(ip) =
MathLib::KelvinVector::kelvinVectorToSymmetricTensor(eps);
}
return cache;
}
unsigned getNumberOfIntegrationPoints() const override
{
return _integration_method.getNumberOfPoints();
}
typename MaterialLib::Solids::MechanicsBase<
DisplacementDim>::MaterialStateVariables const&
getMaterialStateVariablesAt(unsigned integration_point) const override
{
return *_ip_data[integration_point].material_state_variables;
}
private:
SmallDeformationProcessData<DisplacementDim>& _process_data;
std::vector<
IntegrationPointData<BMatricesType, ShapeMatricesType, DisplacementDim>,
Eigen::aligned_allocator<IntegrationPointData<
BMatricesType, ShapeMatricesType, DisplacementDim>>>
_ip_data;
IntegrationMethod _integration_method;
MeshLib::Element const& _element;
SecondaryData<typename ShapeMatrices::ShapeType> _secondary_data;
bool const _is_axially_symmetric;
static const int displacement_size =
ShapeFunction::NPOINTS * DisplacementDim;
};
} // namespace SmallDeformation
} // namespace ProcessLib