OGS
ThermoHydroMechanicsFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <memory>
14#include <vector>
15
24#include "NumLib/DOF/LocalDOF.h"
35
36namespace ProcessLib
37{
38namespace ThermoHydroMechanics
39{
40namespace MPL = MaterialPropertyLib;
41
44template <typename ShapeMatrixType>
46{
47 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>> N_u;
48};
49
50template <typename ShapeFunctionDisplacement, typename ShapeFunctionPressure,
51 int DisplacementDim>
53 : public LocalAssemblerInterface<DisplacementDim>
54{
55public:
58
59 // Types for pressure.
62
65
68
69 static int const KelvinVectorSize =
72
73 using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>;
74
75 static constexpr auto& N_u_op = MathLib::eigenBlockMatrixView<
76 DisplacementDim,
78
82 delete;
83
85 MeshLib::Element const& e,
86 std::size_t const /*local_matrix_size*/,
87 NumLib::GenericIntegrationMethod const& integration_method,
88 bool const is_axially_symmetric,
90
93 std::string_view const name,
94 double const* values,
95 int const integration_order) override;
96
97 void assemble(double const /*t*/, double const /*dt*/,
98 std::vector<double> const& /*local_x*/,
99 std::vector<double> const& /*local_x_prev*/,
100 std::vector<double>& /*local_M_data*/,
101 std::vector<double>& /*local_K_data*/,
102 std::vector<double>& /*local_rhs_data*/) override
103 {
104 OGS_FATAL(
105 "ThermoHydroMechanicsLocalAssembler: assembly without Jacobian is "
106 "not implemented.");
107 }
108
109 void assembleWithJacobian(double const t, double const dt,
110 std::vector<double> const& local_x,
111 std::vector<double> const& local_x_prev,
112 std::vector<double>& local_rhs_data,
113 std::vector<double>& local_Jac_data) override;
114
115 void initializeConcrete() override
116 {
117 unsigned const n_integration_points =
119
120 for (unsigned ip = 0; ip < n_integration_points; ip++)
121 {
122 auto& ip_data = _ip_data[ip];
123
124 ParameterLib::SpatialPosition const x_position{
125 std::nullopt, _element.getID(), ip,
128 ShapeFunctionDisplacement,
130
132 if (_process_data.initial_stress.value)
133 {
134 ip_data.sigma_eff =
136 DisplacementDim>((*_process_data.initial_stress.value)(
137 std::numeric_limits<double>::quiet_NaN() /* time
138 independent
139 */
140 ,
141 x_position));
142 }
143
144 double const t = 0; // TODO (naumov) pass t from top
145 ip_data.solid_material.initializeInternalStateVariables(
146 t, x_position, *ip_data.material_state_variables);
147
148 ip_data.pushBackState();
149 }
150 }
151
152 void setInitialConditionsConcrete(Eigen::VectorXd const local_x,
153 double const t,
154 int const process_id) override;
155
156 void preTimestepConcrete(std::vector<double> const& /*local_x*/,
157 double const /*t*/, double const /*dt*/) override
158 {
159 unsigned const n_integration_points =
161
162 for (unsigned ip = 0; ip < n_integration_points; ip++)
163 {
164 _ip_data_output[ip].velocity.setConstant(
165 DisplacementDim, std::numeric_limits<double>::quiet_NaN());
166 }
167 }
168
169 void postTimestepConcrete(Eigen::VectorXd const& local_x,
170 Eigen::VectorXd const& local_x_prev,
171 double const t, double const dt,
172 int const /*process_id*/) override
173 {
174 unsigned const n_integration_points =
176
177 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
178
179 for (unsigned ip = 0; ip < n_integration_points; ip++)
180 {
181 auto& ip_data = _ip_data[ip];
182 auto const& N_u = ip_data.N_u;
183 auto const& dNdx_u = ip_data.dNdx_u;
184
185 ParameterLib::SpatialPosition const x_position{
186 std::nullopt, _element.getID(), ip,
189 ShapeFunctionDisplacement,
191
192 updateConstitutiveRelations(local_x, local_x_prev, x_position, t,
193 dt, _ip_data[ip], _ip_data_output[ip]);
194
195 auto const x_coord =
196 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
198 _element, N_u);
199 auto const B = LinearBMatrix::computeBMatrix<
200 DisplacementDim, ShapeFunctionDisplacement::NPOINTS,
201 typename BMatricesType::BMatrixType>(dNdx_u, N_u, x_coord,
203
205
207 eps_prev = B * u_prev;
208
209 _ip_data[ip].eps0 =
210 _ip_data[ip].eps0_prev +
211 (1 - _ip_data[ip].phi_fr_prev / _ip_data[ip].porosity) *
212 (eps_prev - _ip_data[ip].eps0_prev);
213 _ip_data[ip].pushBackState();
214 }
215 }
216
218 double const t, double const dt, Eigen::VectorXd const& local_x,
219 Eigen::VectorXd const& local_x_prev) override;
220
221 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
222 const unsigned integration_point) const override
223 {
224 auto const& N_u = _secondary_data.N_u[integration_point];
225
226 // assumes N is stored contiguously in memory
227 return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
228 }
229
230 std::vector<double> const& getIntPtDarcyVelocity(
231 const double t,
232 std::vector<GlobalVector*> const& x,
233 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
234 std::vector<double>& cache) const override;
235
236 // TODO (naumov) This method is same as getIntPtSigma but for arguments and
237 // the ordering of the cache_mat.
238 // There should be only one.
239 std::vector<double> getSigma() const override
240 {
241 constexpr int kelvin_vector_size =
243
245 [this](std::vector<double>& values)
246 { return getIntPtSigma(0, {}, {}, values); });
247 }
248
249 std::vector<double> const& getIntPtFluidDensity(
250 const double t,
251 std::vector<GlobalVector*> const& x,
252 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
253 std::vector<double>& cache) const override;
254
255 std::vector<double> const& getIntPtViscosity(
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
261private:
264 using IpData =
266 ShapeMatricesTypePressure, DisplacementDim,
267 ShapeFunctionDisplacement::NPOINTS>;
268
269private:
271 Eigen::Ref<Eigen::VectorXd const> const local_x,
272 Eigen::Ref<Eigen::VectorXd const> const local_x_prev,
273 ParameterLib::SpatialPosition const& x_position, double const t,
274 double const dt, IpData& ip_data,
276
277 std::size_t setSigma(double const* values)
278 {
280 values, _ip_data, &IpData::sigma_eff);
281 }
282
283 std::vector<double> const& getIntPtSigma(
284 const double /*t*/,
285 std::vector<GlobalVector*> const& /*x*/,
286 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
287 std::vector<double>& cache) const override
288 {
290 _ip_data, &IpData::sigma_eff, cache);
291 }
292
293 std::vector<double> const& getIntPtSigmaIce(
294 const double /*t*/,
295 std::vector<GlobalVector*> const& /*x*/,
296 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
297 std::vector<double>& cache) const override
298 {
301 }
302
303 std::vector<double> getEpsilonM() const override
304 {
305 constexpr int kelvin_vector_size =
307
309 [this](std::vector<double>& values)
310 { return getIntPtEpsilonM(0, {}, {}, values); });
311 }
312
313 virtual std::vector<double> const& getIntPtEpsilonM(
314 const double /*t*/,
315 std::vector<GlobalVector*> const& /*x*/,
316 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
317 std::vector<double>& cache) const override
318 {
320 _ip_data, &IpData::eps_m, cache);
321 }
322
323 std::vector<double> getEpsilon() const override
324 {
325 constexpr int kelvin_vector_size =
327
329 [this](std::vector<double>& values)
330 { return getIntPtEpsilon(0, {}, {}, values); });
331 }
332
333 virtual std::vector<double> const& getIntPtEpsilon(
334 const double /*t*/,
335 std::vector<GlobalVector*> const& /*x*/,
336 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
337 std::vector<double>& cache) const override
338 {
340 _ip_data, &IpData::eps, cache);
341 }
342
343 std::vector<double> const& getIntPtIceVolume(
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::phi_fr, cache);
351 }
352
353 unsigned getNumberOfIntegrationPoints() const override
354 {
356 }
357
358 int getMaterialID() const override
359 {
360 return _process_data.material_ids == nullptr
361 ? 0
362 : (*_process_data.material_ids)[_element.getID()];
363 }
364
366 std::function<std::span<double>(
368 MaterialStateVariables&)> const& get_values_span,
369 int const& n_components) const override
370 {
373 n_components);
374 }
375
377 DisplacementDim>::MaterialStateVariables const&
378 getMaterialStateVariablesAt(unsigned integration_point) const override
379 {
380 return *_ip_data[integration_point].material_state_variables;
381 }
382
383private:
384 template <typename SolutionVector>
385 static constexpr auto localDOF(SolutionVector const& x)
386 {
387 return NumLib::localDOF<
388 ShapeFunctionPressure, ShapeFunctionPressure,
390 }
391
393
394 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
395 std::vector<IntegrationPointDataForOutput<DisplacementDim>,
396 Eigen::aligned_allocator<
399
404 typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType>
406
407 // The shape function of pressure has the same form with the shape function
408 // of temperature
409 static const int temperature_index = 0;
410 static const int temperature_size = ShapeFunctionPressure::NPOINTS;
411 static const int pressure_index = ShapeFunctionPressure::NPOINTS;
412 static const int pressure_size = ShapeFunctionPressure::NPOINTS;
413 static const int displacement_index = ShapeFunctionPressure::NPOINTS * 2;
414 static const int displacement_size =
415 ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
416};
417
418} // namespace ThermoHydroMechanics
419} // namespace ProcessLib
420
#define OGS_FATAL(...)
Definition Error.h:26
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
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
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
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.
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)
double interpolateXCoordinate(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
auto localDOF(ElementDOFVector const &x)
Definition LocalDOF.h:64
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::unique_ptr< typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables > material_state_variables
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u