pvsprimaryvariables.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_PRIMARY_VARIABLES_HH
29#define EWOMS_PVS_PRIMARY_VARIABLES_HH
30
31#include <dune/common/fvector.hh>
32
33#include <opm/common/Exceptions.hpp>
34
35#include <opm/material/common/MathToolbox.hpp>
36#include <opm/material/common/Valgrind.hpp>
37#include <opm/material/constraintsolvers/NcpFlash.hpp>
38#include <opm/material/fluidstates/CompositionalFluidState.hpp>
39
43
44#include <cassert>
45#include <cmath>
46#include <ostream>
47#include <type_traits>
48
49namespace Opm {
50
60template <class TypeTag>
62{
63 using ParentType = FvBasePrimaryVariables<TypeTag>;
66
73
74 // primary variable indices
75 enum { pressure0Idx = Indices::pressure0Idx };
76 enum { switch0Idx = Indices::switch0Idx };
77
78 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
79 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
80 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
81
82 using Toolbox = MathToolbox<Evaluation>;
83 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
85 using NcpFlash = ::Opm::NcpFlash<Scalar, FluidSystem>;
86
87public:
88 PvsPrimaryVariables() : ParentType()
89 { Valgrind::SetDefined(*this); }
90
95 PvsPrimaryVariables(const PvsPrimaryVariables& value) : ParentType(value)
96 {
97 Valgrind::SetDefined(*this);
98
99 phasePresence_ = value.phasePresence_;
100 }
101
105 template <class FluidState>
106 void assignMassConservative(const FluidState& fluidState,
107 const MaterialLawParams& matParams,
108 bool isInEquilibrium = false)
109 {
110#ifndef NDEBUG
111 // make sure the temperature is the same in all fluid phases
112 for (unsigned phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
113 assert(std::abs(fluidState.temperature(0) - fluidState.temperature(phaseIdx)) < 1e-30);
114 }
115#endif // NDEBUG
116
117 // for the equilibrium case, we don't need complicated
118 // computations.
119 if (isInEquilibrium) {
120 assignNaive(fluidState);
121 return;
122 }
123
124 // use a flash calculation to calculate a fluid state in
125 // thermodynamic equilibrium
126 typename FluidSystem::template ParameterCache<Scalar> paramCache;
127 CompositionalFluidState<Scalar, FluidSystem> fsFlash;
128
129 // use the externally given fluid state as initial value for
130 // the flash calculation
131 fsFlash.assign(fluidState);
132
133 // calculate the phase densities
134 paramCache.updateAll(fsFlash);
135 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
136 const Scalar rho = FluidSystem::density(fsFlash, paramCache, phaseIdx);
137 fsFlash.setDensity(phaseIdx, rho);
138 }
139 // calculate the "global molarities"
140 ComponentVector globalMolarities(0.0);
141 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
142 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
143 globalMolarities[compIdx] +=
144 fsFlash.saturation(phaseIdx) * fsFlash.molarity(phaseIdx, compIdx);
145 }
146 }
147
148 // run the flash calculation
149 NcpFlash::template solve<MaterialLaw>(fsFlash, matParams, paramCache, globalMolarities);
150
151 // use the result to assign the primary variables
152 assignNaive(fsFlash);
153 }
154
159 short phasePresence() const
160 { return phasePresence_; }
161
168 void setPhasePresence(short value)
169 { phasePresence_ = value; }
170
178 void setPhasePresent(unsigned phaseIdx, bool yesno)
179 {
180 if (yesno) {
181 setPhasePresence(phasePresence_ | (1 << phaseIdx));
182 }
183 else {
184 setPhasePresence(phasePresence_ & ~(1 << phaseIdx));
185 }
186 }
187
192 unsigned implicitSaturationIdx() const
193 { return lowestPresentPhaseIdx(); }
194
203 static bool phaseIsPresent(unsigned phaseIdx, short phasePresence)
204 { return phasePresence & (1 << phaseIdx); }
205
212 bool phaseIsPresent(unsigned phaseIdx) const
213 { return phasePresence_ & (1 << phaseIdx); }
214
218 unsigned lowestPresentPhaseIdx() const
219 { return static_cast<unsigned>(ffs(phasePresence_) - 1); }
220
224 ThisType& operator=(const Implementation& value)
225 {
227 phasePresence_ = value.phasePresence_;
228
229 return *this;
230 }
231
235 ThisType& operator=(Scalar value)
236 {
238
239 phasePresence_ = 0;
240 return *this;
241 }
242
250 Evaluation explicitSaturationValue(unsigned phaseIdx, unsigned timeIdx) const
251 {
252 if (!phaseIsPresent(phaseIdx) || phaseIdx == lowestPresentPhaseIdx()) {
253 // non-present phases have saturation 0
254 return 0.0;
255 }
256
257 const unsigned varIdx = switch0Idx + phaseIdx - 1;
258 if constexpr (std::is_same_v<Evaluation, Scalar>) {
259 return (*this)[varIdx]; // finite differences
260 }
261 else {
262 // automatic differentiation
263 if (timeIdx != 0) {
264 Toolbox::createConstant((*this)[varIdx]);
265 }
266 return Toolbox::createVariable((*this)[varIdx], varIdx);
267 }
268 }
269
273 template <class FluidState>
274 void assignNaive(const FluidState& fluidState)
275 {
276 using FsToolbox = MathToolbox<typename FluidState::Scalar>;
277
278 // assign the phase temperatures. this is out-sourced to
279 // the energy module
280 EnergyModule::setPriVarTemperatures(*this, fluidState);
281
282 // set the pressure of the first phase
283 (*this)[pressure0Idx] = FsToolbox::value(fluidState.pressure(/*phaseIdx=*/0));
284 Valgrind::CheckDefined((*this)[pressure0Idx]);
285
286 // determine the phase presence.
287 phasePresence_ = 0;
288 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
289 // use a NCP condition to determine if the phase is
290 // present or not
291 Scalar a = 1;
292 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
293 a -= FsToolbox::value(fluidState.moleFraction(phaseIdx, compIdx));
294 }
295 const Scalar b = FsToolbox::value(fluidState.saturation(phaseIdx));
296
297 if (b > a) {
298 phasePresence_ |= (1 << phaseIdx);
299 }
300 }
301
302 // some phase must be present
303 if (phasePresence_ == 0) {
304 throw NumericalProblem("Phase state was 0, i.e., no fluid is present");
305 }
306
307 // set the primary variables which correspond to mole
308 // fractions of the present phase which has the lowest index.
309 const unsigned lowestPhaseIdx = lowestPresentPhaseIdx();
310 for (unsigned switchIdx = 0; switchIdx < numPhases - 1; ++switchIdx) {
311 unsigned phaseIdx = switchIdx;
312 const unsigned compIdx = switchIdx + 1;
313 if (switchIdx >= lowestPhaseIdx) {
314 ++phaseIdx;
315 }
316
317 if (phaseIsPresent(phaseIdx)) {
318 (*this)[switch0Idx + switchIdx] = FsToolbox::value(fluidState.saturation(phaseIdx));
319 Valgrind::CheckDefined((*this)[switch0Idx + switchIdx]);
320 }
321 else {
322 (*this)[switch0Idx + switchIdx] =
323 FsToolbox::value(fluidState.moleFraction(lowestPhaseIdx, compIdx));
324 Valgrind::CheckDefined((*this)[switch0Idx + switchIdx]);
325 }
326 }
327
328 // set the mole fractions in of the remaining components in
329 // the phase with the lowest index
330 for (unsigned compIdx = numPhases - 1; compIdx < numComponents - 1; ++compIdx) {
331 (*this)[switch0Idx + compIdx] =
332 FsToolbox::value(fluidState.moleFraction(lowestPhaseIdx, compIdx + 1));
333 Valgrind::CheckDefined((*this)[switch0Idx + compIdx]);
334 }
335 }
336
340 void print(std::ostream& os) const
341 {
342 os << "(p_" << FluidSystem::phaseName(0) << " = "
343 << (*this)[pressure0Idx];
344 const unsigned lowestPhaseIdx = lowestPresentPhaseIdx();
345 for (unsigned switchIdx = 0; switchIdx < numPhases - 1; ++switchIdx) {
346 unsigned phaseIdx = switchIdx;
347 const unsigned compIdx = switchIdx + 1;
348 if (phaseIdx >= lowestPhaseIdx) {
349 ++phaseIdx; // skip the saturation of the present
350 // phase with the lowest index
351 }
352
353 if (phaseIsPresent(phaseIdx)) {
354 os << ", S_" << FluidSystem::phaseName(phaseIdx) << " = "
355 << (*this)[switch0Idx + switchIdx];
356 }
357 else {
358 os << ", x_" << FluidSystem::phaseName(lowestPhaseIdx) << "^"
359 << FluidSystem::componentName(compIdx) << " = "
360 << (*this)[switch0Idx + switchIdx];
361 }
362 }
363 for (unsigned compIdx = numPhases - 1; compIdx < numComponents - 1; ++compIdx) {
364 os << ", x_" << FluidSystem::phaseName(lowestPhaseIdx) << "^"
365 << FluidSystem::componentName(compIdx + 1) << " = "
366 << (*this)[switch0Idx + compIdx];
367 }
368 os << ")";
369 os << ", phase presence: " << static_cast<int>(phasePresence_) << std::flush;
370 }
371
372private:
373 short phasePresence_{};
374};
375
376} // namespace Opm
377
378#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:54
Represents the primary variables used by the a model.
Definition: fvbaseprimaryvariables.hh:53
FvBasePrimaryVariables & operator=(const FvBasePrimaryVariables &value)=default
Assignment from another primary variables object.
Represents the primary variables used in the primary variable switching compositional model.
Definition: pvsprimaryvariables.hh:62
unsigned implicitSaturationIdx() const
Returns the index of the phase with's its saturation is determined by the closure condition of satura...
Definition: pvsprimaryvariables.hh:192
void setPhasePresence(short value)
Set which fluid phases are present in a given control volume.
Definition: pvsprimaryvariables.hh:168
static bool phaseIsPresent(unsigned phaseIdx, short phasePresence)
Returns true iff a phase is present for a given phase presence.
Definition: pvsprimaryvariables.hh:203
ThisType & operator=(const Implementation &value)
Assignment operator from an other primary variables object.
Definition: pvsprimaryvariables.hh:224
void assignNaive(const FluidState &fluidState)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: pvsprimaryvariables.hh:274
PvsPrimaryVariables(const PvsPrimaryVariables &value)
Definition: pvsprimaryvariables.hh:95
unsigned lowestPresentPhaseIdx() const
Returns the phase with the lowest index that is present.
Definition: pvsprimaryvariables.hh:218
ThisType & operator=(Scalar value)
Assignment operator from a scalar value.
Definition: pvsprimaryvariables.hh:235
void print(std::ostream &os) const
Prints the names of the primary variables and their values.
Definition: pvsprimaryvariables.hh:340
Evaluation explicitSaturationValue(unsigned phaseIdx, unsigned timeIdx) const
Returns an explcitly stored saturation for a given phase.
Definition: pvsprimaryvariables.hh:250
bool phaseIsPresent(unsigned phaseIdx) const
Returns true iff a phase is present for the current phase presence.
Definition: pvsprimaryvariables.hh:212
short phasePresence() const
Return the fluid phases which are present in a given control volume.
Definition: pvsprimaryvariables.hh:159
void assignMassConservative(const FluidState &fluidState, const MaterialLawParams &matParams, bool isInEquilibrium=false)
< Import base class assignment operators.
Definition: pvsprimaryvariables.hh:106
void setPhasePresent(unsigned phaseIdx, bool yesno)
Set whether a given indivividual phase should be present or not.
Definition: pvsprimaryvariables.hh:178
PvsPrimaryVariables()
Definition: pvsprimaryvariables.hh:88
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
Declares the properties required for the compositional multi-phase primary variable switching model.