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
37
38#include <algorithm>
39
40namespace Opm {
41
48template <class TypeTag>
49class ImmiscibleBoundaryRateVector : public GetPropType<TypeTag, Properties::RateVector>
50{
57
58 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
59 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
60 enum { conti0EqIdx = Indices::conti0EqIdx };
61 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
62
63 using Toolbox = MathToolbox<Evaluation>;
65
66public:
68
75 ImmiscibleBoundaryRateVector(const Evaluation& value)
76 : ParentType(value)
77 {}
78
85
87
99 template <class Context, class FluidState>
100 void setFreeFlow(const Context& context, unsigned bfIdx, unsigned timeIdx, const FluidState& fluidState)
101 {
102 ExtensiveQuantities extQuants;
103 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
104 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
105 const unsigned focusDofIdx = context.focusDofIndex();
106 const unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
107
109 // advective fluxes of all components in all phases
111 (*this) = Evaluation(0.0);
112 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
113 const auto& pBoundary = fluidState.pressure(phaseIdx);
114 const Evaluation& pInside = insideIntQuants.fluidState().pressure(phaseIdx);
115
116 // mass conservation
117 Evaluation density;
118 if (pBoundary > pInside) {
119 if (focusDofIdx == interiorDofIdx) {
120 density = fluidState.density(phaseIdx);
121 }
122 else {
123 density = getValue(fluidState.density(phaseIdx));
124 }
125 }
126 else if (focusDofIdx == interiorDofIdx) {
127 density = insideIntQuants.fluidState().density(phaseIdx);
128 }
129 else {
130 density = getValue(insideIntQuants.fluidState().density(phaseIdx));
131 }
132
133 Valgrind::CheckDefined(density);
134 Valgrind::CheckDefined(extQuants.volumeFlux(phaseIdx));
135
136 (*this)[conti0EqIdx + phaseIdx] += extQuants.volumeFlux(phaseIdx) * density;
137
138 // energy conservation
139 if constexpr (enableEnergy) {
140 Evaluation specificEnthalpy;
141 if (pBoundary > pInside) {
142 if (focusDofIdx == interiorDofIdx) {
143 specificEnthalpy = fluidState.enthalpy(phaseIdx);
144 }
145 else {
146 specificEnthalpy = getValue(fluidState.enthalpy(phaseIdx));
147 }
148 }
149 else if (focusDofIdx == interiorDofIdx) {
150 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
151 }
152 else {
153 specificEnthalpy = getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
154 }
155
156 const Evaluation enthalpyRate = density * extQuants.volumeFlux(phaseIdx) * specificEnthalpy;
157 EnergyModule::addToEnthalpyRate(*this, enthalpyRate);
158 }
159 }
160
161 // thermal conduction
162 EnergyModule::addToEnthalpyRate(*this, EnergyModule::thermalConductionRate(extQuants));
163
164#ifndef NDEBUG
165 for (unsigned i = 0; i < numEq; ++i) {
166 Valgrind::CheckDefined((*this)[i]);
167 }
168 Valgrind::CheckDefined(*this);
169#endif
170 }
171
181 template <class Context, class FluidState>
182 void setInFlow(const Context& context,
183 unsigned bfIdx,
184 unsigned timeIdx,
185 const FluidState& fluidState)
186 {
187 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
188
189 // we only allow fluxes in the direction opposite to the outer unit normal
190 std::for_each(this->begin(), this->end(),
191 [](auto& val) { val = Toolbox::min(0.0, val); });
192 }
193
203 template <class Context, class FluidState>
204 void setOutFlow(const Context& context,
205 unsigned bfIdx,
206 unsigned timeIdx,
207 const FluidState& fluidState)
208 {
209 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
210
211 // we only allow fluxes in the same direction as the outer unit normal
212 std::for_each(this->begin(), this->end(),
213 [](auto& val) { val = Toolbox::max(0.0, val); });
214 }
215
220 { (*this) = Evaluation(0.0); }
221};
222
223} // namespace Opm
224
225#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:54
Implements a boundary vector for the fully implicit multi-phase model which assumes immiscibility.
Definition: immiscibleboundaryratevector.hh:50
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: immiscibleboundaryratevector.hh:219
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:100
ImmiscibleBoundaryRateVector(const Evaluation &value)
Constructor that assigns all entries to a scalar value.
Definition: immiscibleboundaryratevector.hh:75
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: immiscibleboundaryratevector.hh:204
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: immiscibleboundaryratevector.hh:182
ImmiscibleBoundaryRateVector(const ImmiscibleBoundaryRateVector &value)=default
Copy constructor.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
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