flashboundaryratevector.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_FLASH_BOUNDARY_RATE_VECTOR_HH
29#define EWOMS_FLASH_BOUNDARY_RATE_VECTOR_HH
30
32#include <opm/material/common/Valgrind.hpp>
33
34namespace Opm {
35
43template <class TypeTag>
44class FlashBoundaryRateVector : public GetPropType<TypeTag, Properties::RateVector>
45{
52
53 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
54 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
55 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
56 enum { conti0EqIdx = Indices::conti0EqIdx };
57 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
58
60 using Toolbox = Opm::MathToolbox<Evaluation>;
61
62public:
63 FlashBoundaryRateVector() : ParentType()
64 {}
65
70 FlashBoundaryRateVector(const Evaluation& value) : ParentType(value)
71 {}
72
79
83 template <class Context, class FluidState>
84 void setFreeFlow(const Context& context,
85 unsigned bfIdx,
86 unsigned timeIdx,
87 const FluidState& fluidState)
88 {
89 ExtensiveQuantities extQuants;
90 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
91 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
92 unsigned focusDofIdx = context.focusDofIndex();
93 unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
94
96 // advective fluxes of all components in all phases
98 (*this) = Evaluation(0.0);
99 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
100 Evaluation density;
101 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
102 if (focusDofIdx == interiorDofIdx)
103 density = fluidState.density(phaseIdx);
104 else
105 density = Opm::getValue(fluidState.density(phaseIdx));
106 }
107 else if (focusDofIdx == interiorDofIdx)
108 density = insideIntQuants.fluidState().density(phaseIdx);
109 else
110 density = Opm::getValue(insideIntQuants.fluidState().density(phaseIdx));
111
112 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
113 Evaluation molarity;
114 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
115 if (focusDofIdx == interiorDofIdx)
116 molarity = fluidState.molarity(phaseIdx, compIdx);
117 else
118 molarity = Opm::getValue(fluidState.molarity(phaseIdx, compIdx));
119 }
120 else if (focusDofIdx == interiorDofIdx)
121 molarity = insideIntQuants.fluidState().molarity(phaseIdx, compIdx);
122 else
123 molarity = Opm::getValue(insideIntQuants.fluidState().molarity(phaseIdx, compIdx));
124
125 // add advective flux of current component in current
126 // phase
127 (*this)[conti0EqIdx + compIdx] += extQuants.volumeFlux(phaseIdx)*molarity;
128 }
129
130 if (enableEnergy) {
131 Evaluation specificEnthalpy;
132 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
133 if (focusDofIdx == interiorDofIdx)
134 specificEnthalpy = fluidState.enthalpy(phaseIdx);
135 else
136 specificEnthalpy = Opm::getValue(fluidState.enthalpy(phaseIdx));
137 }
138 else if (focusDofIdx == interiorDofIdx)
139 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
140 else
141 specificEnthalpy = Opm::getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
142
143 Evaluation enthalpyRate = density*extQuants.volumeFlux(phaseIdx)*specificEnthalpy;
144 EnergyModule::addToEnthalpyRate(*this, enthalpyRate);
145 }
146 }
147
148 // thermal conduction
149 EnergyModule::addToEnthalpyRate(*this, EnergyModule::thermalConductionRate(extQuants));
150
151#ifndef NDEBUG
152 for (unsigned i = 0; i < numEq; ++i) {
153 Opm::Valgrind::CheckDefined((*this)[i]);
154 }
155#endif
156 }
157
161 template <class Context, class FluidState>
162 void setInFlow(const Context& context,
163 unsigned bfIdx,
164 unsigned timeIdx,
165 const FluidState& fluidState)
166 {
167 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
168
169 // we only allow fluxes in the direction opposite to the outer unit normal
170 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
171 Evaluation& val = this->operator[](eqIdx);
172 val = Toolbox::min(0.0, val);
173 }
174 }
175
179 template <class Context, class FluidState>
180 void setOutFlow(const Context& context,
181 unsigned bfIdx,
182 unsigned timeIdx,
183 const FluidState& fluidState)
184 {
185 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
186
187 // we only allow fluxes in the same direction as the outer unit normal
188 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
189 Evaluation& val = this->operator[](eqIdx);
190 val = Toolbox::max(0.0, val);
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 model which is based on...
Definition: flashboundaryratevector.hh:45
FlashBoundaryRateVector(const FlashBoundaryRateVector &value)=default
FlashBoundaryRateVector(const Evaluation &value)
Definition: flashboundaryratevector.hh:70
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: flashboundaryratevector.hh:197
FlashBoundaryRateVector & operator=(const FlashBoundaryRateVector &value)=default
void setFreeFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: flashboundaryratevector.hh:84
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: flashboundaryratevector.hh:162
FlashBoundaryRateVector()
Definition: flashboundaryratevector.hh:63
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: flashboundaryratevector.hh:180
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