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