ncpboundaryratevector.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_NCP_BOUNDARY_RATE_VECTOR_HH
29#define EWOMS_NCP_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
48template <class TypeTag>
49class NcpBoundaryRateVector : public GetPropType<TypeTag, Properties::RateVector>
50{
57
58 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
59 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
60 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
61 enum { conti0EqIdx = Indices::conti0EqIdx };
62 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
63
65
66 using Toolbox = MathToolbox<Evaluation>;
67
68public:
70
75 NcpBoundaryRateVector(const Evaluation& value)
76 : ParentType(value)
77 {}
78
85
89 template <class Context, class FluidState>
90 void setFreeFlow(const Context& context,
91 unsigned bfIdx,
92 unsigned timeIdx,
93 const FluidState& fluidState)
94 {
95 ExtensiveQuantities extQuants;
96 extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
97 const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
98 const unsigned focusDofIdx = context.focusDofIndex();
99 const unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
100
102 // advective fluxes of all components in all phases
104 (*this) = Evaluation(0.0);
105 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
106 Evaluation density;
107 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
108 if (focusDofIdx == interiorDofIdx) {
109 density = fluidState.density(phaseIdx);
110 }
111 else {
112 density = getValue(fluidState.density(phaseIdx));
113 }
114 }
115 else if (focusDofIdx == interiorDofIdx) {
116 density = insideIntQuants.fluidState().density(phaseIdx);
117 }
118 else {
119 density = getValue(insideIntQuants.fluidState().density(phaseIdx));
120 }
121
122 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
123 Evaluation molarity;
124 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
125 if (focusDofIdx == interiorDofIdx) {
126 molarity = fluidState.molarity(phaseIdx, compIdx);
127 }
128 else {
129 molarity = getValue(fluidState.molarity(phaseIdx, compIdx));
130 }
131 }
132 else if (focusDofIdx == interiorDofIdx) {
133 molarity = insideIntQuants.fluidState().molarity(phaseIdx, compIdx);
134 }
135 else {
136 molarity = getValue(insideIntQuants.fluidState().molarity(phaseIdx, compIdx));
137 }
138
139 // add advective flux of current component in current
140 // phase
141 (*this)[conti0EqIdx + compIdx] += extQuants.volumeFlux(phaseIdx) * molarity;
142 }
143
144 if (enableEnergy) {
145 Evaluation specificEnthalpy;
146 if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
147 if (focusDofIdx == interiorDofIdx) {
148 specificEnthalpy = fluidState.enthalpy(phaseIdx);
149 }
150 else {
151 specificEnthalpy = getValue(fluidState.enthalpy(phaseIdx));
152 }
153 }
154 else if (focusDofIdx == interiorDofIdx) {
155 specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
156 }
157 else {
158 specificEnthalpy = getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
159 }
160
161 const Evaluation enthalpyRate = density * extQuants.volumeFlux(phaseIdx) * specificEnthalpy;
162 EnergyModule::addToEnthalpyRate(*this, enthalpyRate);
163 }
164 }
165
166 // thermal conduction
167 EnergyModule::addToEnthalpyRate(*this, EnergyModule::thermalConductionRate(extQuants));
168
169#ifndef NDEBUG
170 for (unsigned i = 0; i < numEq; ++i) {
171 Valgrind::CheckDefined((*this)[i]);
172 }
173#endif
174 }
175
179 template <class Context, class FluidState>
180 void setInFlow(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 direction opposite to the outer unit normal
188 std::for_each(this->begin(), this->end(),
189 [](auto& val) { val = Toolbox::min(Scalar{0.0}, val); });
190 }
191
195 template <class Context, class FluidState>
196 void setOutFlow(const Context& context,
197 unsigned bfIdx,
198 unsigned timeIdx,
199 const FluidState& fluidState)
200 {
201 this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
202
203 // we only allow fluxes in the same direction as the outer unit normal
204 std::for_each(this->begin(), this->end(),
205 [](auto& val) { val = Toolbox::max(Scalar{0.0}, val); });
206 }
207
212 { (*this) = Evaluation(0.0); }
213};
214
215} // namespace Opm
216
217#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 NCP model.
Definition: ncpboundaryratevector.hh:50
NcpBoundaryRateVector(const Evaluation &value)
Definition: ncpboundaryratevector.hh:75
void setInFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: ncpboundaryratevector.hh:180
void setFreeFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: ncpboundaryratevector.hh:90
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: ncpboundaryratevector.hh:211
NcpBoundaryRateVector & operator=(const NcpBoundaryRateVector &value)=default
void setOutFlow(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: ncpboundaryratevector.hh:196
NcpBoundaryRateVector(const NcpBoundaryRateVector &value)=default
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