OGS
StaggeredHTFEM-impl.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 <typeinfo>
7
13#include "StaggeredHTFEM.h"
14
15namespace ProcessLib
16{
17namespace HT
18{
19template <typename ShapeFunction, int GlobalDim>
21 double const t, double const dt, Eigen::VectorXd const& local_x,
22 Eigen::VectorXd const& local_x_prev, int const process_id,
23 std::vector<double>& local_M_data, std::vector<double>& local_K_data,
24 std::vector<double>& local_b_data)
25{
26 if (process_id == this->_process_data.heat_transport_process_id)
27 {
28 assembleHeatTransportEquation(t, dt, local_x, local_M_data,
29 local_K_data);
30 return;
31 }
32
33 assembleHydraulicEquation(t, dt, local_x, local_x_prev, local_M_data,
34 local_K_data, local_b_data);
35}
36
37template <typename ShapeFunction, int GlobalDim>
39 double const t, double const dt, Eigen::VectorXd const& local_x,
40 Eigen::VectorXd const& local_x_prev, std::vector<double>& local_M_data,
41 std::vector<double>& local_K_data, std::vector<double>& local_b_data)
42{
43 auto const local_p =
44 local_x.template segment<pressure_size>(pressure_index);
45 auto const local_T =
46 local_x.template segment<temperature_size>(temperature_index);
47
48 auto const local_T_prev =
49 local_x_prev.template segment<temperature_size>(temperature_index);
50
52 local_M_data, pressure_size, pressure_size);
54 local_K_data, pressure_size, pressure_size);
55 auto local_b = MathLib::createZeroedVector<LocalVectorType>(local_b_data,
57
58 auto const& process_data = this->_process_data;
59 auto const& medium =
60 *this->_process_data.media_map.getMedium(this->_element.getID());
61 auto const& liquid_phase =
63 auto const& solid_phase =
65
66 auto const& b =
67 process_data
68 .projected_specific_body_force_vectors[this->_element.getID()];
69
71
72 unsigned const n_integration_points =
73 this->_integration_method.getNumberOfPoints();
74
75 auto const& Ns =
76 process_data.shape_matrix_cache
77 .template NsHigherOrder<typename ShapeFunction::MeshElement>();
78
79 for (unsigned ip(0); ip < n_integration_points; ip++)
80 {
81 auto const& ip_data = this->_ip_data[ip];
82 auto const& dNdx = ip_data.dNdx;
83 auto const& N = Ns[ip];
84 auto const& w = ip_data.integration_weight;
85
87 std::nullopt, this->_element.getID(),
90 this->_element, N))};
91
92 double p_int_pt = 0.0;
93 double T_int_pt = 0.0;
94 NumLib::shapeFunctionInterpolate(local_p, N, p_int_pt);
95 NumLib::shapeFunctionInterpolate(local_T, N, T_int_pt);
96
97 vars.temperature = T_int_pt;
98 vars.liquid_phase_pressure = p_int_pt;
99
100 vars.liquid_saturation = 1.0;
101
102 auto const porosity =
104 .template value<double>(vars, pos, t, dt);
105 auto const fluid_density =
106 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
107 .template value<double>(vars, pos, t, dt);
108
109 vars.density = fluid_density;
110 const double dfluid_density_dp =
111 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
112 .template dValue<double>(
114 pos, t, dt);
115
116 // Use the viscosity model to compute the viscosity
117 auto const viscosity =
119 .template value<double>(vars, pos, t, dt);
120
121 // \todo the argument to getValue() has to be changed for non
122 // constant storage model
123 auto const specific_storage =
125 .template value<double>(vars, pos, t, dt);
126
127 auto const intrinsic_permeability =
130 .value(vars, pos, t, dt));
131 GlobalDimMatrixType const K_over_mu =
132 intrinsic_permeability / viscosity;
133
134 // matrix assembly
135 local_M.noalias() +=
136 w *
137 (porosity * dfluid_density_dp / fluid_density + specific_storage) *
138 N.transpose() * N;
139
140 local_K.noalias() += w * dNdx.transpose() * K_over_mu * dNdx;
141
142 if (process_data.has_gravity)
143 {
144 local_b.noalias() +=
145 w * fluid_density * dNdx.transpose() * K_over_mu * b;
146 }
147
148 if (!process_data.has_fluid_thermal_expansion)
149 {
150 return;
151 }
152
153 // Add the thermal expansion term
154 {
155 auto const solid_thermal_expansion =
156 process_data.solid_thermal_expansion(t, pos)[0];
157 const double dfluid_density_dT =
158 liquid_phase
160 .template dValue<double>(
162 t, dt);
163 double const Tdot_int_pt = (T_int_pt - local_T_prev.dot(N)) / dt;
164 auto const biot_constant = process_data.biot_constant(t, pos)[0];
165 const double eff_thermal_expansion =
166 3.0 * (biot_constant - porosity) * solid_thermal_expansion -
167 porosity * dfluid_density_dT / fluid_density;
168 local_b.noalias() += eff_thermal_expansion * Tdot_int_pt * w * N;
169 }
170 }
171}
172
173template <typename ShapeFunction, int GlobalDim>
175 double const t,
176 double const dt,
177 Eigen::VectorXd const& local_x,
178 std::vector<double>& local_M_data,
179 std::vector<double>& local_K_data)
180{
181 auto const local_p =
182 local_x.template segment<pressure_size>(pressure_index);
183 auto const local_T =
184 local_x.template segment<temperature_size>(temperature_index);
185
187 local_M_data, temperature_size, temperature_size);
189 local_K_data, temperature_size, temperature_size);
190
191 auto const& process_data = this->_process_data;
192 auto const& medium =
193 *process_data.media_map.getMedium(this->_element.getID());
194 auto const& liquid_phase =
196
197 auto const& b =
198 process_data
199 .projected_specific_body_force_vectors[this->_element.getID()];
200
202
203 unsigned const n_integration_points =
204 this->_integration_method.getNumberOfPoints();
205
206 std::vector<GlobalDimVectorType> ip_flux_vector;
207 double average_velocity_norm = 0.0;
208 ip_flux_vector.reserve(n_integration_points);
209
210 auto const& Ns =
211 process_data.shape_matrix_cache
212 .template NsHigherOrder<typename ShapeFunction::MeshElement>();
213
214 for (unsigned ip(0); ip < n_integration_points; ip++)
215 {
216 auto const& ip_data = this->_ip_data[ip];
217 auto const& dNdx = ip_data.dNdx;
218 auto const& N = Ns[ip];
219 auto const& w = ip_data.integration_weight;
220
222 std::nullopt, this->_element.getID(),
225 this->_element, N))};
226
227 double p_at_xi = 0.;
228 NumLib::shapeFunctionInterpolate(local_p, N, p_at_xi);
229 double T_at_xi = 0.;
230 NumLib::shapeFunctionInterpolate(local_T, N, T_at_xi);
231
232 vars.temperature = T_at_xi;
233 vars.liquid_phase_pressure = p_at_xi;
234
235 vars.liquid_saturation = 1.0;
236
237 auto const porosity =
239 .template value<double>(vars, pos, t, dt);
240 vars.porosity = porosity;
241
242 // Use the fluid density model to compute the density
243 auto const fluid_density =
244 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
245 .template value<double>(vars, pos, t, dt);
246 vars.density = fluid_density;
247 auto const specific_heat_capacity_fluid =
249 .template value<double>(vars, pos, t, dt);
250
251 // Assemble mass matrix
252 local_M.noalias() += w *
254 vars, porosity, fluid_density,
255 specific_heat_capacity_fluid, pos, t, dt) *
256 N.transpose() * N;
257
258 // Assemble Laplace matrix
259 auto const viscosity =
261 .template value<double>(vars, pos, t, dt);
262
263 auto const intrinsic_permeability =
266 .value(vars, pos, t, dt));
267
268 GlobalDimMatrixType const K_over_mu =
269 intrinsic_permeability / viscosity;
270 GlobalDimVectorType const velocity =
271 process_data.has_gravity
272 ? GlobalDimVectorType(-K_over_mu *
273 (dNdx * local_p - fluid_density * b))
274 : GlobalDimVectorType(-K_over_mu * dNdx * local_p);
275
276 GlobalDimMatrixType const thermal_conductivity_dispersivity =
278 vars, fluid_density, specific_heat_capacity_fluid, velocity,
279 pos, t, dt);
280
281 local_K.noalias() +=
282 w * dNdx.transpose() * thermal_conductivity_dispersivity * dNdx;
283
284 ip_flux_vector.emplace_back(velocity * fluid_density *
285 specific_heat_capacity_fluid);
286 average_velocity_norm += velocity.norm();
287 }
288
290 process_data.stabilizer, this->_ip_data,
291 process_data.shape_matrix_cache, ip_flux_vector,
292 average_velocity_norm / static_cast<double>(n_integration_points),
293 local_K);
294}
295
296template <typename ShapeFunction, int GlobalDim>
297std::vector<double> const&
299 const double t,
300 std::vector<GlobalVector*> const& x,
301 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
302 std::vector<double>& cache) const
303{
304 assert(x.size() == dof_table.size());
305 auto const n_processes = dof_table.size();
306
307 std::vector<std::vector<GlobalIndexType>> indices_of_all_coupled_processes;
308 indices_of_all_coupled_processes.reserve(n_processes);
309 for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
310 {
311 auto const indices =
312 NumLib::getIndices(this->_element.getID(), *dof_table[process_id]);
313 assert(!indices.empty());
314 indices_of_all_coupled_processes.push_back(indices);
315 }
316 auto const local_xs =
317 getCoupledLocalSolutions(x, indices_of_all_coupled_processes);
318
319 return this->getIntPtDarcyVelocityLocal(t, local_xs, cache);
320}
321} // namespace HT
322} // namespace ProcessLib
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
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
void assembleForStaggeredScheme(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
void assembleHydraulicEquation(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
void assembleHeatTransportEquation(double const t, double const dt, Eigen::VectorXd const &local_x, std::vector< double > &local_M_data, std::vector< double > &local_K_data)
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
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
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
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)
std::vector< double > getCoupledLocalSolutions(std::vector< GlobalVector * > const &global_solutions, std::vector< std::vector< GlobalIndexType > > const &indices)