blackoilboundaryratevector.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_BLACK_OIL_BOUNDARY_RATE_VECTOR_HH
29#define EWOMS_BLACK_OIL_BOUNDARY_RATE_VECTOR_HH
30
31#include <opm/material/common/Valgrind.hpp>
32#include <opm/material/constraintsolvers/NcpFlash.hpp>
33
35
36#include <algorithm>
37#include <type_traits>
38
39namespace Opm {
40
46template <class TypeTag>
47class BlackOilBoundaryRateVector : public GetPropType<TypeTag, Properties::RateVector>
48{
57
58 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
59 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
60 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
61 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
62 enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
63 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
64 enum { conti0EqIdx = Indices::conti0EqIdx };
65 enum { contiEnergyEqIdx = Indices::contiEnergyEqIdx };
66 enum { enableFoam = getPropValue<TypeTag, Properties::EnableFoam>() };
67 enum { enableMICP = getPropValue<TypeTag, Properties::EnableMICP>() };
68
69 static constexpr bool blackoilConserveSurfaceVolume =
70 getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
71
73
74public:
79
83 BlackOilBoundaryRateVector(Scalar value) : ParentType(value)
84 {}
85
91
95 template <class Context, class FluidState>
96 void setFreeFlow(const Context& context,
97 unsigned bfIdx,
98 unsigned timeIdx,
99 const FluidState& fluidState)
100 {
101 ExtensiveQuantities extQuants;
102 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
103 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
104 const unsigned focusDofIdx = context.focusDofIndex();
105 const unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
106
108 // advective fluxes of all components in all phases
110 (*this) = 0.0;
111 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
112 if (!FluidSystem::phaseIsActive(phaseIdx)) {
113 continue;
114 }
115 const auto& pBoundary = fluidState.pressure(phaseIdx);
116 const Evaluation& pInside = insideIntQuants.fluidState().pressure(phaseIdx);
117
118 RateVector tmp;
119
120 // mass conservation
121 if (pBoundary < pInside) {
122 // outflux
123 LocalResidual::template evalPhaseFluxes_<Evaluation>(tmp,
124 phaseIdx,
125 insideIntQuants.pvtRegionIndex(),
126 extQuants,
127 insideIntQuants.fluidState());
128 }
129 else if (pBoundary > pInside) {
130 using RhsEval = std::conditional_t<std::is_same_v<typename FluidState::Scalar, Evaluation>,
131 Evaluation, Scalar>;
132 // influx
133 LocalResidual::template evalPhaseFluxes_<RhsEval>(tmp,
134 phaseIdx,
135 insideIntQuants.pvtRegionIndex(),
136 extQuants,
137 fluidState);
138 }
139
140 for (unsigned i = 0; i < tmp.size(); ++i) {
141 (*this)[i] += tmp[i];
142 }
143
144 // energy conservation
145 if constexpr (enableEnergy) {
146 Evaluation density;
147 Evaluation specificEnthalpy;
148 if (pBoundary > pInside) {
149 if (focusDofIdx == interiorDofIdx) {
150 density = fluidState.density(phaseIdx);
151 specificEnthalpy = fluidState.enthalpy(phaseIdx);
152 }
153 else {
154 density = getValue(fluidState.density(phaseIdx));
155 specificEnthalpy = getValue(fluidState.enthalpy(phaseIdx));
156 }
157 }
158 else if (focusDofIdx == interiorDofIdx) {
159 density = insideIntQuants.fluidState().density(phaseIdx);
160 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
161 }
162 else {
163 density = getValue(insideIntQuants.fluidState().density(phaseIdx));
164 specificEnthalpy = getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
165 }
166
167 const Evaluation enthalpyRate = density * extQuants.volumeFlux(phaseIdx) * specificEnthalpy;
168 EnergyModule::addToEnthalpyRate(*this, enthalpyRate *
169 getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>());
170 }
171 }
172
173 if constexpr (enableSolvent) {
174 (*this)[Indices::contiSolventEqIdx] = extQuants.solventVolumeFlux();
175 if (blackoilConserveSurfaceVolume) {
176 (*this)[Indices::contiSolventEqIdx] *= insideIntQuants.solventInverseFormationVolumeFactor();
177 }
178 else {
179 (*this)[Indices::contiSolventEqIdx] *= insideIntQuants.solventDensity();
180 }
181 }
182
183 if constexpr (enablePolymer) {
184 (*this)[Indices::contiPolymerEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
185 insideIntQuants.polymerConcentration();
186 }
187
188 if constexpr (enableMICP) {
189 (*this)[Indices::contiMicrobialEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
190 insideIntQuants.microbialConcentration();
191 (*this)[Indices::contiOxygenEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
192 insideIntQuants.oxygenConcentration();
193 (*this)[Indices::contiUreaEqIdx] = extQuants.volumeFlux(FluidSystem::waterPhaseIdx) *
194 insideIntQuants.ureaConcentration();
195 // since the urea concentration can be much larger than 1, then we apply a scaling factor
196 (*this)[Indices::contiUreaEqIdx] *= getPropValue<TypeTag, Properties::BlackOilUreaScalingFactor>();
197 }
198
199 // make sure that the right mass conservation quantities are used
200 LocalResidual::adaptMassConservationQuantities_(*this, insideIntQuants.pvtRegionIndex());
201
202 // heat conduction
203 if constexpr (enableEnergy) {
204 EnergyModule::addToEnthalpyRate(*this, extQuants.energyFlux() *
205 getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>());
206 }
207
208#ifndef NDEBUG
209 for (unsigned i = 0; i < numEq; ++i) {
210 Valgrind::CheckDefined((*this)[i]);
211 }
212 Valgrind::CheckDefined(*this);
213#endif
214 }
215
219 template <class Context, class FluidState>
220 void setInFlow(const Context& context,
221 unsigned bfIdx,
222 unsigned timeIdx,
223 const FluidState& fluidState)
224 {
225 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
226
227 // we only allow fluxes in the direction opposite to the outer
228 // unit normal
229 std::for_each(this->begin(), this->end(),
230 [](auto& val) { val = std::min(Scalar(0), val); });
231 }
232
236 template <class Context, class FluidState>
237 void setOutFlow(const Context& context,
238 unsigned bfIdx,
239 unsigned timeIdx,
240 const FluidState& fluidState)
241 {
242 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
243
244 // we only allow fluxes in the same direction as the outer
245 // unit normal
246 std::for_each(this->begin(), this->end(),
247 [](auto& val) { val = std::max(Scalar(0), val); });
248 }
249
254 { (*this) = Scalar(0); }
255
263 template <class Context, class FluidState>
264 void setThermalFlow([[maybe_unused]] const Context& context,
265 [[maybe_unused]] unsigned bfIdx,
266 [[maybe_unused]] unsigned timeIdx,
267 [[maybe_unused]] const FluidState& boundaryFluidState)
268 {
269 // set the mass no-flow condition
270 setNoFlow();
271
272 // if we do not conserve energy there is nothing we should do in addition
273 if constexpr (enableEnergy) {
274 ExtensiveQuantities extQuants;
275 extQuants.updateBoundary(context, bfIdx, timeIdx, boundaryFluidState);
276
277 (*this)[contiEnergyEqIdx] += extQuants.energyFlux();
278
279#ifndef NDEBUG
280 for (unsigned i = 0; i < numEq; ++i) {
281 Valgrind::CheckDefined((*this)[i]);
282 }
283 Valgrind::CheckDefined(*this);
284#endif
285 }
286 }
287};
288
289} // namespace Opm
290
291#endif
Contains the classes required to extend the black-oil model by energy.
Implements a boundary vector for the fully implicit black-oil model.
Definition: blackoilboundaryratevector.hh:48
BlackOilBoundaryRateVector(Scalar value)
Definition: blackoilboundaryratevector.hh:83
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: blackoilboundaryratevector.hh:237
BlackOilBoundaryRateVector & operator=(const BlackOilBoundaryRateVector &value)=default
BlackOilBoundaryRateVector(const BlackOilBoundaryRateVector &value)=default
Copy constructor.
void setFreeFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: blackoilboundaryratevector.hh:96
BlackOilBoundaryRateVector()=default
Default constructor.
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: blackoilboundaryratevector.hh:253
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: blackoilboundaryratevector.hh:220
void setThermalFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &boundaryFluidState)
an energy flux that corresponds to the thermal conduction from
Definition: blackoilboundaryratevector.hh:264
Contains the high level supplements required to extend the black oil model by energy.
Definition: blackoilenergymodules.hh:58
static void addToEnthalpyRate(RateVector &flux, const Evaluation &hRate)
Definition: blackoilenergymodules.hh:251
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