26template <
typename ShapeFunction,
int GlobalDim>
28 double const t,
double const dt, std::vector<double>
const& local_x,
29 std::vector<double>
const& ,
30 std::vector<double>& local_M_data, std::vector<double>& local_K_data,
31 std::vector<double>& local_b_data)
36 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
40 .template value<double>(vars, pos, t, dt);
49 assert(permeability.rows() == GlobalDim || permeability.rows() == 1);
51 if (permeability.size() == 1)
53 assembleMatrixAndVector<IsotropicCalculator>(
54 t, dt, local_x, local_M_data, local_K_data, local_b_data);
58 assembleMatrixAndVector<AnisotropicCalculator>(
59 t, dt, local_x, local_M_data, local_K_data, local_b_data);
63template <
typename ShapeFunction,
int GlobalDim>
66 std::vector<double>
const& local_x)
const
70 double const dt = std::numeric_limits<double>::quiet_NaN();
74 auto const shape_matrices =
78 std::array{p_local_coords})[0];
84 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
85 auto const& fluid_phase = fluidPhase(medium);
89 double pressure = 0.0;
97 auto const viscosity =
99 .template value<double>(vars, pos, t, dt);
101 Eigen::Vector3d flux(0.0, 0.0, 0.0);
102 flux.head<GlobalDim>() =
103 -intrinsic_permeability / viscosity * shape_matrices.dNdx *
104 Eigen::Map<const NodalVectorType>(local_x.data(), local_x.size());
109template <
typename ShapeFunction,
int GlobalDim>
110template <
typename LaplacianGravityVelocityCalculator>
113 std::vector<double>
const& local_x,
114 std::vector<double>& local_M_data,
115 std::vector<double>& local_K_data,
116 std::vector<double>& local_b_data)
118 auto const local_matrix_size = local_x.size();
119 assert(local_matrix_size == ShapeFunction::NPOINTS);
122 local_M_data, local_matrix_size, local_matrix_size);
124 local_K_data, local_matrix_size, local_matrix_size);
126 local_b_data, local_matrix_size);
127 const auto local_p_vec =
130 unsigned const n_integration_points =
131 _integration_method.getNumberOfPoints();
136 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
137 auto const& fluid_phase = fluidPhase(medium);
145 .template value<double>(vars, pos, t, dt);
148 _process_data.element_rotation_matrices[_element.getID()] *
149 _process_data.element_rotation_matrices[_element.getID()].transpose() *
150 _process_data.specific_body_force;
152 auto const& Ns = _process_data.shape_matrix_cache
153 .NsHigherOrder<
typename ShapeFunction::MeshElement>();
155 for (
unsigned ip = 0; ip < n_integration_points; ip++)
157 auto const& ip_data = _ip_data[ip];
158 auto const& N = Ns[ip];
160 phase_pressure = N.dot(local_p_vec);
162 auto const [fluid_density, viscosity] =
166 auto const porosity =
168 .template value<double>(vars, pos, t, dt);
169 auto const specific_storage =
171 .template value<double>(vars, pos, t, dt);
173 auto const get_drho_dp = [&]()
176 .template dValue<double>(vars, _process_data.phase_variable,
182 : specific_storage + porosity * get_drho_dp() / fluid_density;
184 double const scaling_factor =
189 local_M.noalias() += scaling_factor * storage * N.transpose() * N *
190 ip_data.integration_weight;
198 LaplacianGravityVelocityCalculator::calculateLaplacianAndGravityTerm(
199 local_K, local_b, ip_data, scaling_factor * permeability, viscosity,
200 fluid_density, projected_body_force_vector,
201 _process_data.has_gravity);
205template <
typename ShapeFunction,
int GlobalDim>
206template <
typename VelocityCacheType>
208 bool const is_scalar_permeability,
const double t,
const double dt,
209 std::vector<double>
const& local_x,
211 VelocityCacheType& darcy_velocity_at_ips)
const
213 if (is_scalar_permeability)
215 computeProjectedDarcyVelocity<IsotropicCalculator>(
216 t, dt, local_x, pos, darcy_velocity_at_ips);
220 computeProjectedDarcyVelocity<AnisotropicCalculator>(
221 t, dt, local_x, pos, darcy_velocity_at_ips);
225template <
typename ShapeFunction,
int GlobalDim>
226std::vector<double>
const&
229 std::vector<GlobalVector*>
const& x,
230 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_tables,
231 std::vector<double>& velocity_cache)
const
235 double const dt = std::numeric_limits<double>::quiet_NaN();
237 constexpr int process_id = 0;
240 auto const local_x = x[process_id]->get(indices);
241 auto const n_integration_points = _integration_method.getNumberOfPoints();
242 velocity_cache.clear();
247 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
251 .template value<double>(vars, pos, t, dt);
259 assert(permeability.rows() == GlobalDim || permeability.rows() == 1);
261 bool const is_scalar_permeability = (permeability.size() == 1);
264 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
265 velocity_cache, GlobalDim, n_integration_points);
267 computeDarcyVelocity(is_scalar_permeability, t, dt, local_x, pos,
268 velocity_cache_vectors);
270 return velocity_cache;
273template <
typename ShapeFunction,
int GlobalDim>
274template <
typename LaplacianGravityVelocityCalculator,
275 typename VelocityCacheType>
278 const double t,
const double dt, std::vector<double>
const& local_x,
280 VelocityCacheType& darcy_velocity_at_ips)
const
282 auto const local_matrix_size = local_x.size();
283 assert(local_matrix_size == ShapeFunction::NPOINTS);
285 const auto local_p_vec =
288 unsigned const n_integration_points =
289 _integration_method.getNumberOfPoints();
291 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
292 auto const& fluid_phase = fluidPhase(medium);
300 .template value<double>(vars, pos, t, dt);
303 _process_data.element_rotation_matrices[_element.getID()] *
304 _process_data.element_rotation_matrices[_element.getID()].transpose() *
305 _process_data.specific_body_force;
307 auto const& Ns = _process_data.shape_matrix_cache
308 .NsHigherOrder<
typename ShapeFunction::MeshElement>();
310 for (
unsigned ip = 0; ip < n_integration_points; ip++)
312 auto const& ip_data = _ip_data[ip];
313 auto const& N = Ns[ip];
315 phase_pressure = N.dot(local_p_vec);
317 auto const [fluid_density, viscosity] =
326 darcy_velocity_at_ips.col(ip) =
327 LaplacianGravityVelocityCalculator::calculateVelocity(
328 local_p_vec, ip_data, permeability, viscosity, fluid_density,
329 projected_body_force_vector, _process_data.has_gravity);
333template <
typename ShapeFunction,
int GlobalDim>
336 Eigen::Map<NodalMatrixType>& local_K,
337 Eigen::Map<NodalVectorType>& local_b,
340 double const mu,
double const rho_L,
343 const double K = permeability_with_density_factor(0, 0) / mu;
345 local_K.noalias() += fac * ip_data.
dNdx.transpose() * ip_data.
dNdx;
350 (fac * rho_L) * ip_data.
dNdx.transpose() * specific_body_force;
354template <
typename ShapeFunction,
int GlobalDim>
355Eigen::Matrix<double, GlobalDim, 1>
358 Eigen::Map<const NodalVectorType>
const& local_p,
362 bool const has_gravity)
364 const double K = permeability(0, 0) / mu;
370 velocity += (K * rho_L) * specific_body_force;
375template <
typename ShapeFunction,
int GlobalDim>
378 Eigen::Map<NodalMatrixType>& local_K,
379 Eigen::Map<NodalVectorType>& local_b,
382 double const mu,
double const rho_L,
386 local_K.noalias() += fac * ip_data.
dNdx.transpose() *
387 permeability_with_density_factor * ip_data.
dNdx;
391 local_b.noalias() += (fac * rho_L) * ip_data.
dNdx.transpose() *
392 permeability_with_density_factor *
397template <
typename ShapeFunction,
int GlobalDim>
398Eigen::Matrix<double, GlobalDim, 1>
401 Eigen::Map<const NodalVectorType>
const& local_p,
405 bool const has_gravity)
413 velocity += (rho_L / mu) * permeability * specific_body_force;
double gas_phase_pressure
double liquid_phase_pressure
void setElementID(std::size_t element_id)
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
std::tuple< double, double > getFluidDensityAndViscosity(double const t, double const dt, ParameterLib::SpatialPosition const &pos, MaterialPropertyLib::Phase const &fluid_phase, MaterialPropertyLib::VariableArray &vars)
It computes fluid density and viscosity for single phase flow model.
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)
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)