33template <
typename ShapeFunction,
int GlobalDim>
49 std::size_t
const local_matrix_size,
51 bool const is_axially_symmetric,
53 const unsigned dof_per_node)
61 assert(local_matrix_size == ShapeFunction::NPOINTS * dof_per_node);
62 (void)local_matrix_size;
65 unsigned const n_integration_points =
67 _ip_data.reserve(n_integration_points);
74 auto const shape_matrices =
76 GlobalDim>(element, is_axially_symmetric,
79 for (
unsigned ip = 0; ip < n_integration_points; ip++)
82 shape_matrices[ip].dNdx,
84 shape_matrices[ip].integralMeasure *
85 shape_matrices[ip].detJ * aperture_size);
90 const unsigned integration_point)
const override
93 typename ShapeFunction::MeshElement>()[integration_point];
96 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
103 std::vector<double>
const& local_x)
const override
108 auto const shape_matrices =
112 std::array{pnt_local_coords})[0];
120 double T_int_pt = 0.0;
121 double p_int_pt = 0.0;
130 auto const& liquid_phase = medium.
phase(
"AqueousLiquid");
134 double const dt = std::numeric_limits<double>::quiet_NaN();
138 .value(vars, pos, t, dt));
142 .template value<double>(vars, pos, t, dt);
145 auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
146 &local_x[local_x.size() / 2], ShapeFunction::NPOINTS);
148 -K_over_mu * shape_matrices.dNdx * p_nodal_values;
155 .template value<double>(vars, pos, t, dt);
158 [this->_element.getID()];
159 q += K_over_mu * rho_w * b;
162 Eigen::Vector3d flux;
163 flux.head<GlobalDim>() = q;
172 std::vector<IntegrationPointData<GlobalDimNodalMatrixType>>
_ip_data;
176 const double fluid_density,
const double specific_heat_capacity_fluid,
182 auto const& solid_phase = medium.
phase(
"Solid");
184 auto const specific_heat_capacity_solid =
188 .template value<double>(vars, pos, t, dt);
190 auto const solid_density =
192 .template value<double>(vars, pos, t, dt);
194 return solid_density * specific_heat_capacity_solid * (1 - porosity) +
195 fluid_density * specific_heat_capacity_fluid * porosity;
200 const double fluid_density,
const double specific_heat_capacity_fluid,
208 auto thermal_conductivity =
213 .value(vars, pos, t, dt));
215 auto const thermal_dispersivity_transversal =
218 thermal_transversal_dispersivity)
219 .
template value<double>();
221 auto const thermal_dispersivity_longitudinal =
224 thermal_longitudinal_dispersivity)
225 .
template value<double>();
230 return thermal_conductivity +
231 fluid_density * specific_heat_capacity_fluid *
234 GlobalDimMatrixType::Zero(GlobalDim, GlobalDim),
235 velocity, 0 , thermal_dispersivity_transversal,
236 thermal_dispersivity_longitudinal);
240 const double t, std::vector<double>
const& local_x,
241 std::vector<double>& cache)
const
243 std::vector<double> local_p{
246 std::vector<double> local_T{
250 auto const n_integration_points =
255 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
256 cache, GlobalDim, n_integration_points);
263 auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
264 &local_p[0], ShapeFunction::NPOINTS);
268 auto const& liquid_phase = medium.
phase(
"AqueousLiquid");
274 for (
unsigned ip = 0; ip < n_integration_points; ++ip)
277 auto const& dNdx = ip_data.dNdx;
278 auto const& N = Ns[ip];
280 pos.setIntegrationPoint(ip);
282 double T_int_pt = 0.0;
283 double p_int_pt = 0.0;
292 double const dt = std::numeric_limits<double>::quiet_NaN();
295 .value(vars, pos, t, dt));
300 .template value<double>(vars, pos, t, dt);
303 cache_mat.col(ip).noalias() = -K_over_mu * dNdx * p_nodal_values;
310 .template value<double>(vars, pos, t, dt);
315 cache_mat.col(ip).noalias() += K_over_mu * rho_w * b;
Medium * getMedium(std::size_t element_id)
Phase const & phase(std::size_t index) const
Property const & property(PropertyType const &p) const
double liquid_phase_pressure
std::size_t getID() const
Returns the ID of the element.
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
unsigned getNumberOfPoints() const
auto const & NsHigherOrder() const
void setElementID(std::size_t element_id)
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
std::vector< double > const & getIntPtDarcyVelocityLocal(const double t, std::vector< double > const &local_x, std::vector< double > &cache) const
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
NumLib::GenericIntegrationMethod const & _integration_method
typename ShapeMatricesType::GlobalDimNodalMatrixType GlobalDimNodalMatrixType
static const int temperature_index
typename ShapeMatricesType::NodalVectorType NodalVectorType
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)
HTFEM(MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, HTProcessData const &process_data, const unsigned dof_per_node)
static const int pressure_size
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
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)
Eigen::Vector3d getFlux(MathLib::Point3d const &pnt_local_coords, double const t, std::vector< double > const &local_x) const override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
MeshLib::Element const & _element
HTProcessData const & _process_data
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
static const int pressure_index
typename ShapeMatricesType::NodalRowVectorType NodalRowVectorType
Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > computeShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, PointContainer const &points)
Eigen::MatrixXd computeHydrodynamicDispersion(NumericalStabilization const &stabilizer, std::size_t const element_id, Eigen::MatrixXd const &pore_diffusion_coefficient, Eigen::VectorXd const &velocity, double const porosity, double const solute_dispersivity_transverse, double const solute_dispersivity_longitudinal)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
ParameterLib::Parameter< double > const & aperture_size
MaterialPropertyLib::MaterialSpatialDistributionMap media_map
NumLib::NumericalStabilization stabilizer
NumLib::ShapeMatrixCache shape_matrix_cache
caches for each mesh element type the shape matrix
std::vector< Eigen::VectorXd > const projected_specific_body_force_vectors
Projected specific body force vector: R * R^T * b.