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