25template <
typename ShapeFunction,
int GlobalDim>
27 double const t,
double const dt, std::vector<double>
const& local_x,
28 std::vector<double>
const& ,
29 std::vector<double>& local_M_data, std::vector<double>& local_K_data,
30 std::vector<double>& local_b_data)
35 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
39 .template value<double>(vars, pos, t, dt);
48 assert(permeability.rows() == GlobalDim || permeability.rows() == 1);
50 if (permeability.size() == 1)
52 assembleMatrixAndVector<IsotropicCalculator>(
53 t, dt, local_x, local_M_data, local_K_data, local_b_data);
57 assembleMatrixAndVector<AnisotropicCalculator>(
58 t, dt, local_x, local_M_data, local_K_data, local_b_data);
62template <
typename ShapeFunction,
int GlobalDim>
65 std::vector<double>
const& local_x)
const
69 double const dt = std::numeric_limits<double>::quiet_NaN();
73 auto const shape_matrices =
77 std::array{p_local_coords})[0];
83 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
84 auto const& fluid_phase = fluidPhase(medium);
88 double pressure = 0.0;
96 auto const viscosity =
98 .template value<double>(vars, pos, t, dt);
100 Eigen::Vector3d flux(0.0, 0.0, 0.0);
101 flux.head<GlobalDim>() =
102 -intrinsic_permeability / viscosity * shape_matrices.dNdx *
103 Eigen::Map<const NodalVectorType>(local_x.data(), local_x.size());
108template <
typename ShapeFunction,
int GlobalDim>
109template <
typename LaplacianGravityVelocityCalculator>
112 std::vector<double>
const& local_x,
113 std::vector<double>& local_M_data,
114 std::vector<double>& local_K_data,
115 std::vector<double>& local_b_data)
117 auto const local_matrix_size = local_x.size();
118 assert(local_matrix_size == ShapeFunction::NPOINTS);
121 local_M_data, local_matrix_size, local_matrix_size);
123 local_K_data, local_matrix_size, local_matrix_size);
125 local_b_data, local_matrix_size);
126 const auto local_p_vec =
129 unsigned const n_integration_points =
130 _integration_method.getNumberOfPoints();
135 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
136 auto const& fluid_phase = fluidPhase(medium);
144 .template value<double>(vars, pos, t, dt);
147 _process_data.element_rotation_matrices[_element.getID()] *
148 _process_data.element_rotation_matrices[_element.getID()].transpose() *
149 _process_data.specific_body_force;
151 auto const& Ns = _process_data.shape_matrix_cache
152 .NsHigherOrder<
typename ShapeFunction::MeshElement>();
154 for (
unsigned ip = 0; ip < n_integration_points; ip++)
156 auto const& ip_data = _ip_data[ip];
157 auto const& N = Ns[ip];
159 phase_pressure = N.dot(local_p_vec);
161 auto const [fluid_density, viscosity] =
164 auto const porosity =
166 .template value<double>(vars, pos, t, dt);
167 auto const specific_storage =
169 .template value<double>(vars, pos, t, dt);
171 auto const get_drho_dp = [&]()
174 .template dValue<double>(vars, _process_data.phase_variable,
180 : specific_storage + porosity * get_drho_dp() / fluid_density;
182 double const scaling_factor =
187 local_M.noalias() += scaling_factor * storage * N.transpose() * N *
188 ip_data.integration_weight;
197 LaplacianGravityVelocityCalculator::calculateLaplacianAndGravityTerm(
198 local_K, local_b, ip_data, scaling_factor * permeability, viscosity,
199 fluid_density, projected_body_force_vector,
200 _process_data.has_gravity);
204template <
typename ShapeFunction,
int GlobalDim>
205template <
typename VelocityCacheType>
207 bool const is_scalar_permeability,
const double t,
const double dt,
208 std::vector<double>
const& local_x,
210 VelocityCacheType& darcy_velocity_at_ips)
const
212 if (is_scalar_permeability)
214 computeProjectedDarcyVelocity<IsotropicCalculator>(
215 t, dt, local_x, pos, darcy_velocity_at_ips);
219 computeProjectedDarcyVelocity<AnisotropicCalculator>(
220 t, dt, local_x, pos, darcy_velocity_at_ips);
224template <
typename ShapeFunction,
int GlobalDim>
225std::vector<double>
const&
228 std::vector<GlobalVector*>
const& x,
229 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_tables,
230 std::vector<double>& velocity_cache)
const
234 double const dt = std::numeric_limits<double>::quiet_NaN();
236 constexpr int process_id = 0;
239 auto const local_x = x[process_id]->get(indices);
240 auto const n_integration_points = _integration_method.getNumberOfPoints();
241 velocity_cache.clear();
246 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
250 .template value<double>(vars, pos, t, dt);
258 assert(permeability.rows() == GlobalDim || permeability.rows() == 1);
260 bool const is_scalar_permeability = (permeability.size() == 1);
263 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
264 velocity_cache, GlobalDim, n_integration_points);
266 computeDarcyVelocity(is_scalar_permeability, t, dt, local_x, pos,
267 velocity_cache_vectors);
269 return velocity_cache;
272template <
typename ShapeFunction,
int GlobalDim>
273template <
typename LaplacianGravityVelocityCalculator,
274 typename VelocityCacheType>
277 const double t,
const double dt, std::vector<double>
const& local_x,
279 VelocityCacheType& darcy_velocity_at_ips)
const
281 auto const local_matrix_size = local_x.size();
282 assert(local_matrix_size == ShapeFunction::NPOINTS);
284 const auto local_p_vec =
287 unsigned const n_integration_points =
288 _integration_method.getNumberOfPoints();
290 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
291 auto const& fluid_phase = fluidPhase(medium);
299 .template value<double>(vars, pos, t, dt);
302 _process_data.element_rotation_matrices[_element.getID()] *
303 _process_data.element_rotation_matrices[_element.getID()].transpose() *
304 _process_data.specific_body_force;
306 auto const& Ns = _process_data.shape_matrix_cache
307 .NsHigherOrder<
typename ShapeFunction::MeshElement>();
309 for (
unsigned ip = 0; ip < n_integration_points; ip++)
311 auto const& ip_data = _ip_data[ip];
312 auto const& N = Ns[ip];
314 phase_pressure = N.dot(local_p_vec);
316 auto const [fluid_density, viscosity] =
324 darcy_velocity_at_ips.col(ip) =
325 LaplacianGravityVelocityCalculator::calculateVelocity(
326 local_p_vec, ip_data, permeability, viscosity, fluid_density,
327 projected_body_force_vector, _process_data.has_gravity);
331template <
typename ShapeFunction,
int GlobalDim>
334 Eigen::Map<NodalMatrixType>& local_K,
335 Eigen::Map<NodalVectorType>& local_b,
338 double const mu,
double const rho_L,
341 const double K = permeability_with_density_factor(0, 0) / mu;
343 local_K.noalias() += fac * ip_data.
dNdx.transpose() * ip_data.
dNdx;
348 (fac * rho_L) * ip_data.
dNdx.transpose() * specific_body_force;
352template <
typename ShapeFunction,
int GlobalDim>
353Eigen::Matrix<double, GlobalDim, 1>
356 Eigen::Map<const NodalVectorType>
const& local_p,
360 bool const has_gravity)
362 const double K = permeability(0, 0) / mu;
368 velocity += (K * rho_L) * specific_body_force;
373template <
typename ShapeFunction,
int GlobalDim>
376 Eigen::Map<NodalMatrixType>& local_K,
377 Eigen::Map<NodalVectorType>& local_b,
380 double const mu,
double const rho_L,
384 local_K.noalias() += fac * ip_data.
dNdx.transpose() *
385 permeability_with_density_factor * ip_data.
dNdx;
389 local_b.noalias() += (fac * rho_L) * ip_data.
dNdx.transpose() *
390 permeability_with_density_factor *
395template <
typename ShapeFunction,
int GlobalDim>
396Eigen::Matrix<double, GlobalDim, 1>
399 Eigen::Map<const NodalVectorType>
const& local_p,
403 bool const has_gravity)
411 velocity += (rho_L / mu) * permeability * specific_body_force;
double gas_phase_pressure
double liquid_phase_pressure
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
Eigen::Vector3d getFlux(MathLib::Point3d const &p_local_coords, double const t, std::vector< double > const &local_x) const override
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
void computeProjectedDarcyVelocity(const double t, const double dt, std::vector< double > const &local_x, ParameterLib::SpatialPosition const &pos, VelocityCacheType &darcy_velocity_at_ips) const
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
void assembleMatrixAndVector(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
void computeDarcyVelocity(bool const is_scalar_permeability, const double t, const double dt, std::vector< double > const &local_x, ParameterLib::SpatialPosition const &pos, VelocityCacheType &darcy_velocity_at_ips) const
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 > &velocity_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< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > computeShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, PointContainer const &points)
std::tuple< double, double > getFluidDensityAndViscosity(double const t, double const dt, ParameterLib::SpatialPosition const &pos, MaterialPropertyLib::Phase const &fluid_phase, MaterialPropertyLib::VariableArray &vars)
GlobalDimNodalMatrixType const dNdx
double const integration_weight
static void calculateLaplacianAndGravityTerm(Eigen::Map< NodalMatrixType > &local_K, Eigen::Map< NodalVectorType > &local_b, IntegrationPointData< GlobalDimNodalMatrixType > const &ip_data, GlobalDimMatrixType const &permeability_with_density_factor, double const mu, double const rho_L, GlobalDimVectorType const &specific_body_force, bool const has_gravity)
static Eigen::Matrix< double, GlobalDim, 1 > calculateVelocity(Eigen::Map< const NodalVectorType > const &local_p, IntegrationPointData< GlobalDimNodalMatrixType > const &ip_data, GlobalDimMatrixType const &permeability, double const mu, double const rho_L, GlobalDimVectorType const &specific_body_force, bool const has_gravity)
static void calculateLaplacianAndGravityTerm(Eigen::Map< NodalMatrixType > &local_K, Eigen::Map< NodalVectorType > &local_b, IntegrationPointData< GlobalDimNodalMatrixType > const &ip_data, GlobalDimMatrixType const &permeability_with_density_factor, double const mu, double const rho_L, GlobalDimVectorType const &specific_body_force, bool const has_gravity)
static Eigen::Matrix< double, GlobalDim, 1 > calculateVelocity(Eigen::Map< const NodalVectorType > const &local_p, IntegrationPointData< GlobalDimNodalMatrixType > const &ip_data, GlobalDimMatrixType const &permeability, double const mu, double const rho_L, GlobalDimVectorType const &specific_body_force, bool const has_gravity)