version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_2P_SCV_SATURATION_RECONSTRUCTION_HH
14#define DUMUX_2P_SCV_SATURATION_RECONSTRUCTION_HH
15
17
18namespace Dumux {
19
29template<class DiscretizationMethod, bool enableReconstruction>
31{
32public:
46 template<class SpatialParams, class Element, class Scv, class ElemSol>
47 static typename ElemSol::PrimaryVariables::value_type
48 reconstructSn(const SpatialParams& spatialParams,
49 const Element& element,
50 const Scv& scv,
51 const ElemSol& elemSol,
52 typename ElemSol::PrimaryVariables::value_type sn)
53 { return sn; }
54};
55
57template<>
58class TwoPScvSaturationReconstruction<DiscretizationMethods::Box, /*enableReconstruction*/true>
59{
60public:
70 template<class SpatialParams, class Element, class Scv, class ElemSol>
71 static typename ElemSol::PrimaryVariables::value_type
72 reconstructSn(const SpatialParams& spatialParams,
73 const Element& element,
74 const Scv& scv,
75 const ElemSol& elemSol,
76 typename ElemSol::PrimaryVariables::value_type sn)
77 {
78 // if this dof doesn't lie on a material interface, simply return Sn
79 const auto& materialInterfaces = spatialParams.materialInterfaces();
80 if (!materialInterfaces.isOnMaterialInterface(scv))
81 return sn;
82
83 // compute capillary pressure using material parameters associated with the dof
84 const auto& interfacePcSw = materialInterfaces.pcSwAtDof(scv);
85 const auto pc = interfacePcSw.pc(/*ww=*/1.0 - sn);
86
87 // reconstruct by inverting the pc-sw curve
88 const auto& pcSw = spatialParams.fluidMatrixInteraction(element, scv, elemSol).pcSwCurve();
89 const auto pcMin = pcSw.endPointPc();
90
91 if (pc < pcMin && pcMin > 0.0)
92 return 0.0;
93 else
94 return 1.0 - pcSw.sw(pc);
95 }
96};
97
98} // end namespace Dumux
99
100#endif
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 nonwetting phase saturation in an scv.
Definition: saturationreconstruction.hh:72
Class that computes the nonwetting saturation in an scv from the saturation at the global degree of f...
Definition: saturationreconstruction.hh:31
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 nonwetting phase saturation in an scv.
Definition: saturationreconstruction.hh:48
The available discretization methods in Dumux.
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
Definition: adapt.hh:17