Loading [MathJax]/extensions/tex2jax.js
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
58 local_M_data, pressure_size, pressure_size);
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 auto const& ip_data = this->_ip_data[ip];
89 auto const& dNdx = ip_data.dNdx;
90 auto const& N = Ns[ip];
91 auto const& w = ip_data.integration_weight;
92
93 double p_int_pt = 0.0;
94 double T_int_pt = 0.0;
95 NumLib::shapeFunctionInterpolate(local_p, N, p_int_pt);
96 NumLib::shapeFunctionInterpolate(local_T, N, T_int_pt);
97
98 vars.temperature = T_int_pt;
99 vars.liquid_phase_pressure = p_int_pt;
100
101 vars.liquid_saturation = 1.0;
102
103 auto const porosity =
105 .template value<double>(vars, pos, t, dt);
106 auto const fluid_density =
107 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
108 .template value<double>(vars, pos, t, dt);
109
110 vars.density = fluid_density;
111 const double dfluid_density_dp =
112 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
113 .template dValue<double>(
115 pos, t, dt);
116
117 // Use the viscosity model to compute the viscosity
118 auto const viscosity =
120 .template value<double>(vars, pos, t, dt);
121
122 // \todo the argument to getValue() has to be changed for non
123 // constant storage model
124 auto const specific_storage =
126 .template value<double>(vars, pos, t, dt);
127
128 auto const intrinsic_permeability =
131 .value(vars, pos, t, dt));
132 GlobalDimMatrixType const K_over_mu =
133 intrinsic_permeability / viscosity;
134
135 // matrix assembly
136 local_M.noalias() +=
137 w *
138 (porosity * dfluid_density_dp / fluid_density + specific_storage) *
139 N.transpose() * N;
140
141 local_K.noalias() += w * dNdx.transpose() * K_over_mu * dNdx;
142
143 if (process_data.has_gravity)
144 {
145 local_b.noalias() +=
146 w * fluid_density * dNdx.transpose() * K_over_mu * b;
147 }
148
149 if (!process_data.has_fluid_thermal_expansion)
150 {
151 return;
152 }
153
154 // Add the thermal expansion term
155 {
156 auto const solid_thermal_expansion =
157 process_data.solid_thermal_expansion(t, pos)[0];
158 const double dfluid_density_dT =
159 liquid_phase
161 .template dValue<double>(
163 t, dt);
164 double const Tdot_int_pt = (T_int_pt - local_T_prev.dot(N)) / dt;
165 auto const biot_constant = process_data.biot_constant(t, pos)[0];
166 const double eff_thermal_expansion =
167 3.0 * (biot_constant - porosity) * solid_thermal_expansion -
168 porosity * dfluid_density_dT / fluid_density;
169 local_b.noalias() += eff_thermal_expansion * Tdot_int_pt * w * N;
170 }
171 }
172}
173
174template <typename ShapeFunction, int GlobalDim>
176 double const t,
177 double const dt,
178 Eigen::VectorXd const& local_x,
179 std::vector<double>& local_M_data,
180 std::vector<double>& local_K_data)
181{
182 auto const local_p =
183 local_x.template segment<pressure_size>(pressure_index);
184 auto const local_T =
185 local_x.template segment<temperature_size>(temperature_index);
186
188 local_M_data, temperature_size, temperature_size);
190 local_K_data, temperature_size, temperature_size);
191
193 pos.setElementID(this->_element.getID());
194
195 auto const& process_data = this->_process_data;
196 auto const& medium =
197 *process_data.media_map.getMedium(this->_element.getID());
198 auto const& liquid_phase = medium.phase("AqueousLiquid");
199
200 auto const& b =
201 process_data
202 .projected_specific_body_force_vectors[this->_element.getID()];
203
205
206 unsigned const n_integration_points =
207 this->_integration_method.getNumberOfPoints();
208
209 std::vector<GlobalDimVectorType> ip_flux_vector;
210 double average_velocity_norm = 0.0;
211 ip_flux_vector.reserve(n_integration_points);
212
213 auto const& Ns =
214 process_data.shape_matrix_cache
215 .template NsHigherOrder<typename ShapeFunction::MeshElement>();
216
217 for (unsigned ip(0); ip < n_integration_points; ip++)
218 {
219 auto const& ip_data = this->_ip_data[ip];
220 auto const& dNdx = ip_data.dNdx;
221 auto const& N = Ns[ip];
222 auto const& w = ip_data.integration_weight;
223
224 double p_at_xi = 0.;
225 NumLib::shapeFunctionInterpolate(local_p, N, p_at_xi);
226 double T_at_xi = 0.;
227 NumLib::shapeFunctionInterpolate(local_T, N, T_at_xi);
228
229 vars.temperature = T_at_xi;
230 vars.liquid_phase_pressure = p_at_xi;
231
232 vars.liquid_saturation = 1.0;
233
234 auto const porosity =
236 .template value<double>(vars, pos, t, dt);
237 vars.porosity = porosity;
238
239 // Use the fluid density model to compute the density
240 auto const fluid_density =
241 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
242 .template value<double>(vars, pos, t, dt);
243 vars.density = fluid_density;
244 auto const specific_heat_capacity_fluid =
246 .template value<double>(vars, pos, t, dt);
247
248 // Assemble mass matrix
249 local_M.noalias() += w *
250 this->getHeatEnergyCoefficient(
251 vars, porosity, fluid_density,
252 specific_heat_capacity_fluid, pos, t, dt) *
253 N.transpose() * N;
254
255 // Assemble Laplace matrix
256 auto const viscosity =
258 .template value<double>(vars, pos, t, dt);
259
260 auto const intrinsic_permeability =
263 .value(vars, pos, t, dt));
264
265 GlobalDimMatrixType const K_over_mu =
266 intrinsic_permeability / viscosity;
267 GlobalDimVectorType const velocity =
268 process_data.has_gravity
269 ? GlobalDimVectorType(-K_over_mu *
270 (dNdx * local_p - fluid_density * b))
271 : GlobalDimVectorType(-K_over_mu * dNdx * local_p);
272
273 GlobalDimMatrixType const thermal_conductivity_dispersivity =
274 this->getThermalConductivityDispersivity(
275 vars, fluid_density, specific_heat_capacity_fluid, velocity,
276 pos, t, dt);
277
278 local_K.noalias() +=
279 w * dNdx.transpose() * thermal_conductivity_dispersivity * dNdx;
280
281 ip_flux_vector.emplace_back(velocity * fluid_density *
282 specific_heat_capacity_fluid);
283 average_velocity_norm += velocity.norm();
284 }
285
287 process_data.stabilizer, this->_ip_data,
288 process_data.shape_matrix_cache, ip_flux_vector,
289 average_velocity_norm / static_cast<double>(n_integration_points),
290 local_K);
291}
292
293template <typename ShapeFunction, int GlobalDim>
294std::vector<double> const&
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
300{
301 assert(x.size() == dof_table.size());
302 auto const n_processes = dof_table.size();
303
304 std::vector<std::vector<GlobalIndexType>> indices_of_all_coupled_processes;
305 indices_of_all_coupled_processes.reserve(n_processes);
306 for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
307 {
308 auto const indices =
309 NumLib::getIndices(this->_element.getID(), *dof_table[process_id]);
310 assert(!indices.empty());
311 indices_of_all_coupled_processes.push_back(indices);
312 }
313 auto const local_xs =
314 getCoupledLocalSolutions(x, indices_of_all_coupled_processes);
315
316 return this->getIntPtDarcyVelocityLocal(t, local_xs, cache);
317}
318} // namespace HT
319} // namespace ProcessLib
void setElementID(std::size_t element_id)
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
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::vector< double > getCoupledLocalSolutions(std::vector< GlobalVector * > const &global_solutions, std::vector< std::vector< GlobalIndexType > > const &indices)