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/Valgrind.hpp>
32#include <opm/material/constraintsolvers/NcpFlash.hpp>
33
35
36namespace Opm {
37
43template <class TypeTag>
44class RichardsBoundaryRateVector : public GetPropType<TypeTag, Properties::RateVector>
45{
52
53 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
54 enum { contiEqIdx = Indices::contiEqIdx };
55 enum { liquidPhaseIdx = getPropValue<TypeTag, Properties::LiquidPhaseIndex>() };
56
57 using Toolbox = Opm::MathToolbox<Evaluation>;
58
59public:
61 {}
62
67 RichardsBoundaryRateVector(const Evaluation& value)
68 : ParentType(value)
69 {}
70
77
81 template <class Context, class FluidState>
82 void setFreeFlow(const Context& context, unsigned bfIdx, unsigned timeIdx, const FluidState& fluidState)
83 {
84 ExtensiveQuantities extQuants;
85 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
86 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
87 unsigned focusDofIdx = context.focusDofIndex();
88 unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
89
91 // advective fluxes of all components in all phases
93 (*this) = Evaluation(0.0);
94
95 unsigned phaseIdx = liquidPhaseIdx;
96 Evaluation density;
97 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
98 if (focusDofIdx == interiorDofIdx)
99 density = fluidState.density(phaseIdx);
100 else
101 density = Opm::getValue(fluidState.density(phaseIdx));
102 }
103 else if (focusDofIdx == interiorDofIdx)
104 density = insideIntQuants.fluidState().density(phaseIdx);
105 else
106 density = Opm::getValue(insideIntQuants.fluidState().density(phaseIdx));
107
108 // add advective flux of current component in current
109 // phase
110 (*this)[contiEqIdx] += extQuants.volumeFlux(phaseIdx) * density;
111
112#ifndef NDEBUG
113 for (unsigned i = 0; i < numEq; ++i) {
114 Opm::Valgrind::CheckDefined((*this)[i]);
115 }
116 Opm::Valgrind::CheckDefined(*this);
117#endif
118 }
119
123 template <class Context, class FluidState>
124 void setInFlow(const Context& context,
125 unsigned bfIdx,
126 unsigned timeIdx,
127 const FluidState& fluidState)
128 {
129 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
130
131 // we only allow fluxes in the direction opposite to the outer unit normal
132 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
133 Evaluation& val = this->operator[](eqIdx);
134 val = Toolbox::min(0.0, val);
135 }
136 }
137
141 template <class Context, class FluidState>
142 void setOutFlow(const Context& context,
143 unsigned bfIdx,
144 unsigned timeIdx,
145 const FluidState& fluidState)
146 {
147 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
148
149 // we only allow fluxes in the same direction as the outer unit normal
150 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
151 Evaluation& val = this->operator[](eqIdx);
152 val = Toolbox::max(0.0, val);
153 }
154 }
155
160 { (*this) = Evaluation(0.0); }
161};
162
163} // namespace Opm
164
165#endif
Implements a boundary vector for the fully implicit Richards model.
Definition: richardsboundaryratevector.hh:45
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:82
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: richardsboundaryratevector.hh:142
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: richardsboundaryratevector.hh:159
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: richardsboundaryratevector.hh:124
RichardsBoundaryRateVector()
Definition: richardsboundaryratevector.hh:60
RichardsBoundaryRateVector & operator=(const RichardsBoundaryRateVector &value)=default
RichardsBoundaryRateVector(const Evaluation &value)
Definition: richardsboundaryratevector.hh:67
Definition: blackoilboundaryratevector.hh:37
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:242