immiscibleboundaryratevector.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_IMMISCIBLE_BOUNDARY_RATE_VECTOR_HH
29#define EWOMS_IMMISCIBLE_BOUNDARY_RATE_VECTOR_HH
30
31#include <opm/material/common/Valgrind.hpp>
32#include <opm/material/constraintsolvers/NcpFlash.hpp>
33
35
36namespace Opm {
37
44template <class TypeTag>
45class ImmiscibleBoundaryRateVector : public GetPropType<TypeTag, Properties::RateVector>
46{
53
54 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
55 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
56 enum { conti0EqIdx = Indices::conti0EqIdx };
57 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
58
59 using Toolbox = Opm::MathToolbox<Evaluation>;
61
62public:
64 : ParentType()
65 {}
66
73 ImmiscibleBoundaryRateVector(const Evaluation& value)
74 : ParentType(value)
75 {}
76
83
85
97 template <class Context, class FluidState>
98 void setFreeFlow(const Context& context, unsigned bfIdx, unsigned timeIdx, const FluidState& fluidState)
99 {
100 ExtensiveQuantities extQuants;
101 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
102 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
103 unsigned focusDofIdx = context.focusDofIndex();
104 unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
105
107 // advective fluxes of all components in all phases
109 (*this) = Evaluation(0.0);
110 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
111 const auto& pBoundary = fluidState.pressure(phaseIdx);
112 const Evaluation& pInside = insideIntQuants.fluidState().pressure(phaseIdx);
113
114 // mass conservation
115 Evaluation density;
116 if (pBoundary > pInside) {
117 if (focusDofIdx == interiorDofIdx)
118 density = fluidState.density(phaseIdx);
119 else
120 density = Opm::getValue(fluidState.density(phaseIdx));
121 }
122 else if (focusDofIdx == interiorDofIdx)
123 density = insideIntQuants.fluidState().density(phaseIdx);
124 else
125 density = Opm::getValue(insideIntQuants.fluidState().density(phaseIdx));
126
127 Opm::Valgrind::CheckDefined(density);
128 Opm::Valgrind::CheckDefined(extQuants.volumeFlux(phaseIdx));
129
130 (*this)[conti0EqIdx + phaseIdx] += extQuants.volumeFlux(phaseIdx)*density;
131
132 // energy conservation
133 if (enableEnergy) {
134 Evaluation specificEnthalpy;
135 if (pBoundary > pInside) {
136 if (focusDofIdx == interiorDofIdx)
137 specificEnthalpy = fluidState.enthalpy(phaseIdx);
138 else
139 specificEnthalpy = Opm::getValue(fluidState.enthalpy(phaseIdx));
140 }
141 else if (focusDofIdx == interiorDofIdx)
142 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
143 else
144 specificEnthalpy = Opm::getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
145
146 Evaluation enthalpyRate = density*extQuants.volumeFlux(phaseIdx)*specificEnthalpy;
147 EnergyModule::addToEnthalpyRate(*this, enthalpyRate);
148 }
149 }
150
151 // thermal conduction
152 EnergyModule::addToEnthalpyRate(*this, EnergyModule::thermalConductionRate(extQuants));
153
154#ifndef NDEBUG
155 for (unsigned i = 0; i < numEq; ++i)
156 Opm::Valgrind::CheckDefined((*this)[i]);
157 Opm::Valgrind::CheckDefined(*this);
158#endif
159 }
160
170 template <class Context, class FluidState>
171 void setInFlow(const Context& context,
172 unsigned bfIdx,
173 unsigned timeIdx,
174 const FluidState& fluidState)
175 {
176 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
177
178 // we only allow fluxes in the direction opposite to the outer unit normal
179 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
180 Evaluation& val = this->operator[](eqIdx);
181 val = Toolbox::min(0.0, val);
182 }
183 }
184
194 template <class Context, class FluidState>
195 void setOutFlow(const Context& context,
196 unsigned bfIdx,
197 unsigned timeIdx,
198 const FluidState& fluidState)
199 {
200 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
201
202 // we only allow fluxes in the same direction as the outer unit normal
203 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
204 Evaluation& val = this->operator[](eqIdx);
205 val = Toolbox::max(0.0, val);
206 }
207 }
208
213 { (*this) = Evaluation(0.0); }
214};
215
216} // namespace Opm
217
218#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:50
Implements a boundary vector for the fully implicit multi-phase model which assumes immiscibility.
Definition: immiscibleboundaryratevector.hh:46
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: immiscibleboundaryratevector.hh:212
ImmiscibleBoundaryRateVector & operator=(const ImmiscibleBoundaryRateVector &value)=default
void setFreeFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: immiscibleboundaryratevector.hh:98
ImmiscibleBoundaryRateVector(const Evaluation &value)
Constructor that assigns all entries to a scalar value.
Definition: immiscibleboundaryratevector.hh:73
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: immiscibleboundaryratevector.hh:195
ImmiscibleBoundaryRateVector()
Definition: immiscibleboundaryratevector.hh:63
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: immiscibleboundaryratevector.hh:171
ImmiscibleBoundaryRateVector(const ImmiscibleBoundaryRateVector &value)=default
Copy constructor.
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:235