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)
34 unsigned const ip = 0;
36 auto const& Ns = _process_data.shape_matrix_cache
37 .NsHigherOrder<
typename ShapeFunction::MeshElement>();
38 auto const& N = Ns[ip];
41 std::nullopt, _element.getID(),
46 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
50 .template value<double>(vars, pos, t, dt);
59 assert(permeability.rows() == GlobalDim || permeability.rows() == 1);
61 if (permeability.size() == 1)
63 assembleMatrixAndVector<IsotropicCalculator>(
64 t, dt, local_x, local_M_data, local_K_data, local_b_data);
68 assembleMatrixAndVector<AnisotropicCalculator>(
69 t, dt, local_x, local_M_data, local_K_data, local_b_data);
73template <
typename ShapeFunction,
int GlobalDim>
76 std::vector<double>
const& local_x)
const
80 double const dt = std::numeric_limits<double>::quiet_NaN();
84 auto const shape_matrices =
88 std::array{p_local_coords})[0];
94 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
95 auto const& fluid_phase = fluidPhase(medium);
99 double pressure = 0.0;
107 auto const viscosity =
109 .template value<double>(vars, pos, t, dt);
111 Eigen::Vector3d flux(0.0, 0.0, 0.0);
112 flux.head<GlobalDim>() =
113 -intrinsic_permeability / viscosity * shape_matrices.dNdx *
114 Eigen::Map<const NodalVectorType>(local_x.data(), local_x.size());
119template <
typename ShapeFunction,
int GlobalDim>
120template <
typename LaplacianGravityVelocityCalculator>
123 std::vector<double>
const& local_x,
124 std::vector<double>& local_M_data,
125 std::vector<double>& local_K_data,
126 std::vector<double>& local_b_data)
128 auto const local_matrix_size = local_x.size();
129 assert(local_matrix_size == ShapeFunction::NPOINTS);
132 local_M_data, local_matrix_size, local_matrix_size);
134 local_K_data, local_matrix_size, local_matrix_size);
136 local_b_data, local_matrix_size);
137 const auto local_p_vec =
140 unsigned const n_integration_points =
141 _integration_method.getNumberOfPoints();
146 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
147 auto const& fluid_phase = fluidPhase(medium);
154 _process_data.element_rotation_matrices[_element.getID()] *
155 _process_data.element_rotation_matrices[_element.getID()].transpose() *
156 _process_data.specific_body_force;
158 auto const& Ns = _process_data.shape_matrix_cache
159 .NsHigherOrder<
typename ShapeFunction::MeshElement>();
161 for (
unsigned ip = 0; ip < n_integration_points; ip++)
163 auto const& ip_data = _ip_data[ip];
164 auto const& N = Ns[ip];
166 phase_pressure = N.dot(local_p_vec);
168 std::nullopt, _element.getID(),
175 .template value<double>(vars, pos, t, dt);
177 auto const [fluid_density, viscosity] =
181 auto const porosity =
183 .template value<double>(vars, pos, t, dt);
184 auto const specific_storage =
186 .template value<double>(vars, pos, t, dt);
188 auto const get_drho_dp = [&]()
191 .template dValue<double>(vars, _process_data.phase_variable,
197 : specific_storage + porosity * get_drho_dp() / fluid_density;
199 double const scaling_factor =
204 local_M.noalias() += scaling_factor * storage * N.transpose() * N *
205 ip_data.integration_weight;
213 LaplacianGravityVelocityCalculator::calculateLaplacianAndGravityTerm(
214 local_K, local_b, ip_data, scaling_factor * permeability, viscosity,
215 fluid_density, projected_body_force_vector,
216 _process_data.has_gravity);
220template <
typename ShapeFunction,
int GlobalDim>
221template <
typename VelocityCacheType>
223 bool const is_scalar_permeability,
const double t,
const double dt,
224 std::vector<double>
const& local_x,
226 VelocityCacheType& darcy_velocity_at_ips)
const
228 if (is_scalar_permeability)
230 computeProjectedDarcyVelocity<IsotropicCalculator>(
231 t, dt, local_x, pos, darcy_velocity_at_ips);
235 computeProjectedDarcyVelocity<AnisotropicCalculator>(
236 t, dt, local_x, pos, darcy_velocity_at_ips);
240template <
typename ShapeFunction,
int GlobalDim>
241std::vector<double>
const&
244 std::vector<GlobalVector*>
const& x,
245 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_tables,
246 std::vector<double>& velocity_cache)
const
250 double const dt = std::numeric_limits<double>::quiet_NaN();
252 constexpr int process_id = 0;
255 auto const local_x = x[process_id]->get(indices);
256 auto const n_integration_points = _integration_method.getNumberOfPoints();
257 velocity_cache.clear();
260 unsigned const ip = 0;
262 auto const& Ns = _process_data.shape_matrix_cache
263 .NsHigherOrder<
typename ShapeFunction::MeshElement>();
264 auto const& N = Ns[ip];
267 std::nullopt, _element.getID(),
272 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
276 .template value<double>(vars, pos, t, dt);
284 assert(permeability.rows() == GlobalDim || permeability.rows() == 1);
286 bool const is_scalar_permeability = (permeability.size() == 1);
289 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
290 velocity_cache, GlobalDim, n_integration_points);
292 computeDarcyVelocity(is_scalar_permeability, t, dt, local_x, pos,
293 velocity_cache_vectors);
295 return velocity_cache;
298template <
typename ShapeFunction,
int GlobalDim>
299template <
typename LaplacianGravityVelocityCalculator,
300 typename VelocityCacheType>
303 const double t,
const double dt, std::vector<double>
const& local_x,
305 VelocityCacheType& darcy_velocity_at_ips)
const
307 auto const local_matrix_size = local_x.size();
308 assert(local_matrix_size == ShapeFunction::NPOINTS);
310 const auto local_p_vec =
313 unsigned const n_integration_points =
314 _integration_method.getNumberOfPoints();
316 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
317 auto const& fluid_phase = fluidPhase(medium);
324 _process_data.element_rotation_matrices[_element.getID()] *
325 _process_data.element_rotation_matrices[_element.getID()].transpose() *
326 _process_data.specific_body_force;
328 auto const& Ns = _process_data.shape_matrix_cache
329 .NsHigherOrder<
typename ShapeFunction::MeshElement>();
331 for (
unsigned ip = 0; ip < n_integration_points; ip++)
333 auto const& ip_data = _ip_data[ip];
334 auto const& N = Ns[ip];
336 phase_pressure = N.dot(local_p_vec);
338 std::nullopt, _element.getID(),
345 .template value<double>(vars, pos, t, dt);
346 auto const [fluid_density, viscosity] =
355 darcy_velocity_at_ips.col(ip) =
356 LaplacianGravityVelocityCalculator::calculateVelocity(
357 local_p_vec, ip_data, permeability, viscosity, fluid_density,
358 projected_body_force_vector, _process_data.has_gravity);
362template <
typename ShapeFunction,
int GlobalDim>
365 Eigen::Map<NodalMatrixType>& local_K,
366 Eigen::Map<NodalVectorType>& local_b,
369 double const mu,
double const rho_L,
372 const double K = permeability_with_density_factor(0, 0) / mu;
374 local_K.noalias() += fac * ip_data.
dNdx.transpose() * ip_data.
dNdx;
379 (fac * rho_L) * ip_data.
dNdx.transpose() * specific_body_force;
383template <
typename ShapeFunction,
int GlobalDim>
384Eigen::Matrix<double, GlobalDim, 1>
387 Eigen::Map<const NodalVectorType>
const& local_p,
391 bool const has_gravity)
393 const double K = permeability(0, 0) / mu;
399 velocity += (K * rho_L) * specific_body_force;
404template <
typename ShapeFunction,
int GlobalDim>
407 Eigen::Map<NodalMatrixType>& local_K,
408 Eigen::Map<NodalVectorType>& local_b,
411 double const mu,
double const rho_L,
415 local_K.noalias() += fac * ip_data.
dNdx.transpose() *
416 permeability_with_density_factor * ip_data.
dNdx;
420 local_b.noalias() += (fac * rho_L) * ip_data.
dNdx.transpose() *
421 permeability_with_density_factor *
426template <
typename ShapeFunction,
int GlobalDim>
427Eigen::Matrix<double, GlobalDim, 1>
430 Eigen::Map<const NodalVectorType>
const& local_p,
434 bool const has_gravity)
442 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)
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)