OGS
MonolithicHTFEM.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 <Eigen/Core>
7#include <typeinfo>
8#include <vector>
9
10#include "HTFEM.h"
11#include "HTProcessData.h"
22
23namespace ProcessLib
24{
25namespace HT
26{
27const unsigned NUM_NODAL_DOF = 2;
28
29template <typename ShapeFunction, int GlobalDim>
30class MonolithicHTFEM : public HTFEM<ShapeFunction, GlobalDim>
31{
34
35 using LocalMatrixType = typename ShapeMatricesType::template MatrixType<
36 NUM_NODAL_DOF * ShapeFunction::NPOINTS,
37 NUM_NODAL_DOF * ShapeFunction::NPOINTS>;
39 typename ShapeMatricesType::template VectorType<NUM_NODAL_DOF *
40 ShapeFunction::NPOINTS>;
41
45
48
49public:
51 std::size_t const local_matrix_size,
52 NumLib::GenericIntegrationMethod const& integration_method,
53 bool is_axially_symmetric,
54 HTProcessData const& process_data)
55 : HTFEM<ShapeFunction, GlobalDim>(
56 element, local_matrix_size, integration_method,
57 is_axially_symmetric, process_data, NUM_NODAL_DOF)
58 {
59 }
60
61 void assemble(double const t, double const dt,
62 std::vector<double> const& local_x,
63 std::vector<double> const& local_x_prev,
64 std::vector<double>& local_M_data,
65 std::vector<double>& local_K_data,
66 std::vector<double>& local_b_data) override
67 {
68 auto const local_matrix_size = local_x.size();
69 // This assertion is valid only if all nodal d.o.f. use the same shape
70 // matrices.
71 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
72
74 local_M_data, local_matrix_size, local_matrix_size);
76 local_K_data, local_matrix_size, local_matrix_size);
78 local_b_data, local_matrix_size);
79
80 auto KTT = local_K.template block<temperature_size, temperature_size>(
82 auto MTT = local_M.template block<temperature_size, temperature_size>(
84 auto Kpp = local_K.template block<pressure_size, pressure_size>(
86 auto Mpp = local_M.template block<pressure_size, pressure_size>(
88 auto Bp = local_b.template block<pressure_size, 1>(pressure_index, 0);
89
90 auto const& process_data = this->_process_data;
91
92 auto p_nodal_values = Eigen::Map<const NodalVectorType>(
93 &local_x[pressure_index], pressure_size);
94
95 auto const& medium =
96 *process_data.media_map.getMedium(this->_element.getID());
97 auto const& liquid_phase =
99 auto const& solid_phase =
101
102 auto const& b =
103 process_data
104 .projected_specific_body_force_vectors[this->_element.getID()];
105
107
108 unsigned const n_integration_points =
109 this->_integration_method.getNumberOfPoints();
110
111 std::vector<GlobalDimVectorType> ip_flux_vector;
112 double average_velocity_norm = 0.0;
113 ip_flux_vector.reserve(n_integration_points);
114
115 auto const& Ns =
116 process_data.shape_matrix_cache
117 .template NsHigherOrder<typename ShapeFunction::MeshElement>();
118
119 for (unsigned ip(0); ip < n_integration_points; ip++)
120 {
121 auto const& ip_data = this->_ip_data[ip];
122 auto const& dNdx = ip_data.dNdx;
123 auto const& N = Ns[ip];
124 auto const& w = ip_data.integration_weight;
125
127 std::nullopt, this->_element.getID(),
131 this->_element, N))};
132
133 double T_int_pt = 0.0;
134 double p_int_pt = 0.0;
135 // Order matters: First T, then P!
136 NumLib::shapeFunctionInterpolate(local_x, N, T_int_pt, p_int_pt);
137
138 vars.temperature = T_int_pt;
139 vars.liquid_phase_pressure = p_int_pt;
140
141 vars.liquid_saturation = 1.0;
142 // \todo the argument to getValue() has to be changed for non
143 // constant storage model
144 auto const specific_storage =
146 .template value<double>(vars, pos, t, dt);
147
148 auto const porosity =
150 .template value<double>(vars, pos, t, dt);
151 vars.porosity = porosity;
152
153 auto const intrinsic_permeability =
155 medium
156 .property(
158 .value(vars, pos, t, dt));
159
160 auto const specific_heat_capacity_fluid =
161 liquid_phase
163 .template value<double>(vars, pos, t, dt);
164
165 // Use the fluid density model to compute the density
166 auto const fluid_density =
167 liquid_phase
169 .template value<double>(vars, pos, t, dt);
170
171 vars.density = fluid_density;
172 const double dfluid_density_dp =
173 liquid_phase
175 .template dValue<double>(
176 vars,
178 pos, t, dt);
179
180 // Use the viscosity model to compute the viscosity
181 auto const viscosity =
182 liquid_phase
184 .template value<double>(vars, pos, t, dt);
185 GlobalDimMatrixType K_over_mu = intrinsic_permeability / viscosity;
186
187 GlobalDimVectorType const velocity =
188 process_data.has_gravity
189 ? GlobalDimVectorType(-K_over_mu * (dNdx * p_nodal_values -
190 fluid_density * b))
191 : GlobalDimVectorType(-K_over_mu * dNdx * p_nodal_values);
192
193 // matrix assembly
194 GlobalDimMatrixType const thermal_conductivity_dispersivity =
196 vars, fluid_density, specific_heat_capacity_fluid, velocity,
197 pos, t, dt);
198
199 KTT.noalias() +=
200 dNdx.transpose() * thermal_conductivity_dispersivity * dNdx * w;
201
202 ip_flux_vector.emplace_back(velocity * fluid_density *
203 specific_heat_capacity_fluid);
204 average_velocity_norm += velocity.norm();
205
206 Kpp.noalias() += w * dNdx.transpose() * K_over_mu * dNdx;
207 MTT.noalias() += w *
209 vars, porosity, fluid_density,
210 specific_heat_capacity_fluid, pos, t, dt) *
211 N.transpose() * N;
212 Mpp.noalias() += w *
213 (porosity * dfluid_density_dp / fluid_density +
214 specific_storage) *
215 N.transpose() * N;
216 if (process_data.has_gravity)
217 {
218 Bp += w * fluid_density * dNdx.transpose() * K_over_mu * b;
219 }
220
221 if (process_data.has_fluid_thermal_expansion)
222 {
223 double const solid_thermal_expansion =
224 process_data.solid_thermal_expansion(t, pos)[0];
225 double const dfluid_density_dT =
226 liquid_phase
228 .template dValue<double>(
230 pos, t, dt);
231 double const Tdot_int_pt =
232 (T_int_pt -
233 Eigen::Map<const NodalVectorType>(
234 &local_x_prev[temperature_index], temperature_size)
235 .dot(N)) /
236 dt;
237 double const biot_constant =
238 process_data.biot_constant(t, pos)[0];
239 double const eff_thermal_expansion =
240 3.0 * (biot_constant - porosity) * solid_thermal_expansion -
241 porosity * dfluid_density_dT / fluid_density;
242 Bp.noalias() += eff_thermal_expansion * Tdot_int_pt * w * N;
243 }
244 }
245
247 process_data.stabilizer, this->_ip_data,
248 process_data.shape_matrix_cache, ip_flux_vector,
249 average_velocity_norm / static_cast<double>(n_integration_points),
250 KTT);
251 }
252
253 std::vector<double> const& getIntPtDarcyVelocity(
254 const double t,
255 std::vector<GlobalVector*> const& x,
256 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
257 std::vector<double>& cache) const override
258 {
259 int const process_id = 0; // monolithic case.
260 auto const indices =
261 NumLib::getIndices(this->_element.getID(), *dof_table[process_id]);
262 assert(!indices.empty());
263 auto const& local_x = x[process_id]->get(indices);
264
265 return this->getIntPtDarcyVelocityLocal(t, local_x, cache);
266 }
267
268private:
269 using HTFEM<ShapeFunction, GlobalDim>::pressure_index;
270 using HTFEM<ShapeFunction, GlobalDim>::pressure_size;
271 using HTFEM<ShapeFunction, GlobalDim>::temperature_index;
272 using HTFEM<ShapeFunction, GlobalDim>::temperature_size;
273};
274
275} // namespace HT
276} // namespace ProcessLib
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
std::vector< double > const & getIntPtDarcyVelocityLocal(const double t, std::vector< double > const &local_x, std::vector< double > &cache) const
Definition HTFEM.h:234
NumLib::GenericIntegrationMethod const & _integration_method
Definition HTFEM.h:165
static const int temperature_index
Definition HTFEM.h:319
static const int temperature_size
Definition HTFEM.h:320
double getHeatEnergyCoefficient(MaterialPropertyLib::VariableArray const &vars, const double porosity, const double fluid_density, const double specific_heat_capacity_fluid, ParameterLib::SpatialPosition const &pos, double const t, double const dt)
Definition HTFEM.h:168
HTFEM(MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, HTProcessData const &process_data, const unsigned dof_per_node)
Definition HTFEM.h:41
static const int pressure_size
Definition HTFEM.h:318
GlobalDimMatrixType getThermalConductivityDispersivity(MaterialPropertyLib::VariableArray const &vars, const double fluid_density, const double specific_heat_capacity_fluid, const GlobalDimVectorType &velocity, ParameterLib::SpatialPosition const &pos, double const t, double const dt)
Definition HTFEM.h:193
MeshLib::Element const & _element
Definition HTFEM.h:162
HTProcessData const & _process_data
Definition HTFEM.h:163
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
Definition HTFEM.h:166
static const int pressure_index
Definition HTFEM.h:317
typename ShapeMatricesType::NodalMatrixType NodalMatrixType
typename ShapeMatricesType::template MatrixType< NUM_NODAL_DOF *ShapeFunction::NPOINTS, NUM_NODAL_DOF *ShapeFunction::NPOINTS > LocalMatrixType
typename ShapeMatricesType::template VectorType< NUM_NODAL_DOF * ShapeFunction::NPOINTS > LocalVectorType
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
typename ShapeMatricesType::NodalRowVectorType NodalRowVectorType
MonolithicHTFEM(MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool is_axially_symmetric, HTProcessData const &process_data)
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
typename ShapeMatricesType::NodalVectorType NodalVectorType
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
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
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
constexpr Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
void assembleAdvectionMatrix(IPData const &ip_data_vector, NumLib::ShapeMatrixCache const &shape_matrix_cache, std::vector< FluxVectorType > const &ip_flux_vector, Eigen::MatrixBase< Derived > &laplacian_matrix)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
const unsigned NUM_NODAL_DOF
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType