My Project
Loading...
Searching...
No Matches
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
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 <iostream>
44
45namespace Opm {
46
56template <class TypeTag>
58{
62
69
70 // primary variable indices
71 enum { pressure0Idx = Indices::pressure0Idx };
72 enum { switch0Idx = Indices::switch0Idx };
73
75 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
77
78 using Toolbox = MathToolbox<Evaluation>;
79 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
82
83public:
85 { Valgrind::SetDefined(*this); }
86
92 {
93 Valgrind::SetDefined(*this);
94
95 phasePresence_ = value.phasePresence_;
96 }
97
101 template <class FluidState>
102 void assignMassConservative(const FluidState& fluidState,
103 const MaterialLawParams& matParams,
104 bool isInEquilibrium = false)
105 {
106#ifndef NDEBUG
107 // make sure the temperature is the same in all fluid phases
108 for (unsigned phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
109 assert(std::abs(fluidState.temperature(0) - fluidState.temperature(phaseIdx)) < 1e-30);
110 }
111#endif // NDEBUG
112
113 // for the equilibrium case, we don't need complicated
114 // computations.
115 if (isInEquilibrium) {
116 assignNaive(fluidState);
117 return;
118 }
119
120 // use a flash calculation to calculate a fluid state in
121 // thermodynamic equilibrium
122 typename FluidSystem::template ParameterCache<Scalar> paramCache;
124
125 // use the externally given fluid state as initial value for
126 // the flash calculation
127 fsFlash.assign(fluidState);
128
129 // calculate the phase densities
130 paramCache.updateAll(fsFlash);
131 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
132 Scalar rho = FluidSystem::density(fsFlash, paramCache, phaseIdx);
133 fsFlash.setDensity(phaseIdx, rho);
134 }
135 // calculate the "global molarities"
136 ComponentVector globalMolarities(0.0);
137 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
138 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
140 fsFlash.saturation(phaseIdx) * fsFlash.molarity(phaseIdx, compIdx);
141 }
142 }
143
144 // run the flash calculation
146
147 // use the result to assign the primary variables
149 }
150
155 short phasePresence() const
156 { return phasePresence_; }
157
164 void setPhasePresence(short value)
165 { phasePresence_ = value; }
166
174 void setPhasePresent(unsigned phaseIdx, bool yesno = true)
175 {
176 if (yesno)
177 setPhasePresence(phasePresence_ | (1 << phaseIdx));
178 else
179 setPhasePresence(phasePresence_& ~(1 << phaseIdx));
180 }
181
186 unsigned implicitSaturationIdx() const
187 { return lowestPresentPhaseIdx(); }
188
197 static bool phaseIsPresent(unsigned phaseIdx, short phasePresence)
198 { return phasePresence& (1 << phaseIdx); }
199
206 bool phaseIsPresent(unsigned phaseIdx) const
207 { return phasePresence_& (1 << phaseIdx); }
208
212 unsigned lowestPresentPhaseIdx() const
213 { return static_cast<unsigned>(ffs(phasePresence_) - 1); }
214
218 ThisType& operator=(const Implementation& value)
219 {
221 phasePresence_ = value.phasePresence_;
222
223 return *this;
224 }
225
229 ThisType& operator=(Scalar value)
230 {
232
233 phasePresence_ = 0;
234 return *this;
235 }
236
244 Evaluation explicitSaturationValue(unsigned phaseIdx, unsigned timeIdx) const
245 {
247 // non-present phases have saturation 0
248 return 0.0;
249
250 unsigned varIdx = switch0Idx + phaseIdx - 1;
251 if (std::is_same<Evaluation, Scalar>::value)
252 return (*this)[varIdx]; // finite differences
253 else {
254 // automatic differentiation
255 if (timeIdx != 0)
256 Toolbox::createConstant((*this)[varIdx]);
257 return Toolbox::createVariable((*this)[varIdx], varIdx);
258 }
259 }
260
264 template <class FluidState>
265 void assignNaive(const FluidState& fluidState)
266 {
268
269 // assign the phase temperatures. this is out-sourced to
270 // the energy module
271 EnergyModule::setPriVarTemperatures(*this, fluidState);
272
273 // set the pressure of the first phase
274 (*this)[pressure0Idx] = FsToolbox::value(fluidState.pressure(/*phaseIdx=*/0));
275 Valgrind::CheckDefined((*this)[pressure0Idx]);
276
277 // determine the phase presence.
278 phasePresence_ = 0;
279 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
280 // use a NCP condition to determine if the phase is
281 // present or not
282 Scalar a = 1;
283 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
284 a -= FsToolbox::value(fluidState.moleFraction(phaseIdx, compIdx));
285 }
286 Scalar b = FsToolbox::value(fluidState.saturation(phaseIdx));
287
288 if (b > a)
289 phasePresence_ |= (1 << phaseIdx);
290 }
291
292 // some phase must be present
293 if (phasePresence_ == 0)
294 throw NumericalProblem("Phase state was 0, i.e., no fluid is present");
295
296 // set the primary variables which correspond to mole
297 // fractions of the present phase which has the lowest index.
299 for (unsigned switchIdx = 0; switchIdx < numPhases - 1; ++switchIdx) {
300 unsigned phaseIdx = switchIdx;
301 unsigned compIdx = switchIdx + 1;
303 ++phaseIdx;
304
306 (*this)[switch0Idx + switchIdx] = FsToolbox::value(fluidState.saturation(phaseIdx));
307 Valgrind::CheckDefined((*this)[switch0Idx + switchIdx]);
308 }
309 else {
310 (*this)[switch0Idx + switchIdx] =
311 FsToolbox::value(fluidState.moleFraction(lowestPhaseIdx, compIdx));
312 Valgrind::CheckDefined((*this)[switch0Idx + switchIdx]);
313 }
314 }
315
316 // set the mole fractions in of the remaining components in
317 // the phase with the lowest index
318 for (unsigned compIdx = numPhases - 1; compIdx < numComponents - 1; ++compIdx) {
319 (*this)[switch0Idx + compIdx] =
320 FsToolbox::value(fluidState.moleFraction(lowestPhaseIdx, compIdx + 1));
321 Valgrind::CheckDefined((*this)[switch0Idx + compIdx]);
322 }
323 }
324
328 void print(std::ostream& os = std::cout) const
329 {
330 os << "(p_" << FluidSystem::phaseName(0) << " = "
331 << this->operator[](pressure0Idx);
333 for (unsigned switchIdx = 0; switchIdx < numPhases - 1; ++switchIdx) {
334 unsigned phaseIdx = switchIdx;
335 unsigned compIdx = switchIdx + 1;
337 ++phaseIdx; // skip the saturation of the present
338 // phase with the lowest index
339
341 os << ", S_" << FluidSystem::phaseName(phaseIdx) << " = "
342 << (*this)[switch0Idx + switchIdx];
343 }
344 else {
345 os << ", x_" << FluidSystem::phaseName(lowestPhaseIdx) << "^"
346 << FluidSystem::componentName(compIdx) << " = "
347 << (*this)[switch0Idx + switchIdx];
348 }
349 }
350 for (unsigned compIdx = numPhases - 1; compIdx < numComponents - 1;
351 ++compIdx) {
352 os << ", x_" << FluidSystem::phaseName(lowestPhaseIdx) << "^"
353 << FluidSystem::componentName(compIdx + 1) << " = "
354 << (*this)[switch0Idx + compIdx];
355 }
356 os << ")";
357 os << ", phase presence: " << static_cast<int>(phasePresence_) << std::flush;
358 }
359
360private:
361 short phasePresence_;
362};
363
364} // namespace Opm
365
366#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:58
void print(std::ostream &os=std::cout) const
Prints the names of the primary variables and their values.
Definition pvsprimaryvariables.hh:328
unsigned implicitSaturationIdx() const
Returns the index of the phase with's its saturation is determined by the closure condition of satura...
Definition pvsprimaryvariables.hh:186
void setPhasePresence(short value)
Set which fluid phases are present in a given control volume.
Definition pvsprimaryvariables.hh:164
static bool phaseIsPresent(unsigned phaseIdx, short phasePresence)
Returns true iff a phase is present for a given phase presence.
Definition pvsprimaryvariables.hh:197
void setPhasePresent(unsigned phaseIdx, bool yesno=true)
Set whether a given indivividual phase should be present or not.
Definition pvsprimaryvariables.hh:174
ThisType & operator=(const Implementation &value)
Assignment operator from an other primary variables object.
Definition pvsprimaryvariables.hh:218
void assignNaive(const FluidState &fluidState)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition pvsprimaryvariables.hh:265
PvsPrimaryVariables(const PvsPrimaryVariables &value)
Definition pvsprimaryvariables.hh:91
unsigned lowestPresentPhaseIdx() const
Returns the phase with the lowest index that is present.
Definition pvsprimaryvariables.hh:212
ThisType & operator=(Scalar value)
Assignment operator from a scalar value.
Definition pvsprimaryvariables.hh:229
Evaluation explicitSaturationValue(unsigned phaseIdx, unsigned timeIdx) const
Returns an explcitly stored saturation for a given phase.
Definition pvsprimaryvariables.hh:244
bool phaseIsPresent(unsigned phaseIdx) const
Returns true iff a phase is present for the current phase presence.
Definition pvsprimaryvariables.hh:206
short phasePresence() const
Return the fluid phases which are present in a given control volume.
Definition pvsprimaryvariables.hh:155
void assignMassConservative(const FluidState &fluidState, const MaterialLawParams &matParams, bool isInEquilibrium=false)
< Import base class assignment operators.
Definition pvsprimaryvariables.hh:102
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Represents the primary variables used by the a model.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilboundaryratevector.hh:37
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242
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.