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