OGS
ThermoHydroMechanicsFEM.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <memory>
7#include <vector>
8
17#include "NumLib/DOF/LocalDOF.h"
28
29namespace ProcessLib
30{
32{
33namespace MPL = MaterialPropertyLib;
34
37template <typename ShapeMatrixType>
39{
40 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N_u;
41};
42
43template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
44 int DisplacementDim>
46 : public LocalAssemblerInterface<DisplacementDim>
47{
48public:
51
52 // Types for pressure.
55
58
61
62 static int const KelvinVectorSize =
65
66 using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>;
67
68 static constexpr auto& N_u_op = MathLib::eigenBlockMatrixView<
69 DisplacementDim,
71
75 delete;
76
78 MeshLib::Element const& e,
79 std::size_t const /*local_matrix_size*/,
80 NumLib::GenericIntegrationMethod const& integration_method,
81 bool const is_axially_symmetric,
83
86 std::string_view const name,
87 double const* values,
88 int const integration_order) override;
89
90 void assemble(double const /*t*/, double const /*dt*/,
91 std::vector<double> const& /*local_x*/,
92 std::vector<double> const& /*local_x_prev*/,
93 std::vector<double>& /*local_M_data*/,
94 std::vector<double>& /*local_K_data*/,
95 std::vector<double>& /*local_rhs_data*/) override
96 {
98 "ThermoHydroMechanicsLocalAssembler: assembly without Jacobian is "
99 "not implemented.");
100 }
101
102 void assembleWithJacobian(double const t, double const dt,
103 std::vector<double> const& local_x,
104 std::vector<double> const& local_x_prev,
105 std::vector<double>& local_rhs_data,
106 std::vector<double>& local_Jac_data) override;
107
108 void initializeConcrete() override
109 {
110 unsigned const n_integration_points =
111 _integration_method.getNumberOfPoints();
112
113 for (unsigned ip = 0; ip < n_integration_points; ip++)
114 {
115 auto& ip_data = _ip_data[ip];
116
117 ParameterLib::SpatialPosition const x_position{
118 std::nullopt, _element.getID(),
121 ShapeFunctionDisplacement,
123
125 if (_process_data.initial_stress.value)
126 {
127 ip_data.sigma_eff =
129 DisplacementDim>((*_process_data.initial_stress.value)(
130 std::numeric_limits<double>::quiet_NaN() /* time
131 independent
132 */
133 ,
134 x_position));
135 }
136
137 double const t = 0; // TODO (naumov) pass t from top
138 ip_data.solid_material.initializeInternalStateVariables(
139 t, x_position, *ip_data.material_state_variables);
140
141 ip_data.pushBackState();
142 }
143 }
144
145 void setInitialConditionsConcrete(Eigen::VectorXd const local_x,
146 double const t,
147 int const process_id) override;
148
149 void preTimestepConcrete(std::vector<double> const& /*local_x*/,
150 double const /*t*/, double const /*dt*/) override
151 {
152 unsigned const n_integration_points =
153 _integration_method.getNumberOfPoints();
154
155 for (unsigned ip = 0; ip < n_integration_points; ip++)
156 {
157 _ip_data_output[ip].velocity.setConstant(
158 DisplacementDim, std::numeric_limits<double>::quiet_NaN());
159 }
160 }
161
162 void postTimestepConcrete(Eigen::VectorXd const& local_x,
163 Eigen::VectorXd const& local_x_prev,
164 double const t, double const dt,
165 int const /*process_id*/) override
166 {
167 unsigned const n_integration_points =
168 _integration_method.getNumberOfPoints();
169
170 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
171
172 for (unsigned ip = 0; ip < n_integration_points; ip++)
173 {
174 auto& ip_data = _ip_data[ip];
175 auto const& N_u = ip_data.N_u;
176 auto const& dNdx_u = ip_data.dNdx_u;
177
178 ParameterLib::SpatialPosition const x_position{
179 std::nullopt, _element.getID(),
182 ShapeFunctionDisplacement,
184
185 updateConstitutiveRelations(local_x, local_x_prev, x_position, t,
186 dt, _ip_data[ip], _ip_data_output[ip]);
187
188 auto const x_coord =
189 x_position.getCoordinates().value()[0]; // r for axisymmetry
190 auto const B = LinearBMatrix::computeBMatrix<
191 DisplacementDim, ShapeFunctionDisplacement::NPOINTS,
192 typename BMatricesType::BMatrixType>(dNdx_u, N_u, x_coord,
194
196 eps_prev = B * u_prev;
197
198 _ip_data[ip].eps0 =
199 _ip_data[ip].eps0_prev +
200 (1 - _ip_data[ip].phi_fr_prev / _ip_data[ip].porosity) *
201 (eps_prev - _ip_data[ip].eps0_prev);
202 _ip_data[ip].pushBackState();
203 }
204 }
205
207 double const t, double const dt, Eigen::VectorXd const& local_x,
208 Eigen::VectorXd const& local_x_prev) override;
209
210 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
211 const unsigned integration_point) const override
212 {
213 auto const& N_u = _secondary_data.N_u[integration_point];
214
215 // assumes N is stored contiguously in memory
216 return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
217 }
218
219 std::vector<double> const& getIntPtDarcyVelocity(
220 const double t,
221 std::vector<GlobalVector*> const& x,
222 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
223 std::vector<double>& cache) const override;
224
225 // TODO (naumov) This method is same as getIntPtSigma but for arguments and
226 // the ordering of the cache_mat.
227 // There should be only one.
228 std::vector<double> getSigma() const override
229 {
230 constexpr int kelvin_vector_size =
232
234 [this](std::vector<double>& values)
235 { return getIntPtSigma(0, {}, {}, values); });
236 }
237
238 std::vector<double> getSigmaIce() const override
239 {
240 constexpr int kelvin_vector_size =
242
244 [this](std::vector<double>& values)
245 { return getIntPtSigmaIce(0, {}, {}, values); });
246 }
247
248 std::vector<double> getIceVolumeFraction() const override
249 {
250 std::vector<double> result;
251 getIntPtIceVolume(0, {}, {}, result);
252 return result;
253 }
254
255 std::vector<double> const& getIntPtFluidDensity(
256 const double t,
257 std::vector<GlobalVector*> const& x,
258 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
259 std::vector<double>& cache) const override;
260
261 std::vector<double> const& getIntPtViscosity(
262 const double t,
263 std::vector<GlobalVector*> const& x,
264 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
265 std::vector<double>& cache) const override;
266
268 {
269 return displacement_size;
270 }
271
272private:
275 using IpData =
277 ShapeMatricesTypePressure, DisplacementDim,
278 ShapeFunctionDisplacement::NPOINTS>;
279
280private:
282 Eigen::Ref<Eigen::VectorXd const> const local_x,
283 Eigen::Ref<Eigen::VectorXd const> const local_x_prev,
284 ParameterLib::SpatialPosition const& x_position, double const t,
285 double const dt, IpData& ip_data,
287
288 std::size_t setSigma(double const* values)
289 {
291 values, _ip_data, &IpData::sigma_eff);
292 }
293
294 std::vector<double> const& getIntPtSigma(
295 const double /*t*/,
296 std::vector<GlobalVector*> const& /*x*/,
297 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
298 std::vector<double>& cache) const override
299 {
301 _ip_data, &IpData::sigma_eff, cache);
302 }
303
304 std::vector<double> const& getIntPtSigmaIce(
305 const double /*t*/,
306 std::vector<GlobalVector*> const& /*x*/,
307 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
308 std::vector<double>& cache) const override
309 {
312 }
313
314 std::vector<double> getEpsilon0() const override
315 {
316 constexpr int kelvin_vector_size =
318
320 [this](std::vector<double>& values)
321 { return getIntPtEpsilon0(0, {}, {}, values); });
322 }
323
324 virtual std::vector<double> const& getIntPtEpsilon0(
325 const double /*t*/,
326 std::vector<GlobalVector*> const& /*x*/,
327 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
328 std::vector<double>& cache) const override
329 {
331 _ip_data, &IpData::eps0, cache);
332 }
333 std::vector<double> getEpsilonM() const override
334 {
335 constexpr int kelvin_vector_size =
337
339 [this](std::vector<double>& values)
340 { return getIntPtEpsilonM(0, {}, {}, values); });
341 }
342
343 virtual std::vector<double> const& getIntPtEpsilonM(
344 const double /*t*/,
345 std::vector<GlobalVector*> const& /*x*/,
346 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
347 std::vector<double>& cache) const override
348 {
350 _ip_data, &IpData::eps_m, cache);
351 }
352
353 std::vector<double> getEpsilon() const override
354 {
355 constexpr int kelvin_vector_size =
357
359 [this](std::vector<double>& values)
360 { return getIntPtEpsilon(0, {}, {}, values); });
361 }
362
363 virtual std::vector<double> const& getIntPtEpsilon(
364 const double /*t*/,
365 std::vector<GlobalVector*> const& /*x*/,
366 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
367 std::vector<double>& cache) const override
368 {
370 _ip_data, &IpData::eps, cache);
371 }
372
373 std::vector<double> const& getIntPtIceVolume(
374 const double /*t*/,
375 std::vector<GlobalVector*> const& /*x*/,
376 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
377 std::vector<double>& cache) const override
378 {
380 _ip_data, &IpData::phi_fr, cache);
381 }
382
383 unsigned getNumberOfIntegrationPoints() const override
384 {
385 return _integration_method.getNumberOfPoints();
386 }
387
388 int getMaterialID() const override
389 {
390 return _process_data.material_ids == nullptr
391 ? 0
392 : (*_process_data.material_ids)[_element.getID()];
393 }
394
396 std::function<std::span<double>(
398 MaterialStateVariables&)> const& get_values_span,
399 int const& n_components) const override
400 {
403 n_components);
404 }
405
407 DisplacementDim>::MaterialStateVariables const&
408 getMaterialStateVariablesAt(unsigned integration_point) const override
409 {
410 return *_ip_data[integration_point].material_state_variables;
411 }
412
413private:
414 template <typename SolutionVector>
415 static constexpr auto localDOF(SolutionVector const& x)
416 {
417 return NumLib::localDOF<
418 ShapeFunctionPressure, ShapeFunctionPressure,
420 }
421
423
424 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
425 std::vector<IntegrationPointDataForOutput<DisplacementDim>,
426 Eigen::aligned_allocator<
429
434 typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType>
436
437 // The shape function of pressure has the same form with the shape function
438 // of temperature
439 static const int temperature_index = 0;
440 static const int temperature_size = ShapeFunctionPressure::NPOINTS;
441 static const int pressure_index = ShapeFunctionPressure::NPOINTS;
442 static const int pressure_size = ShapeFunctionPressure::NPOINTS;
443 static const int displacement_index = ShapeFunctionPressure::NPOINTS * 2;
444 static const int displacement_size =
445 ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
446};
447
448} // namespace ThermoHydroMechanics
449} // namespace ProcessLib
450
#define OGS_FATAL(...)
Definition Error.h:19
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
std::optional< MathLib::Point3d > const getCoordinates() const
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
std::vector< double > const & getIntPtViscosity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
std::vector< double > const & getIntPtFluidDensity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
ThermoHydroMechanicsLocalAssembler(ThermoHydroMechanicsLocalAssembler const &)=delete
std::vector< double > getMaterialStateVariableInternalState(std::function< std::span< double >(typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables &)> const &get_values_span, int const &n_components) const override
std::vector< double > const & getIntPtSigmaIce(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
typename ShapeMatricesTypePressure::GlobalDimMatrixType GlobalDimMatrixType
ConstitutiveRelationsValues< DisplacementDim > updateConstitutiveRelations(Eigen::Ref< Eigen::VectorXd const > const local_x, Eigen::Ref< Eigen::VectorXd const > const local_x_prev, ParameterLib::SpatialPosition const &x_position, double const t, double const dt, IpData &ip_data, IntegrationPointDataForOutput< DisplacementDim > &ip_data_output) const
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
void preTimestepConcrete(std::vector< double > const &, double const, double const) override
BMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > BMatricesType
typename ShapeMatricesTypePressure::GlobalDimVectorType GlobalDimVectorType
virtual std::vector< double > const & getIntPtEpsilon(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) override
IntegrationPointData< BMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, DisplacementDim, ShapeFunctionDisplacement::NPOINTS > IpData
ThermoHydroMechanicsProcessData< DisplacementDim > & _process_data
std::vector< IntegrationPointDataForOutput< DisplacementDim >, Eigen::aligned_allocator< IntegrationPointDataForOutput< DisplacementDim > > > _ip_data_output
void postTimestepConcrete(Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, double const t, double const dt, int const) override
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
virtual std::vector< double > const & getIntPtEpsilonM(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< double > const & getIntPtIceVolume(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
MathLib::KelvinVector::Invariants< KelvinVectorSize > Invariants
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt(unsigned integration_point) const override
std::size_t setIPDataInitialConditions(std::string_view const name, double const *values, int const integration_order) override
Returns number of read integration points.
virtual std::vector< double > const & getIntPtEpsilon0(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
ThermoHydroMechanicsLocalAssembler(ThermoHydroMechanicsLocalAssembler &&)=delete
void computeSecondaryVariableConcrete(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
std::vector< double > const & getIntPtSigma(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
constexpr Eigen::CwiseNullaryOp< EigenBlockMatrixViewFunctor< D, M >, typename EigenBlockMatrixViewFunctor< D, M >::Matrix > eigenBlockMatrixView(const Eigen::MatrixBase< M > &matrix)
auto localDOF(ElementDOFVector const &x)
Definition LocalDOF.h:56
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, const double radius, const bool is_axially_symmetric)
Fills a B-matrix based on given shape function dN/dx values.
std::vector< double > const & getIntegrationPointScalarData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)
std::vector< double > transposeInPlace(StoreValuesFunction const &store_values_function)
std::vector< double > getIntegrationPointDataMaterialStateVariables(IntegrationPointDataVector const &ip_data_vector, MemberType member, std::function< std::span< double >(MaterialStateVariables &)> get_values_span, int const n_components)
std::vector< double > const & getIntegrationPointKelvinVectorData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)
std::size_t setIntegrationPointKelvinVectorData(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType IpData::*const member)
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u