3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
saturationreconstruction.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_2P_SCV_SATURATION_RECONSTRUCTION_HH
26#define DUMUX_2P_SCV_SATURATION_RECONSTRUCTION_HH
27
29
30namespace Dumux {
31
41template<DiscretizationMethod M, bool enableReconstruction>
43{
44public:
58 template<class SpatialParams, class Element, class Scv, class ElemSol>
59 static typename ElemSol::PrimaryVariables::value_type
60 reconstructSn(const SpatialParams& spatialParams,
61 const Element& element,
62 const Scv& scv,
63 const ElemSol& elemSol,
64 typename ElemSol::PrimaryVariables::value_type Sn)
65 { return Sn; }
66};
67
69template<>
70class TwoPScvSaturationReconstruction<DiscretizationMethod::box, /*enableReconstruction*/true>
71{
72public:
82 template<class SpatialParams, class Element, class Scv, class ElemSol>
83 static typename ElemSol::PrimaryVariables::value_type
84 reconstructSn(const SpatialParams& spatialParams,
85 const Element& element,
86 const Scv& scv,
87 const ElemSol& elemSol,
88 typename ElemSol::PrimaryVariables::value_type Sn)
89 {
90 // if this dof doesn't lie on a material interface, simply return Sn
91 if (!spatialParams.materialInterfaceParams().isOnMaterialInterface(scv))
92 return Sn;
93
94 using MaterialLaw = typename SpatialParams::MaterialLaw;
95 // compute capillary pressure using material parameters associated with the dof
96 const auto& ifMaterialParams = spatialParams.materialInterfaceParams().getDofParams(scv);
97 const auto pc = MaterialLaw::pc(ifMaterialParams, /*Sw=*/1.0 - Sn);
98
99 // reconstruct by inverting the pc-sw curve
100 const auto& materialLawParams = spatialParams.materialLawParams(element, scv, elemSol);
101 const auto pcMin = MaterialLaw::endPointPc(materialLawParams);
102
103 if (pc < pcMin && pcMin > 0.0) return 0.0;
104 else return 1.0 - MaterialLaw::sw(materialLawParams, pc);
105 }
106};
107
108} // end namespace Dumux
109
110#endif
The available discretization methods in Dumux.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Definition: adapt.hh:29
Class that computes the non-wetting saturation in an scv from the saturation at the global degree of ...
Definition: saturationreconstruction.hh:43
static ElemSol::PrimaryVariables::value_type reconstructSn(const SpatialParams &spatialParams, const Element &element, const Scv &scv, const ElemSol &elemSol, typename ElemSol::PrimaryVariables::value_type Sn)
Compute the non-wetting phase saturation in an scv.
Definition: saturationreconstruction.hh:60
static ElemSol::PrimaryVariables::value_type reconstructSn(const SpatialParams &spatialParams, const Element &element, const Scv &scv, const ElemSol &elemSol, typename ElemSol::PrimaryVariables::value_type Sn)
Compute the non-wetting phase saturation in an scv.
Definition: saturationreconstruction.hh:84