ncpboundaryratevector.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_NCP_BOUNDARY_RATE_VECTOR_HH
29#define EWOMS_NCP_BOUNDARY_RATE_VECTOR_HH
30
31#include "ncpproperties.hh"
32
34#include <opm/material/common/Valgrind.hpp>
35
36namespace Opm {
37
44template <class TypeTag>
45class NcpBoundaryRateVector : public GetPropType<TypeTag, Properties::RateVector>
46{
53
54 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
55 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
56 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
57 enum { conti0EqIdx = Indices::conti0EqIdx };
58 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
59
61
62 using Toolbox = Opm::MathToolbox<Evaluation>;
63
64public:
65 NcpBoundaryRateVector() : ParentType()
66 {}
67
72 NcpBoundaryRateVector(const Evaluation& value) : ParentType(value)
73 {}
74
81
85 template <class Context, class FluidState>
86 void setFreeFlow(const Context& context,
87 unsigned bfIdx,
88 unsigned timeIdx,
89 const FluidState& fluidState)
90 {
91 ExtensiveQuantities extQuants;
92 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
93 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
94 unsigned focusDofIdx = context.focusDofIndex();
95 unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
96
98 // advective fluxes of all components in all phases
100 (*this) = Evaluation(0.0);
101 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
102 Evaluation density;
103 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
104 if (focusDofIdx == interiorDofIdx)
105 density = fluidState.density(phaseIdx);
106 else
107 density = Opm::getValue(fluidState.density(phaseIdx));
108 }
109 else if (focusDofIdx == interiorDofIdx)
110 density = insideIntQuants.fluidState().density(phaseIdx);
111 else
112 density = Opm::getValue(insideIntQuants.fluidState().density(phaseIdx));
113
114 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
115 Evaluation molarity;
116 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
117 if (focusDofIdx == interiorDofIdx)
118 molarity = fluidState.molarity(phaseIdx, compIdx);
119 else
120 molarity = Opm::getValue(fluidState.molarity(phaseIdx, compIdx));
121 }
122 else if (focusDofIdx == interiorDofIdx)
123 molarity = insideIntQuants.fluidState().molarity(phaseIdx, compIdx);
124 else
125 molarity = Opm::getValue(insideIntQuants.fluidState().molarity(phaseIdx, compIdx));
126
127 // add advective flux of current component in current
128 // phase
129 (*this)[conti0EqIdx + compIdx] += extQuants.volumeFlux(phaseIdx)*molarity;
130 }
131
132 if (enableEnergy) {
133 Evaluation specificEnthalpy;
134 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
135 if (focusDofIdx == interiorDofIdx)
136 specificEnthalpy = fluidState.enthalpy(phaseIdx);
137 else
138 specificEnthalpy = Opm::getValue(fluidState.enthalpy(phaseIdx));
139 }
140 else if (focusDofIdx == interiorDofIdx)
141 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
142 else
143 specificEnthalpy = Opm::getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
144
145 Evaluation enthalpyRate = density*extQuants.volumeFlux(phaseIdx)*specificEnthalpy;
146 EnergyModule::addToEnthalpyRate(*this, enthalpyRate);
147 }
148 }
149
150 // thermal 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 }
157#endif
158 }
159
163 template <class Context, class FluidState>
164 void setInFlow(const Context& context,
165 unsigned bfIdx,
166 unsigned timeIdx,
167 const FluidState& fluidState)
168 {
169 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
170
171 // we only allow fluxes in the direction opposite to the outer unit normal
172 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
173 this->operator[](eqIdx) = Toolbox::min(0.0, this->operator[](eqIdx));
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 this->operator[](eqIdx) = Toolbox::max(0.0, this->operator[](eqIdx));
191 }
192 }
193
198 { (*this) = Evaluation(0.0); }
199};
200
201} // namespace Opm
202
203#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:50
Implements a boundary vector for the fully implicit compositional multi-phase NCP model.
Definition: ncpboundaryratevector.hh:46
NcpBoundaryRateVector(const Evaluation &value)
Definition: ncpboundaryratevector.hh:72
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: ncpboundaryratevector.hh:164
void setFreeFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: ncpboundaryratevector.hh:86
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: ncpboundaryratevector.hh:197
NcpBoundaryRateVector & operator=(const NcpBoundaryRateVector &value)=default
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: ncpboundaryratevector.hh:181
NcpBoundaryRateVector(const NcpBoundaryRateVector &value)=default
NcpBoundaryRateVector()
Definition: ncpboundaryratevector.hh:65
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:235
Declares the properties required for the NCP compositional multi-phase model.