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)
24 if (process_id == this->
_process_data.heat_transport_process_id)
32 local_K_data, local_b_data);
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)
46 auto const local_T_prev =
59 auto const& liquid_phase = medium.phase(
"AqueousLiquid");
60 auto const& solid_phase = medium.phase(
"Solid");
64 .projected_specific_body_force_vectors[this->
_element.getID()];
68 unsigned const n_integration_points =
72 process_data.shape_matrix_cache
73 .template NsHigherOrder<typename ShapeFunction::MeshElement>();
75 for (
unsigned ip(0); ip < n_integration_points; ip++)
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;
83 std::nullopt, this->
_element.getID(),
88 double p_int_pt = 0.0;
89 double T_int_pt = 0.0;
100 .template value<double>(vars, pos, t, dt);
101 auto const fluid_density =
103 .template value<double>(vars, pos, t, dt);
106 const double dfluid_density_dp =
108 .template dValue<double>(
113 auto const viscosity =
115 .template value<double>(vars, pos, t, dt);
119 auto const specific_storage =
121 .template value<double>(vars, pos, t, dt);
123 auto const intrinsic_permeability =
126 .value(vars, pos, t, dt));
128 intrinsic_permeability / viscosity;
133 (porosity * dfluid_density_dp / fluid_density + specific_storage) *
136 local_K.noalias() += w * dNdx.transpose() * K_over_mu * dNdx;
138 if (process_data.has_gravity)
141 w * fluid_density * dNdx.transpose() * K_over_mu * b;
144 if (!process_data.has_fluid_thermal_expansion)
151 auto const solid_thermal_expansion =
152 process_data.solid_thermal_expansion(t, pos)[0];
153 const double dfluid_density_dT =
156 .template dValue<double>(
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;
169template <
typename ShapeFunction,
int GlobalDim>
173 Eigen::VectorXd
const& local_x,
174 std::vector<double>& local_M_data,
175 std::vector<double>& local_K_data)
189 *process_data.media_map.getMedium(this->
_element.getID());
190 auto const& liquid_phase = medium.phase(
"AqueousLiquid");
194 .projected_specific_body_force_vectors[this->
_element.getID()];
198 unsigned const n_integration_points =
201 std::vector<GlobalDimVectorType> ip_flux_vector;
202 double average_velocity_norm = 0.0;
203 ip_flux_vector.reserve(n_integration_points);
206 process_data.shape_matrix_cache
207 .template NsHigherOrder<typename ShapeFunction::MeshElement>();
209 for (
unsigned ip(0); ip < n_integration_points; ip++)
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;
217 std::nullopt, this->
_element.getID(),
232 auto const porosity =
234 .template value<double>(vars, pos, t, dt);
238 auto const fluid_density =
240 .template value<double>(vars, pos, t, dt);
242 auto const specific_heat_capacity_fluid =
244 .template value<double>(vars, pos, t, dt);
247 local_M.noalias() += w *
249 vars, porosity, fluid_density,
250 specific_heat_capacity_fluid, pos, t, dt) *
254 auto const viscosity =
256 .template value<double>(vars, pos, t, dt);
258 auto const intrinsic_permeability =
261 .value(vars, pos, t, dt));
264 intrinsic_permeability / viscosity;
266 process_data.has_gravity
268 (dNdx * local_p - fluid_density * b))
273 vars, fluid_density, specific_heat_capacity_fluid, velocity,
277 w * dNdx.transpose() * thermal_conductivity_dispersivity * dNdx;
279 ip_flux_vector.emplace_back(velocity * fluid_density *
280 specific_heat_capacity_fluid);
281 average_velocity_norm += velocity.norm();
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),
291template <
typename ShapeFunction,
int GlobalDim>
292std::vector<double>
const&
295 std::vector<GlobalVector*>
const& x,
296 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
297 std::vector<double>& cache)
const
299 assert(x.size() == dof_table.size());
300 auto const n_processes = dof_table.size();
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)
308 assert(!indices.empty());
309 indices_of_all_coupled_processes.push_back(indices);
311 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
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)