OGS
HydrodynamicDispersion.h
Go to the documentation of this file.
1
10#pragma once
11
12#include <Eigen/Core>
13
15
16namespace NumLib
17{
18namespace detail
19{
20
21inline Eigen::MatrixXd getHydrodynamicDispersion(
22 Eigen::MatrixXd const& pore_diffusion_coefficient,
23 Eigen::VectorXd const& velocity,
24 double const porosity,
25 double const solute_dispersivity_transverse,
26 double const solute_dispersivity_longitudinal)
27{
28 double const velocity_magnitude = velocity.norm();
29 if (velocity_magnitude == 0.0)
30 {
31 return porosity * pore_diffusion_coefficient;
32 }
33
34 auto const dim = velocity.size();
35 Eigen::MatrixXd const& I(Eigen::MatrixXd::Identity(dim, dim));
36 return porosity * pore_diffusion_coefficient +
37 solute_dispersivity_transverse * velocity_magnitude * I +
38 (solute_dispersivity_longitudinal - solute_dispersivity_transverse) /
39 velocity_magnitude * velocity * velocity.transpose();
40}
41
43 IsotropicDiffusionStabilization const& stabilizer,
44 std::size_t const element_id,
45 Eigen::MatrixXd const& pore_diffusion_coefficient,
46 Eigen::VectorXd const& velocity,
47 double const porosity,
48 double const solute_dispersivity_transverse,
49 double const solute_dispersivity_longitudinal)
50{
51 double const velocity_magnitude = velocity.norm();
52 if (velocity_magnitude == 0.0)
53 {
54 return porosity * pore_diffusion_coefficient;
55 }
56
57 double const artificial_diffusion =
58 stabilizer.computeArtificialDiffusion(element_id, velocity_magnitude);
59
60 auto const dim = velocity.size();
61 Eigen::MatrixXd const& I(Eigen::MatrixXd::Identity(dim, dim));
62 return porosity * pore_diffusion_coefficient +
63 (solute_dispersivity_transverse * velocity_magnitude +
64 artificial_diffusion) *
65 I +
66 (solute_dispersivity_longitudinal - solute_dispersivity_transverse) /
67 velocity_magnitude * velocity * velocity.transpose();
68}
69} // namespace detail
70
71inline Eigen::MatrixXd computeHydrodynamicDispersion(
72 NumericalStabilization const& stabilizer,
73 std::size_t const element_id,
74 Eigen::MatrixXd const& pore_diffusion_coefficient,
75 Eigen::VectorXd const& velocity,
76 double const porosity,
77 double const solute_dispersivity_transverse,
78 double const solute_dispersivity_longitudinal)
79{
80 return std::visit(
81 [&](auto&& stabilizer)
82 {
83 using Stabilizer = std::decay_t<decltype(stabilizer)>;
84 if constexpr (std::is_same_v<Stabilizer,
86 {
88 stabilizer,
89 element_id,
90 pore_diffusion_coefficient,
91 velocity,
92 porosity,
93 solute_dispersivity_transverse,
94 solute_dispersivity_longitudinal);
95 }
96
98 pore_diffusion_coefficient,
99 velocity,
100 porosity,
101 solute_dispersivity_transverse,
102 solute_dispersivity_longitudinal);
103 },
104 stabilizer);
105}
106
107} // namespace NumLib
double computeArtificialDiffusion(std::size_t const element_id, double const velocity_norm) const
Eigen::MatrixXd getHydrodynamicDispersion(Eigen::MatrixXd const &pore_diffusion_coefficient, Eigen::VectorXd const &velocity, double const porosity, double const solute_dispersivity_transverse, double const solute_dispersivity_longitudinal)
Eigen::MatrixXd getHydrodynamicDispersionWithArtificialDiffusion(IsotropicDiffusionStabilization 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)
std::variant< NoStabilization, IsotropicDiffusionStabilization, FullUpwind, FluxCorrectedTransport > NumericalStabilization
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)