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)
32 if (process_id == this->
_process_data.heat_transport_process_id)
40 local_K_data, local_b_data);
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)
54 auto const local_T_prev =
67 auto const& liquid_phase = medium.phase(
"AqueousLiquid");
68 auto const& solid_phase = medium.phase(
"Solid");
72 .projected_specific_body_force_vectors[this->
_element.getID()];
76 unsigned const n_integration_points =
80 process_data.shape_matrix_cache
81 .template NsHigherOrder<typename ShapeFunction::MeshElement>();
83 for (
unsigned ip(0); ip < n_integration_points; ip++)
85 auto const& ip_data = this->
_ip_data[ip];
86 auto const& dNdx = ip_data.dNdx;
87 auto const& N = Ns[ip];
88 auto const& w = ip_data.integration_weight;
91 std::nullopt, this->
_element.getID(),
96 double p_int_pt = 0.0;
97 double T_int_pt = 0.0;
106 auto const porosity =
108 .template value<double>(vars, pos, t, dt);
109 auto const fluid_density =
111 .template value<double>(vars, pos, t, dt);
114 const double dfluid_density_dp =
116 .template dValue<double>(
121 auto const viscosity =
123 .template value<double>(vars, pos, t, dt);
127 auto const specific_storage =
129 .template value<double>(vars, pos, t, dt);
131 auto const intrinsic_permeability =
134 .value(vars, pos, t, dt));
136 intrinsic_permeability / viscosity;
141 (porosity * dfluid_density_dp / fluid_density + specific_storage) *
144 local_K.noalias() += w * dNdx.transpose() * K_over_mu * dNdx;
146 if (process_data.has_gravity)
149 w * fluid_density * dNdx.transpose() * K_over_mu * b;
152 if (!process_data.has_fluid_thermal_expansion)
159 auto const solid_thermal_expansion =
160 process_data.solid_thermal_expansion(t, pos)[0];
161 const double dfluid_density_dT =
164 .template dValue<double>(
167 double const Tdot_int_pt = (T_int_pt - local_T_prev.dot(N)) / dt;
168 auto const biot_constant = process_data.biot_constant(t, pos)[0];
169 const double eff_thermal_expansion =
170 3.0 * (biot_constant - porosity) * solid_thermal_expansion -
171 porosity * dfluid_density_dT / fluid_density;
172 local_b.noalias() += eff_thermal_expansion * Tdot_int_pt * w * N;
177template <
typename ShapeFunction,
int GlobalDim>
181 Eigen::VectorXd
const& local_x,
182 std::vector<double>& local_M_data,
183 std::vector<double>& local_K_data)
197 *process_data.media_map.getMedium(this->
_element.getID());
198 auto const& liquid_phase = medium.phase(
"AqueousLiquid");
202 .projected_specific_body_force_vectors[this->
_element.getID()];
206 unsigned const n_integration_points =
209 std::vector<GlobalDimVectorType> ip_flux_vector;
210 double average_velocity_norm = 0.0;
211 ip_flux_vector.reserve(n_integration_points);
214 process_data.shape_matrix_cache
215 .template NsHigherOrder<typename ShapeFunction::MeshElement>();
217 for (
unsigned ip(0); ip < n_integration_points; ip++)
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;
225 std::nullopt, this->
_element.getID(),
240 auto const porosity =
242 .template value<double>(vars, pos, t, dt);
246 auto const fluid_density =
248 .template value<double>(vars, pos, t, dt);
250 auto const specific_heat_capacity_fluid =
252 .template value<double>(vars, pos, t, dt);
255 local_M.noalias() += w *
257 vars, porosity, fluid_density,
258 specific_heat_capacity_fluid, pos, t, dt) *
262 auto const viscosity =
264 .template value<double>(vars, pos, t, dt);
266 auto const intrinsic_permeability =
269 .value(vars, pos, t, dt));
272 intrinsic_permeability / viscosity;
274 process_data.has_gravity
276 (dNdx * local_p - fluid_density * b))
281 vars, fluid_density, specific_heat_capacity_fluid, velocity,
285 w * dNdx.transpose() * thermal_conductivity_dispersivity * dNdx;
287 ip_flux_vector.emplace_back(velocity * fluid_density *
288 specific_heat_capacity_fluid);
289 average_velocity_norm += velocity.norm();
293 process_data.stabilizer, this->_ip_data,
294 process_data.shape_matrix_cache, ip_flux_vector,
295 average_velocity_norm /
static_cast<double>(n_integration_points),
299template <
typename ShapeFunction,
int GlobalDim>
300std::vector<double>
const&
303 std::vector<GlobalVector*>
const& x,
304 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
305 std::vector<double>& cache)
const
307 assert(x.size() == dof_table.size());
308 auto const n_processes = dof_table.size();
310 std::vector<std::vector<GlobalIndexType>> indices_of_all_coupled_processes;
311 indices_of_all_coupled_processes.reserve(n_processes);
312 for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
316 assert(!indices.empty());
317 indices_of_all_coupled_processes.push_back(indices);
319 auto const local_xs =
double liquid_phase_pressure
std::vector< double > const & getIntPtDarcyVelocityLocal(const double t, std::vector< double > const &local_x, std::vector< double > &cache) const
NumLib::GenericIntegrationMethod const & _integration_method
static const int temperature_index
static const int temperature_size
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)
static const int pressure_size
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)
MeshLib::Element const & _element
HTProcessData const & _process_data
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
static const int pressure_index
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
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)