OGS
TwoPhaseFlowWithPPMaterialProperties.cpp
Go to the documentation of this file.
1 
11 
12 #include <boost/math/special_functions/pow.hpp>
13 #include <utility>
14 
15 #include "BaseLib/Logging.h"
22 #include "MeshLib/Mesh.h"
23 #include "MeshLib/PropertyVector.h"
24 #include "ParameterLib/Parameter.h"
26 
27 namespace MaterialLib
28 {
29 namespace TwoPhaseFlowWithPP
30 {
32  MeshLib::PropertyVector<int> const* material_ids,
33  std::unique_ptr<MaterialLib::Fluid::FluidProperty>&& liquid_density,
34  std::unique_ptr<MaterialLib::Fluid::FluidProperty>&& liquid_viscosity,
35  std::unique_ptr<MaterialLib::Fluid::FluidProperty>&& gas_density,
36  std::unique_ptr<MaterialLib::Fluid::FluidProperty>&& gas_viscosity,
37  std::vector<std::unique_ptr<MaterialLib::PorousMedium::Permeability>>&&
38  intrinsic_permeability_models,
39  std::vector<std::unique_ptr<MaterialLib::PorousMedium::Porosity>>&&
40  porosity_models,
41  std::vector<std::unique_ptr<MaterialLib::PorousMedium::Storage>>&&
42  storage_models,
43  std::vector<std::unique_ptr<
45  capillary_pressure_models,
46  std::vector<
47  std::unique_ptr<MaterialLib::PorousMedium::RelativePermeability>>&&
48  relative_permeability_models)
49  : _material_ids(material_ids),
50  _liquid_density(std::move(liquid_density)),
51  _liquid_viscosity(std::move(liquid_viscosity)),
52  _gas_density(std::move(gas_density)),
53  _gas_viscosity(std::move(gas_viscosity)),
54  _intrinsic_permeability_models(std::move(intrinsic_permeability_models)),
55  _porosity_models(std::move(porosity_models)),
56  _storage_models(std::move(storage_models)),
57  _capillary_pressure_models(std::move(capillary_pressure_models)),
58  _relative_permeability_models(std::move(relative_permeability_models))
59 {
60  DBUG("Create material properties for Two-Phase flow with PP model.");
61 }
62 
64  const std::size_t element_id) const
65 {
66  if (_material_ids == nullptr || _material_ids->empty())
67  {
68  return 0;
69  }
70 
71  assert(element_id < _material_ids->size());
72  return (*_material_ids)[element_id];
73 }
74 
76  const double p, const double T) const
77 {
78  ArrayType vars;
79  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::T)] = T;
80  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::p)] = p;
81  return _liquid_density->getValue(vars);
82 }
83 
85  const double T) const
86 {
87  ArrayType vars;
88  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::T)] = T;
89  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::p)] = p;
90  return _gas_density->getValue(vars);
91 }
92 
94  const double p, const double T) const
95 {
96  ArrayType vars;
97  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::T)] = T;
98  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::p)] = p;
99  return _liquid_viscosity->getValue(vars);
100 }
101 
103  const double p, const double T) const
104 {
105  ArrayType vars;
106  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::T)] = T;
107  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::p)] = p;
108  return _gas_viscosity->getValue(vars);
109 }
110 
112  const int material_id, const double t,
113  const ParameterLib::SpatialPosition& pos, const int /*dim*/) const
114 {
115  return _intrinsic_permeability_models[material_id]->getValue(t, pos, 0, 0);
116 }
117 
119  const int material_id, const double t,
120  const ParameterLib::SpatialPosition& pos, const double /*p*/,
121  const double T, const double porosity_variable) const
122 {
123  return _porosity_models[material_id]->getValue(t, pos, porosity_variable,
124  T);
125 }
126 
128  const int material_id, const double /*t*/,
129  const ParameterLib::SpatialPosition& /*pos*/, const double /*p*/,
130  const double /*T*/, const double pc) const
131 {
132  return _capillary_pressure_models[material_id]->getSaturation(pc);
133 }
134 
136  const int material_id, const double /*t*/,
137  const ParameterLib::SpatialPosition& /*pos*/, const double /*p*/,
138  const double /*T*/, const double saturation) const
139 {
140  return _capillary_pressure_models[material_id]->getCapillaryPressure(
141  saturation);
142 }
143 
145  const int material_id, const double /*t*/,
146  const ParameterLib::SpatialPosition& /*pos*/, const double /*p*/,
147  const double /*T*/, const double saturation) const
148 {
149  const double dpcdsw =
150  _capillary_pressure_models[material_id]->getdPcdS(saturation);
151  return 1 / dpcdsw;
152 }
153 
155  const double /*t*/, const ParameterLib::SpatialPosition& /*pos*/,
156  const double /*p*/, const double /*T*/, const double saturation) const
157 {
158  if (saturation < 0.)
159  {
160  return 1.0;
161  }
162  if (saturation > 1)
163  {
164  return 0.0;
165  }
166  return boost::math::pow<3>(1 - saturation);
167 }
168 
170  const double /*t*/, const ParameterLib::SpatialPosition& /*pos*/,
171  const double /*p*/, const double /*T*/, const double saturation) const
172 {
173  if (saturation < 0)
174  {
175  return 0.0;
176  }
177  if (saturation > 1)
178  {
179  return 1.0;
180  }
181  return boost::math::pow<3>(saturation);
182 }
183 } // namespace TwoPhaseFlowWithPP
184 } // namespace MaterialLib
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
Definition of the Mesh class.
Definition of the PiecewiseLinearInterpolation class.
double getWetRelativePermeability(const double t, const ParameterLib::SpatialPosition &pos, const double p, const double T, const double saturation) const
double getPorosity(const int material_id, const double t, const ParameterLib::SpatialPosition &pos, const double p, const double T, const double porosity_variable) const
double getSaturationDerivative(const int material_id, const double t, const ParameterLib::SpatialPosition &pos, const double p, const double T, const double saturation) const
double getSaturation(const int material_id, const double t, const ParameterLib::SpatialPosition &pos, const double p, const double T, const double pc) const
std::vector< std::unique_ptr< MaterialLib::PorousMedium::Porosity > > const _porosity_models
TwoPhaseFlowWithPPMaterialProperties(MeshLib::PropertyVector< int > const *material_ids, std::unique_ptr< MaterialLib::Fluid::FluidProperty > &&liquid_density, std::unique_ptr< MaterialLib::Fluid::FluidProperty > &&liquid_viscosity, std::unique_ptr< MaterialLib::Fluid::FluidProperty > &&gas_density, std::unique_ptr< MaterialLib::Fluid::FluidProperty > &&gas_viscosity, std::vector< std::unique_ptr< MaterialLib::PorousMedium::Permeability >> &&intrinsic_permeability_models, std::vector< std::unique_ptr< MaterialLib::PorousMedium::Porosity >> &&porosity_models, std::vector< std::unique_ptr< MaterialLib::PorousMedium::Storage >> &&storage_models, std::vector< std::unique_ptr< MaterialLib::PorousMedium::CapillaryPressureSaturation >> &&capillary_pressure_models, std::vector< std::unique_ptr< MaterialLib::PorousMedium::RelativePermeability >> &&relative_permeability_models)
Eigen::MatrixXd getPermeability(const int material_id, const double t, const ParameterLib::SpatialPosition &pos, const int dim) const
double getCapillaryPressure(const int material_id, const double t, const ParameterLib::SpatialPosition &pos, const double p, const double T, const double saturation) const
double getNonwetRelativePermeability(const double t, const ParameterLib::SpatialPosition &pos, const double p, const double T, const double saturation) const
std::unique_ptr< MaterialLib::Fluid::FluidProperty > const _liquid_density
std::unique_ptr< MaterialLib::Fluid::FluidProperty > const _liquid_viscosity
std::vector< std::unique_ptr< MaterialLib::PorousMedium::CapillaryPressureSaturation > > const _capillary_pressure_models
std::vector< std::unique_ptr< MaterialLib::PorousMedium::Permeability > > const _intrinsic_permeability_models
std::unique_ptr< MaterialLib::Fluid::FluidProperty > const _gas_viscosity