19template <
typename ShapeFunction,
int GlobalDim>
21 double const t,
double const dt, Eigen::VectorXd
const& local_x,
22 Eigen::VectorXd
const& local_x_prev,
int const process_id,
23 std::vector<double>& local_M_data, std::vector<double>& local_K_data,
24 std::vector<double>& local_b_data)
26 if (process_id == this->
_process_data.heat_transport_process_id)
34 local_K_data, local_b_data);
37template <
typename ShapeFunction,
int GlobalDim>
39 double const t,
double const dt, Eigen::VectorXd
const& local_x,
40 Eigen::VectorXd
const& local_x_prev, std::vector<double>& local_M_data,
41 std::vector<double>& local_K_data, std::vector<double>& local_b_data)
48 auto const local_T_prev =
61 auto const& liquid_phase =
63 auto const& solid_phase =
68 .projected_specific_body_force_vectors[this->
_element.getID()];
72 unsigned const n_integration_points =
76 process_data.shape_matrix_cache
77 .template NsHigherOrder<typename ShapeFunction::MeshElement>();
79 for (
unsigned ip(0); ip < n_integration_points; ip++)
81 auto const& ip_data = this->
_ip_data[ip];
82 auto const& dNdx = ip_data.dNdx;
83 auto const& N = Ns[ip];
84 auto const& w = ip_data.integration_weight;
87 std::nullopt, this->
_element.getID(),
92 double p_int_pt = 0.0;
93 double T_int_pt = 0.0;
102 auto const porosity =
104 .template value<double>(vars, pos, t, dt);
105 auto const fluid_density =
107 .template value<double>(vars, pos, t, dt);
110 const double dfluid_density_dp =
112 .template dValue<double>(
117 auto const viscosity =
119 .template value<double>(vars, pos, t, dt);
123 auto const specific_storage =
125 .template value<double>(vars, pos, t, dt);
127 auto const intrinsic_permeability =
130 .value(vars, pos, t, dt));
132 intrinsic_permeability / viscosity;
137 (porosity * dfluid_density_dp / fluid_density + specific_storage) *
140 local_K.noalias() += w * dNdx.transpose() * K_over_mu * dNdx;
142 if (process_data.has_gravity)
145 w * fluid_density * dNdx.transpose() * K_over_mu * b;
148 if (!process_data.has_fluid_thermal_expansion)
155 auto const solid_thermal_expansion =
156 process_data.solid_thermal_expansion(t, pos)[0];
157 const double dfluid_density_dT =
160 .template dValue<double>(
163 double const Tdot_int_pt = (T_int_pt - local_T_prev.dot(N)) / dt;
164 auto const biot_constant = process_data.biot_constant(t, pos)[0];
165 const double eff_thermal_expansion =
166 3.0 * (biot_constant - porosity) * solid_thermal_expansion -
167 porosity * dfluid_density_dT / fluid_density;
168 local_b.noalias() += eff_thermal_expansion * Tdot_int_pt * w * N;
173template <
typename ShapeFunction,
int GlobalDim>
177 Eigen::VectorXd
const& local_x,
178 std::vector<double>& local_M_data,
179 std::vector<double>& local_K_data)
193 *process_data.media_map.getMedium(this->
_element.getID());
194 auto const& liquid_phase =
199 .projected_specific_body_force_vectors[this->
_element.getID()];
203 unsigned const n_integration_points =
206 std::vector<GlobalDimVectorType> ip_flux_vector;
207 double average_velocity_norm = 0.0;
208 ip_flux_vector.reserve(n_integration_points);
211 process_data.shape_matrix_cache
212 .template NsHigherOrder<typename ShapeFunction::MeshElement>();
214 for (
unsigned ip(0); ip < n_integration_points; ip++)
216 auto const& ip_data = this->
_ip_data[ip];
217 auto const& dNdx = ip_data.dNdx;
218 auto const& N = Ns[ip];
219 auto const& w = ip_data.integration_weight;
222 std::nullopt, this->
_element.getID(),
237 auto const porosity =
239 .template value<double>(vars, pos, t, dt);
243 auto const fluid_density =
245 .template value<double>(vars, pos, t, dt);
247 auto const specific_heat_capacity_fluid =
249 .template value<double>(vars, pos, t, dt);
252 local_M.noalias() += w *
254 vars, porosity, fluid_density,
255 specific_heat_capacity_fluid, pos, t, dt) *
259 auto const viscosity =
261 .template value<double>(vars, pos, t, dt);
263 auto const intrinsic_permeability =
266 .value(vars, pos, t, dt));
269 intrinsic_permeability / viscosity;
271 process_data.has_gravity
273 (dNdx * local_p - fluid_density * b))
278 vars, fluid_density, specific_heat_capacity_fluid, velocity,
282 w * dNdx.transpose() * thermal_conductivity_dispersivity * dNdx;
284 ip_flux_vector.emplace_back(velocity * fluid_density *
285 specific_heat_capacity_fluid);
286 average_velocity_norm += velocity.norm();
290 process_data.stabilizer, this->_ip_data,
291 process_data.shape_matrix_cache, ip_flux_vector,
292 average_velocity_norm /
static_cast<double>(n_integration_points),
296template <
typename ShapeFunction,
int GlobalDim>
297std::vector<double>
const&
300 std::vector<GlobalVector*>
const& x,
301 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
302 std::vector<double>& cache)
const
304 assert(x.size() == dof_table.size());
305 auto const n_processes = dof_table.size();
307 std::vector<std::vector<GlobalIndexType>> indices_of_all_coupled_processes;
308 indices_of_all_coupled_processes.reserve(n_processes);
309 for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
313 assert(!indices.empty());
314 indices_of_all_coupled_processes.push_back(indices);
316 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)