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_M_data*/,
113 std::vector<double>& /*local_K_data*/,
114 std::vector<double>& local_rhs_data,
115 std::vector<double>& local_Jac_data) override;
116
117 void initializeConcrete() override
118 {
119 unsigned const n_integration_points =
121
122 for (unsigned ip = 0; ip < n_integration_points; ip++)
123 {
124 auto& ip_data = _ip_data[ip];
125
126 ParameterLib::SpatialPosition const x_position{
127 std::nullopt, _element.getID(), ip,
130 ShapeFunctionDisplacement,
132
134 if (_process_data.initial_stress.value)
135 {
136 ip_data.sigma_eff =
138 DisplacementDim>((*_process_data.initial_stress.value)(
139 std::numeric_limits<double>::quiet_NaN() /* time
140 independent
141 */
142 ,
143 x_position));
144 }
145
146 double const t = 0; // TODO (naumov) pass t from top
147 ip_data.solid_material.initializeInternalStateVariables(
148 t, x_position, *ip_data.material_state_variables);
149
150 ip_data.pushBackState();
151 }
152 }
153
154 void setInitialConditionsConcrete(Eigen::VectorXd const local_x,
155 double const t,
156 int const process_id) override;
157
158 void preTimestepConcrete(std::vector<double> const& /*local_x*/,
159 double const /*t*/, double const /*dt*/) override
160 {
161 unsigned const n_integration_points =
163
164 for (unsigned ip = 0; ip < n_integration_points; ip++)
165 {
166 _ip_data_output[ip].velocity.setConstant(
167 DisplacementDim, std::numeric_limits<double>::quiet_NaN());
168 }
169 }
170
171 void postTimestepConcrete(Eigen::VectorXd const& local_x,
172 Eigen::VectorXd const& local_x_prev,
173 double const t, double const dt,
174 int const /*process_id*/) override
175 {
176 unsigned const n_integration_points =
178
179 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
180
181 for (unsigned ip = 0; ip < n_integration_points; ip++)
182 {
183 auto& ip_data = _ip_data[ip];
184 auto const& N_u = ip_data.N_u;
185 auto const& dNdx_u = ip_data.dNdx_u;
186
187 ParameterLib::SpatialPosition const x_position{
188 std::nullopt, _element.getID(), ip,
191 ShapeFunctionDisplacement,
193
194 updateConstitutiveRelations(local_x, local_x_prev, x_position, t,
195 dt, _ip_data[ip], _ip_data_output[ip]);
196
197 auto const x_coord =
198 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
200 _element, N_u);
201 auto const B = LinearBMatrix::computeBMatrix<
202 DisplacementDim, ShapeFunctionDisplacement::NPOINTS,
203 typename BMatricesType::BMatrixType>(dNdx_u, N_u, x_coord,
205
207
209 eps_prev = B * u_prev;
210
211 _ip_data[ip].eps0 =
212 _ip_data[ip].eps0_prev +
213 (1 - _ip_data[ip].phi_fr_prev / _ip_data[ip].porosity) *
214 (eps_prev - _ip_data[ip].eps0_prev);
215 _ip_data[ip].pushBackState();
216 }
217 }
218
220 double const t, double const dt, Eigen::VectorXd const& local_x,
221 Eigen::VectorXd const& local_x_prev) override;
222
223 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
224 const unsigned integration_point) const override
225 {
226 auto const& N_u = _secondary_data.N_u[integration_point];
227
228 // assumes N is stored contiguously in memory
229 return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
230 }
231
232 std::vector<double> const& getIntPtDarcyVelocity(
233 const double t,
234 std::vector<GlobalVector*> const& x,
235 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
236 std::vector<double>& cache) const override;
237
238 // TODO (naumov) This method is same as getIntPtSigma but for arguments and
239 // the ordering of the cache_mat.
240 // There should be only one.
241 std::vector<double> getSigma() const override
242 {
243 constexpr int kelvin_vector_size =
245
246 return transposeInPlace<kelvin_vector_size>(
247 [this](std::vector<double>& values)
248 { return getIntPtSigma(0, {}, {}, values); });
249 }
250
251 std::vector<double> const& getIntPtFluidDensity(
252 const double t,
253 std::vector<GlobalVector*> const& x,
254 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
255 std::vector<double>& cache) const override;
256
257 std::vector<double> const& getIntPtViscosity(
258 const double t,
259 std::vector<GlobalVector*> const& x,
260 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
261 std::vector<double>& cache) const override;
262
263private:
266 using IpData =
268 ShapeMatricesTypePressure, DisplacementDim,
269 ShapeFunctionDisplacement::NPOINTS>;
270
271private:
273 Eigen::Ref<Eigen::VectorXd const> const local_x,
274 Eigen::Ref<Eigen::VectorXd const> const local_x_prev,
275 ParameterLib::SpatialPosition const& x_position, double const t,
276 double const dt, IpData& ip_data,
278
279 std::size_t setSigma(double const* values)
280 {
281 return ProcessLib::setIntegrationPointKelvinVectorData<DisplacementDim>(
282 values, _ip_data, &IpData::sigma_eff);
283 }
284
285 std::vector<double> const& getIntPtSigma(
286 const double /*t*/,
287 std::vector<GlobalVector*> const& /*x*/,
288 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
289 std::vector<double>& cache) const override
290 {
291 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
292 _ip_data, &IpData::sigma_eff, cache);
293 }
294
295 std::vector<double> const& getIntPtSigmaIce(
296 const double /*t*/,
297 std::vector<GlobalVector*> const& /*x*/,
298 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
299 std::vector<double>& cache) const override
300 {
301 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
303 }
304
305 std::vector<double> getEpsilonM() const override
306 {
307 constexpr int kelvin_vector_size =
309
310 return transposeInPlace<kelvin_vector_size>(
311 [this](std::vector<double>& values)
312 { return getIntPtEpsilonM(0, {}, {}, values); });
313 }
314
315 virtual std::vector<double> const& getIntPtEpsilonM(
316 const double /*t*/,
317 std::vector<GlobalVector*> const& /*x*/,
318 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
319 std::vector<double>& cache) const override
320 {
321 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
322 _ip_data, &IpData::eps_m, cache);
323 }
324
325 std::vector<double> getEpsilon() const override
326 {
327 constexpr int kelvin_vector_size =
329
330 return transposeInPlace<kelvin_vector_size>(
331 [this](std::vector<double>& values)
332 { return getIntPtEpsilon(0, {}, {}, values); });
333 }
334
335 virtual std::vector<double> const& getIntPtEpsilon(
336 const double /*t*/,
337 std::vector<GlobalVector*> const& /*x*/,
338 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
339 std::vector<double>& cache) const override
340 {
341 return ProcessLib::getIntegrationPointKelvinVectorData<DisplacementDim>(
342 _ip_data, &IpData::eps, cache);
343 }
344
345 std::vector<double> const& getIntPtIceVolume(
346 const double /*t*/,
347 std::vector<GlobalVector*> const& /*x*/,
348 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
349 std::vector<double>& cache) const override
350 {
352 _ip_data, &IpData::phi_fr, cache);
353 }
354
355 unsigned getNumberOfIntegrationPoints() const override
356 {
358 }
359
360 int getMaterialID() const override
361 {
362 return _process_data.material_ids == nullptr
363 ? 0
364 : (*_process_data.material_ids)[_element.getID()];
365 }
366
368 std::function<std::span<double>(
370 MaterialStateVariables&)> const& get_values_span,
371 int const& n_components) const override
372 {
375 n_components);
376 }
377
379 DisplacementDim>::MaterialStateVariables const&
380 getMaterialStateVariablesAt(unsigned integration_point) const override
381 {
382 return *_ip_data[integration_point].material_state_variables;
383 }
384
385private:
386 template <typename SolutionVector>
387 static constexpr auto localDOF(SolutionVector const& x)
388 {
389 return NumLib::localDOF<
390 ShapeFunctionPressure, ShapeFunctionPressure,
392 }
393
395
396 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
397 std::vector<IntegrationPointDataForOutput<DisplacementDim>,
398 Eigen::aligned_allocator<
401
406 typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType>
408
409 // The shape function of pressure has the same form with the shape function
410 // of temperature
411 static const int temperature_index = 0;
412 static const int temperature_size = ShapeFunctionPressure::NPOINTS;
413 static const int pressure_index = ShapeFunctionPressure::NPOINTS;
414 static const int pressure_size = ShapeFunctionPressure::NPOINTS;
415 static const int displacement_index = ShapeFunctionPressure::NPOINTS * 2;
416 static const int displacement_size =
417 ShapeFunctionDisplacement::NPOINTS * DisplacementDim;
418};
419
420} // namespace ThermoHydroMechanics
421} // namespace ProcessLib
422
#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
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.
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &, std::vector< double > &, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) 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)
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 > getIntegrationPointDataMaterialStateVariables(IntegrationPointDataVector const &ip_data_vector, MemberType member, std::function< std::span< double >(MaterialStateVariables &)> get_values_span, int const n_components)
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u
std::unique_ptr< typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables > material_state_variables
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N_u