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