pvsboundaryratevector.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_PVS_BOUNDARY_RATE_VECTOR_HH
29#define EWOMS_PVS_BOUNDARY_RATE_VECTOR_HH
30
31#include "pvsproperties.hh"
32
34#include <opm/material/common/Valgrind.hpp>
35
36namespace Opm {
37
45template <class TypeTag>
46class PvsBoundaryRateVector : public GetPropType<TypeTag, Properties::RateVector>
47{
54
55 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
56 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
57 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
58 enum { conti0EqIdx = Indices::conti0EqIdx };
59 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
60
62 using Toolbox = Opm::MathToolbox<Evaluation>;
63
64public:
66 : ParentType()
67 {}
68
73 PvsBoundaryRateVector(const Evaluation& value)
74 : ParentType(value)
75 {}
76
83
87 template <class Context, class FluidState>
88 void setFreeFlow(const Context& context, unsigned bfIdx, unsigned timeIdx, const FluidState& fluidState)
89 {
90 ExtensiveQuantities extQuants;
91 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
92 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
93 unsigned focusDofIdx = context.focusDofIndex();
94 unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
95
97 // advective fluxes of all components in all phases
99 (*this) = Evaluation(0.0);
100 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
101 Evaluation density;
102 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
103 if (focusDofIdx == interiorDofIdx)
104 density = fluidState.density(phaseIdx);
105 else
106 density = Opm::getValue(fluidState.density(phaseIdx));
107 }
108 else if (focusDofIdx == interiorDofIdx)
109 density = insideIntQuants.fluidState().density(phaseIdx);
110 else
111 density = Opm::getValue(insideIntQuants.fluidState().density(phaseIdx));
112
113 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
114 Evaluation molarity;
115 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
116 if (focusDofIdx == interiorDofIdx)
117 molarity = fluidState.molarity(phaseIdx, compIdx);
118 else
119 molarity = Opm::getValue(fluidState.molarity(phaseIdx, compIdx));
120 }
121 else if (focusDofIdx == interiorDofIdx)
122 molarity = insideIntQuants.fluidState().molarity(phaseIdx, compIdx);
123 else
124 molarity = Opm::getValue(insideIntQuants.fluidState().molarity(phaseIdx, compIdx));
125
126 // add advective flux of current component in current
127 // phase
128 (*this)[conti0EqIdx + compIdx] += extQuants.volumeFlux(phaseIdx)*molarity;
129 }
130
131 if (enableEnergy) {
132 Evaluation specificEnthalpy;
133 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
134 if (focusDofIdx == interiorDofIdx)
135 specificEnthalpy = fluidState.enthalpy(phaseIdx);
136 else
137 specificEnthalpy = Opm::getValue(fluidState.enthalpy(phaseIdx));
138 }
139 else if (focusDofIdx == interiorDofIdx)
140 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
141 else
142 specificEnthalpy = Opm::getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
143
144 Evaluation enthalpyRate = density*extQuants.volumeFlux(phaseIdx)*specificEnthalpy;
145 EnergyModule::addToEnthalpyRate(*this, enthalpyRate);
146 }
147 }
148
149 if (enableEnergy)
150 // heat conduction
151 EnergyModule::addToEnthalpyRate(*this, EnergyModule::thermalConductionRate(extQuants));
152
153#ifndef NDEBUG
154 for (unsigned i = 0; i < numEq; ++i)
155 Opm::Valgrind::CheckDefined((*this)[i]);
156#endif
157 }
158
162 template <class Context, class FluidState>
163 void setInFlow(const Context& context,
164 unsigned bfIdx,
165 unsigned timeIdx,
166 const FluidState& fluidState)
167 {
168 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
169
170 // we only allow fluxes in the direction opposite to the outer unit normal
171 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
172 Evaluation& val = this->operator[](eqIdx);
173 val = Toolbox::min(0.0, val);
174 }
175 }
176
180 template <class Context, class FluidState>
181 void setOutFlow(const Context& context,
182 unsigned bfIdx,
183 unsigned timeIdx,
184 const FluidState& fluidState)
185 {
186 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
187
188 // we only allow fluxes in the same direction as the outer unit normal
189 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
190 Evaluation& val = this->operator[](eqIdx);
191 val = Toolbox::max(0.0, val);
192 }
193 }
194
199 { (*this) = Evaluation(0.0); }
200};
201
202} // namespace Opm
203
204#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:48
Implements a rate vector on the boundary for the fully implicit compositional multi-phase primary var...
Definition: pvsboundaryratevector.hh:47
PvsBoundaryRateVector & operator=(const PvsBoundaryRateVector &value)=default
PvsBoundaryRateVector(const PvsBoundaryRateVector &value)=default
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: pvsboundaryratevector.hh:181
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: pvsboundaryratevector.hh:163
PvsBoundaryRateVector()
Definition: pvsboundaryratevector.hh:65
PvsBoundaryRateVector(const Evaluation &value)
Definition: pvsboundaryratevector.hh:73
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: pvsboundaryratevector.hh:198
void setFreeFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: pvsboundaryratevector.hh:88
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
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
Declares the properties required for the compositional multi-phase primary variable switching model.