17template <
typename ShapeFunction,
int GlobalDim>
19 double const t,
double const dt, std::vector<double>
const& local_x,
20 std::vector<double>
const& ,
21 std::vector<double>& local_M_data, std::vector<double>& local_K_data,
22 std::vector<double>& local_b_data)
25 unsigned const ip = 0;
28 .NsHigherOrder<
typename ShapeFunction::MeshElement>();
29 auto const& N = Ns[ip];
41 .template value<double>(vars, pos, t, dt);
50 assert(permeability.rows() == GlobalDim || permeability.rows() == 1);
52 if (permeability.size() == 1)
55 t, dt, local_x, local_M_data, local_K_data, local_b_data);
60 t, dt, local_x, local_M_data, local_K_data, local_b_data);
64template <
typename ShapeFunction,
int GlobalDim>
67 std::vector<double>
const& local_x)
const
71 double const dt = std::numeric_limits<double>::quiet_NaN();
75 auto const shape_matrices =
79 std::array{p_local_coords})[0];
86 auto const& fluid_phase = fluidPhase(medium);
90 double pressure = 0.0;
98 auto const viscosity =
100 .template value<double>(vars, pos, t, dt);
102 Eigen::Vector3d flux(0.0, 0.0, 0.0);
103 flux.head<GlobalDim>() =
104 -intrinsic_permeability / viscosity * shape_matrices.dNdx *
105 Eigen::Map<const NodalVectorType>(local_x.data(), local_x.size());
110template <
typename ShapeFunction,
int GlobalDim>
111template <
typename LaplacianGravityVelocityCalculator>
114 std::vector<double>
const& local_x,
115 std::vector<double>& local_M_data,
116 std::vector<double>& local_K_data,
117 std::vector<double>& local_b_data)
119 auto const local_matrix_size = local_x.size();
120 assert(local_matrix_size == ShapeFunction::NPOINTS);
123 local_M_data, local_matrix_size, local_matrix_size);
125 local_K_data, local_matrix_size, local_matrix_size);
127 local_b_data, local_matrix_size);
128 const auto local_p_vec =
131 unsigned const n_integration_points =
138 auto const& fluid_phase = fluidPhase(medium);
150 .NsHigherOrder<
typename ShapeFunction::MeshElement>();
152 for (
unsigned ip = 0; ip < n_integration_points; ip++)
155 auto const& N = Ns[ip];
157 phase_pressure = N.dot(local_p_vec);
166 .template value<double>(vars, pos, t, dt);
168 auto const [fluid_density, viscosity] =
172 auto const porosity =
174 .template value<double>(vars, pos, t, dt);
175 auto const specific_storage =
177 .template value<double>(vars, pos, t, dt);
179 auto const get_drho_dp = [&]()
188 : specific_storage + porosity * get_drho_dp() / fluid_density;
190 double const scaling_factor =
195 local_M.noalias() += scaling_factor * storage * N.transpose() * N *
196 ip_data.integration_weight;
204 LaplacianGravityVelocityCalculator::calculateLaplacianAndGravityTerm(
205 local_K, local_b, ip_data, scaling_factor * permeability, viscosity,
206 fluid_density, projected_body_force_vector,
211template <
typename ShapeFunction,
int GlobalDim>
212template <
typename VelocityCacheType>
214 bool const is_scalar_permeability,
const double t,
const double dt,
215 std::vector<double>
const& local_x,
217 VelocityCacheType& darcy_velocity_at_ips)
const
219 if (is_scalar_permeability)
222 t, dt, local_x, pos, darcy_velocity_at_ips);
227 t, dt, local_x, pos, darcy_velocity_at_ips);
231template <
typename ShapeFunction,
int GlobalDim>
232std::vector<double>
const&
235 std::vector<GlobalVector*>
const& x,
236 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_tables,
237 std::vector<double>& velocity_cache)
const
241 double const dt = std::numeric_limits<double>::quiet_NaN();
243 constexpr int process_id = 0;
246 auto const local_x = x[process_id]->get(indices);
248 velocity_cache.clear();
251 unsigned const ip = 0;
254 .NsHigherOrder<
typename ShapeFunction::MeshElement>();
255 auto const& N = Ns[ip];
267 .template value<double>(vars, pos, t, dt);
275 assert(permeability.rows() == GlobalDim || permeability.rows() == 1);
277 bool const is_scalar_permeability = (permeability.size() == 1);
280 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
281 velocity_cache, GlobalDim, n_integration_points);
284 velocity_cache_vectors);
286 return velocity_cache;
289template <
typename ShapeFunction,
int GlobalDim>
290template <
typename LaplacianGravityVelocityCalculator,
291 typename VelocityCacheType>
294 const double t,
const double dt, std::vector<double>
const& local_x,
296 VelocityCacheType& darcy_velocity_at_ips)
const
298 auto const local_matrix_size = local_x.size();
299 assert(local_matrix_size == ShapeFunction::NPOINTS);
301 const auto local_p_vec =
304 unsigned const n_integration_points =
308 auto const& fluid_phase = fluidPhase(medium);
320 .NsHigherOrder<
typename ShapeFunction::MeshElement>();
322 for (
unsigned ip = 0; ip < n_integration_points; ip++)
325 auto const& N = Ns[ip];
327 phase_pressure = N.dot(local_p_vec);
336 .template value<double>(vars, pos, t, dt);
337 auto const [fluid_density, viscosity] =
346 darcy_velocity_at_ips.col(ip) =
347 LaplacianGravityVelocityCalculator::calculateVelocity(
348 local_p_vec, ip_data, permeability, viscosity, fluid_density,
353template <
typename ShapeFunction,
int GlobalDim>
356 Eigen::Map<NodalMatrixType>& local_K,
357 Eigen::Map<NodalVectorType>& local_b,
360 double const mu,
double const rho_L,
363 const double K = permeability_with_density_factor(0, 0) / mu;
365 local_K.noalias() += fac * ip_data.
dNdx.transpose() * ip_data.
dNdx;
370 (fac * rho_L) * ip_data.
dNdx.transpose() * specific_body_force;
374template <
typename ShapeFunction,
int GlobalDim>
375Eigen::Matrix<double, GlobalDim, 1>
378 Eigen::Map<const NodalVectorType>
const& local_p,
382 bool const has_gravity)
384 const double K = permeability(0, 0) / mu;
390 velocity += (K * rho_L) * specific_body_force;
395template <
typename ShapeFunction,
int GlobalDim>
398 Eigen::Map<NodalMatrixType>& local_K,
399 Eigen::Map<NodalVectorType>& local_b,
402 double const mu,
double const rho_L,
406 local_K.noalias() += fac * ip_data.
dNdx.transpose() *
407 permeability_with_density_factor * ip_data.
dNdx;
411 local_b.noalias() += (fac * rho_L) * ip_data.
dNdx.transpose() *
412 permeability_with_density_factor *
417template <
typename ShapeFunction,
int GlobalDim>
418Eigen::Matrix<double, GlobalDim, 1>
421 Eigen::Map<const NodalVectorType>
const& local_p,
425 bool const has_gravity)
433 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
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
const LiquidFlowData & _process_data
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
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
MeshLib::Element const & _element
NumLib::GenericIntegrationMethod const & _integration_method
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.
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< 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::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
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)