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 "pvsindices.hh"
32#include "pvsproperties.hh"
33
34#include <opm/common/Exceptions.hpp>
35
38
39#include <opm/material/constraintsolvers/NcpFlash.hpp>
40#include <opm/material/fluidstates/CompositionalFluidState.hpp>
41#include <opm/material/common/Valgrind.hpp>
42
43#include <dune/common/fvector.hh>
44
45#include <iostream>
46
47namespace Opm {
48
58template <class TypeTag>
60{
61 using ParentType = FvBasePrimaryVariables<TypeTag>;
64
71
72 // primary variable indices
73 enum { pressure0Idx = Indices::pressure0Idx };
74 enum { switch0Idx = Indices::switch0Idx };
75
76 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
77 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
78 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
79
80 using Toolbox = MathToolbox<Evaluation>;
81 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
83 using NcpFlash = Opm::NcpFlash<Scalar, FluidSystem>;
84
85public:
86 PvsPrimaryVariables() : ParentType()
87 { Valgrind::SetDefined(*this); }
88
92 explicit PvsPrimaryVariables(Scalar value) : ParentType(value)
93 {
94 Valgrind::CheckDefined(value);
95 Valgrind::SetDefined(*this);
96
97 phasePresence_ = 0;
98 }
99
104 PvsPrimaryVariables(const PvsPrimaryVariables& value) : ParentType(value)
105 {
106 Valgrind::SetDefined(*this);
107
108 phasePresence_ = value.phasePresence_;
109 }
110
114 template <class FluidState>
115 void assignMassConservative(const FluidState& fluidState,
116 const MaterialLawParams& matParams,
117 bool isInEquilibrium = false)
118 {
119#ifndef NDEBUG
120 // make sure the temperature is the same in all fluid phases
121 for (unsigned phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
122 assert(std::abs(fluidState.temperature(0) - fluidState.temperature(phaseIdx)) < 1e-30);
123 }
124#endif // NDEBUG
125
126 // for the equilibrium case, we don't need complicated
127 // computations.
128 if (isInEquilibrium) {
129 assignNaive(fluidState);
130 return;
131 }
132
133 // use a flash calculation to calculate a fluid state in
134 // thermodynamic equilibrium
135 typename FluidSystem::template ParameterCache<Scalar> paramCache;
136 CompositionalFluidState<Scalar, FluidSystem> fsFlash;
137
138 // use the externally given fluid state as initial value for
139 // the flash calculation
140 fsFlash.assign(fluidState);
141
142 // calculate the phase densities
143 paramCache.updateAll(fsFlash);
144 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
145 Scalar rho = FluidSystem::density(fsFlash, paramCache, phaseIdx);
146 fsFlash.setDensity(phaseIdx, rho);
147 }
148 // calculate the "global molarities"
149 ComponentVector globalMolarities(0.0);
150 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
151 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
152 globalMolarities[compIdx] +=
153 fsFlash.saturation(phaseIdx) * fsFlash.molarity(phaseIdx, compIdx);
154 }
155 }
156
157 // run the flash calculation
158 NcpFlash::template solve<MaterialLaw>(fsFlash, matParams, paramCache, globalMolarities);
159
160 // use the result to assign the primary variables
161 assignNaive(fsFlash);
162 }
163
168 short phasePresence() const
169 { return phasePresence_; }
170
177 void setPhasePresence(short value)
178 { phasePresence_ = value; }
179
187 void setPhasePresent(unsigned phaseIdx, bool yesno = true)
188 {
189 if (yesno)
190 setPhasePresence(phasePresence_ | (1 << phaseIdx));
191 else
192 setPhasePresence(phasePresence_& ~(1 << phaseIdx));
193 }
194
199 unsigned implicitSaturationIdx() const
200 { return lowestPresentPhaseIdx(); }
201
210 static bool phaseIsPresent(unsigned phaseIdx, short phasePresence)
211 { return phasePresence& (1 << phaseIdx); }
212
219 bool phaseIsPresent(unsigned phaseIdx) const
220 { return phasePresence_& (1 << phaseIdx); }
221
225 unsigned lowestPresentPhaseIdx() const
226 { return static_cast<unsigned>(ffs(phasePresence_) - 1); }
227
231 ThisType& operator=(const Implementation& value)
232 {
234 phasePresence_ = value.phasePresence_;
235
236 return *this;
237 }
238
242 ThisType& operator=(Scalar value)
243 {
245
246 phasePresence_ = 0;
247 return *this;
248 }
249
257 Evaluation explicitSaturationValue(unsigned phaseIdx, unsigned timeIdx) const
258 {
259 if (!phaseIsPresent(phaseIdx) || phaseIdx == lowestPresentPhaseIdx())
260 // non-present phases have saturation 0
261 return 0.0;
262
263 unsigned varIdx = switch0Idx + phaseIdx - 1;
264 if (std::is_same<Evaluation, Scalar>::value)
265 return (*this)[varIdx]; // finite differences
266 else {
267 // automatic differentiation
268 if (timeIdx != 0)
269 Toolbox::createConstant((*this)[varIdx]);
270 return Toolbox::createVariable((*this)[varIdx], varIdx);
271 }
272 }
273
277 template <class FluidState>
278 void assignNaive(const FluidState& fluidState)
279 {
280 using FsToolbox = MathToolbox<typename FluidState::Scalar>;
281
282 // assign the phase temperatures. this is out-sourced to
283 // the energy module
284 EnergyModule::setPriVarTemperatures(*this, fluidState);
285
286 // set the pressure of the first phase
287 (*this)[pressure0Idx] = FsToolbox::value(fluidState.pressure(/*phaseIdx=*/0));
288 Valgrind::CheckDefined((*this)[pressure0Idx]);
289
290 // determine the phase presence.
291 phasePresence_ = 0;
292 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
293 // use a NCP condition to determine if the phase is
294 // present or not
295 Scalar a = 1;
296 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
297 a -= FsToolbox::value(fluidState.moleFraction(phaseIdx, compIdx));
298 }
299 Scalar b = FsToolbox::value(fluidState.saturation(phaseIdx));
300
301 if (b > a)
302 phasePresence_ |= (1 << phaseIdx);
303 }
304
305 // some phase must be present
306 if (phasePresence_ == 0)
307 throw NumericalProblem("Phase state was 0, i.e., no fluid is present");
308
309 // set the primary variables which correspond to mole
310 // fractions of the present phase which has the lowest index.
311 unsigned lowestPhaseIdx = lowestPresentPhaseIdx();
312 for (unsigned switchIdx = 0; switchIdx < numPhases - 1; ++switchIdx) {
313 unsigned phaseIdx = switchIdx;
314 unsigned compIdx = switchIdx + 1;
315 if (switchIdx >= lowestPhaseIdx)
316 ++phaseIdx;
317
318 if (phaseIsPresent(phaseIdx)) {
319 (*this)[switch0Idx + switchIdx] = FsToolbox::value(fluidState.saturation(phaseIdx));
320 Valgrind::CheckDefined((*this)[switch0Idx + switchIdx]);
321 }
322 else {
323 (*this)[switch0Idx + switchIdx] =
324 FsToolbox::value(fluidState.moleFraction(lowestPhaseIdx, compIdx));
325 Valgrind::CheckDefined((*this)[switch0Idx + switchIdx]);
326 }
327 }
328
329 // set the mole fractions in of the remaining components in
330 // the phase with the lowest index
331 for (unsigned compIdx = numPhases - 1; compIdx < numComponents - 1; ++compIdx) {
332 (*this)[switch0Idx + compIdx] =
333 FsToolbox::value(fluidState.moleFraction(lowestPhaseIdx, compIdx + 1));
334 Valgrind::CheckDefined((*this)[switch0Idx + compIdx]);
335 }
336 }
337
341 void print(std::ostream& os = std::cout) const
342 {
343 os << "(p_" << FluidSystem::phaseName(0) << " = "
344 << this->operator[](pressure0Idx);
345 unsigned lowestPhaseIdx = lowestPresentPhaseIdx();
346 for (unsigned switchIdx = 0; switchIdx < numPhases - 1; ++switchIdx) {
347 unsigned phaseIdx = switchIdx;
348 unsigned compIdx = switchIdx + 1;
349 if (phaseIdx >= lowestPhaseIdx)
350 ++phaseIdx; // skip the saturation of the present
351 // phase with the lowest index
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;
364 ++compIdx) {
365 os << ", x_" << FluidSystem::phaseName(lowestPhaseIdx) << "^"
366 << FluidSystem::componentName(compIdx + 1) << " = "
367 << (*this)[switch0Idx + compIdx];
368 }
369 os << ")";
370 os << ", phase presence: " << static_cast<int>(phasePresence_) << std::flush;
371 }
372
373private:
374 short phasePresence_;
375};
376
377} // namespace Opm
378
379#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:50
Represents the primary variables used by the a model.
Definition: fvbaseprimaryvariables.hh:52
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:60
void print(std::ostream &os=std::cout) const
Prints the names of the primary variables and their values.
Definition: pvsprimaryvariables.hh:341
unsigned implicitSaturationIdx() const
Returns the index of the phase with's its saturation is determined by the closure condition of satura...
Definition: pvsprimaryvariables.hh:199
void setPhasePresence(short value)
Set which fluid phases are present in a given control volume.
Definition: pvsprimaryvariables.hh:177
static bool phaseIsPresent(unsigned phaseIdx, short phasePresence)
Returns true iff a phase is present for a given phase presence.
Definition: pvsprimaryvariables.hh:210
void setPhasePresent(unsigned phaseIdx, bool yesno=true)
Set whether a given indivividual phase should be present or not.
Definition: pvsprimaryvariables.hh:187
ThisType & operator=(const Implementation &value)
Assignment operator from an other primary variables object.
Definition: pvsprimaryvariables.hh:231
void assignNaive(const FluidState &fluidState)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: pvsprimaryvariables.hh:278
PvsPrimaryVariables(const PvsPrimaryVariables &value)
Definition: pvsprimaryvariables.hh:104
PvsPrimaryVariables(Scalar value)
Constructor with assignment from scalar.
Definition: pvsprimaryvariables.hh:92
unsigned lowestPresentPhaseIdx() const
Returns the phase with the lowest index that is present.
Definition: pvsprimaryvariables.hh:225
ThisType & operator=(Scalar value)
Assignment operator from a scalar value.
Definition: pvsprimaryvariables.hh:242
Evaluation explicitSaturationValue(unsigned phaseIdx, unsigned timeIdx) const
Returns an explcitly stored saturation for a given phase.
Definition: pvsprimaryvariables.hh:257
bool phaseIsPresent(unsigned phaseIdx) const
Returns true iff a phase is present for the current phase presence.
Definition: pvsprimaryvariables.hh:219
short phasePresence() const
Return the fluid phases which are present in a given control volume.
Definition: pvsprimaryvariables.hh:168
void assignMassConservative(const FluidState &fluidState, const MaterialLawParams &matParams, bool isInEquilibrium=false)
Set the primary variables from an arbitrary fluid state in a mass conservative way.
Definition: pvsprimaryvariables.hh:115
PvsPrimaryVariables()
Definition: pvsprimaryvariables.hh:86
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Definition: blackoilboundaryratevector.hh:37
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:235
Declares the properties required for the compositional multi-phase primary variable switching model.