richardsboundaryratevector.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 This file is part of the Open Porous Media project (OPM).
5
6 OPM 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 2 of the License, or
9 (at your option) any later version.
10
11 OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_RICHARDS_BOUNDARY_RATE_VECTOR_HH
29#define EWOMS_RICHARDS_BOUNDARY_RATE_VECTOR_HH
30
31#include <opm/material/common/MathToolbox.hpp>
32#include <opm/material/common/Valgrind.hpp>
33
37
38#include <algorithm>
39
40namespace Opm {
41
47template <class TypeTag>
48class RichardsBoundaryRateVector : public GetPropType<TypeTag, Properties::RateVector>
49{
56
57 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
58 enum { contiEqIdx = Indices::contiEqIdx };
59 enum { liquidPhaseIdx = getPropValue<TypeTag, Properties::LiquidPhaseIndex>() };
60
61 using Toolbox = MathToolbox<Evaluation>;
62
63public:
69 RichardsBoundaryRateVector(const Evaluation& value)
70 : ParentType(value)
71 {}
72
79
83 template <class Context, class FluidState>
84 void setFreeFlow(const Context& context, unsigned bfIdx,
85 unsigned timeIdx, const FluidState& fluidState)
86 {
87 ExtensiveQuantities extQuants;
88 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
89 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
90 const unsigned focusDofIdx = context.focusDofIndex();
91 const unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
92
94 // advective fluxes of all components in all phases
96 (*this) = Evaluation(0.0);
97
98 const unsigned phaseIdx = liquidPhaseIdx;
99 Evaluation density;
100 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
101 if (focusDofIdx == interiorDofIdx) {
102 density = fluidState.density(phaseIdx);
103 }
104 else {
105 density = getValue(fluidState.density(phaseIdx));
106 }
107 }
108 else if (focusDofIdx == interiorDofIdx) {
109 density = insideIntQuants.fluidState().density(phaseIdx);
110 }
111 else {
112 density = getValue(insideIntQuants.fluidState().density(phaseIdx));
113 }
114
115 // add advective flux of current component in current
116 // phase
117 (*this)[contiEqIdx] += extQuants.volumeFlux(phaseIdx) * density;
118
119#ifndef NDEBUG
120 for (unsigned i = 0; i < numEq; ++i) {
121 Valgrind::CheckDefined((*this)[i]);
122 }
123 Valgrind::CheckDefined(*this);
124#endif
125 }
126
130 template <class Context, class FluidState>
131 void setInFlow(const Context& context,
132 unsigned bfIdx,
133 unsigned timeIdx,
134 const FluidState& fluidState)
135 {
136 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
137
138 // we only allow fluxes in the direction opposite to the outer unit normal
139 std::for_each(this->begin(), this->end(),
140 [](auto& val) { val = Toolbox::min(0.0, val); });
141 }
142
146 template <class Context, class FluidState>
147 void setOutFlow(const Context& context,
148 unsigned bfIdx,
149 unsigned timeIdx,
150 const FluidState& fluidState)
151 {
152 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
153
154 // we only allow fluxes in the same direction as the outer unit normal
155 std::for_each(this->begin(), this->end(),
156 [](auto& val) { val = Toolbox::max(0.0, val); });
157 }
158
163 { (*this) = Evaluation(0.0); }
164};
165
166} // namespace Opm
167
168#endif
Implements a boundary vector for the fully implicit Richards model.
Definition: richardsboundaryratevector.hh:49
RichardsBoundaryRateVector(const RichardsBoundaryRateVector &value)=default
void setFreeFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: richardsboundaryratevector.hh:84
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: richardsboundaryratevector.hh:147
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: richardsboundaryratevector.hh:162
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: richardsboundaryratevector.hh:131
RichardsBoundaryRateVector & operator=(const RichardsBoundaryRateVector &value)=default
RichardsBoundaryRateVector(const Evaluation &value)
Definition: richardsboundaryratevector.hh:69
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common properties required by the porous medium multi-phase models.
Definition: blackoilboundaryratevector.hh:39
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
Contains the property declarations for the Richards model.