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