2015-06-18 06:47:07 -05:00
|
|
|
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
|
|
|
|
// vi: set et ts=4 sw=4 sts=4:
|
2013-12-02 09:31:53 -06:00
|
|
|
/*
|
|
|
|
This file is part of the Open Porous Media project (OPM).
|
|
|
|
|
|
|
|
OPM is free software: you can redistribute it and/or modify
|
|
|
|
it under the terms of the GNU General Public License as published by
|
|
|
|
the Free Software Foundation, either version 2 of the License, or
|
|
|
|
(at your option) any later version.
|
|
|
|
|
|
|
|
OPM is distributed in the hope that it will be useful,
|
|
|
|
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
GNU General Public License for more details.
|
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
2016-03-14 08:11:22 -05:00
|
|
|
|
|
|
|
Consult the COPYING file in the top-level source directory of this
|
|
|
|
module for the precise wording of the license and the list of
|
|
|
|
copyright holders.
|
2013-12-02 09:31:53 -06:00
|
|
|
*/
|
2013-05-30 11:44:10 -05:00
|
|
|
/*!
|
|
|
|
* \file
|
|
|
|
*
|
|
|
|
* \brief This is a program to test the flash calculation which uses
|
2013-12-03 04:56:45 -06:00
|
|
|
* non-linear complementarity problems (NCP)
|
2013-05-30 11:44:10 -05:00
|
|
|
*
|
|
|
|
* A flash calculation determines the pressures, saturations and
|
|
|
|
* composition of all phases given the total mass (or, as in this case
|
|
|
|
* the total number of moles) in a given amount of pore space.
|
|
|
|
*/
|
|
|
|
#include "config.h"
|
|
|
|
|
2016-06-02 11:20:16 -05:00
|
|
|
#include <opm/material/densead/Evaluation.hpp>
|
2013-09-10 07:55:59 -05:00
|
|
|
#include <opm/material/constraintsolvers/MiscibleMultiPhaseComposition.hpp>
|
|
|
|
#include <opm/material/constraintsolvers/ComputeFromReferencePhase.hpp>
|
|
|
|
#include <opm/material/constraintsolvers/ImmiscibleFlash.hpp>
|
2013-05-30 11:44:10 -05:00
|
|
|
|
2013-09-10 07:55:59 -05:00
|
|
|
#include <opm/material/fluidstates/ImmiscibleFluidState.hpp>
|
2013-05-30 11:44:10 -05:00
|
|
|
|
2013-09-10 07:55:59 -05:00
|
|
|
#include <opm/material/fluidsystems/H2ON2FluidSystem.hpp>
|
2013-05-30 11:44:10 -05:00
|
|
|
|
2013-11-06 07:32:22 -06:00
|
|
|
#include <opm/material/fluidmatrixinteractions/LinearMaterial.hpp>
|
|
|
|
#include <opm/material/fluidmatrixinteractions/RegularizedBrooksCorey.hpp>
|
|
|
|
#include <opm/material/fluidmatrixinteractions/EffToAbsLaw.hpp>
|
|
|
|
#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
|
2013-05-30 11:44:10 -05:00
|
|
|
|
2016-04-17 04:28:06 -05:00
|
|
|
#include <opm/common/utility/platform_dependent/disable_warnings.h>
|
|
|
|
#include <dune/common/parallel/mpihelper.hh>
|
|
|
|
#include <opm/common/utility/platform_dependent/reenable_warnings.h>
|
|
|
|
|
2013-05-30 11:44:10 -05:00
|
|
|
template <class Scalar, class FluidState>
|
2017-01-18 10:04:59 -06:00
|
|
|
void checkSame(const FluidState& fsRef, const FluidState& fsFlash)
|
2013-05-30 11:44:10 -05:00
|
|
|
{
|
|
|
|
enum { numPhases = FluidState::numPhases };
|
|
|
|
enum { numComponents = FluidState::numComponents };
|
|
|
|
|
2016-04-17 04:28:21 -05:00
|
|
|
Scalar tol = std::max(std::numeric_limits<Scalar>::epsilon()*1e4, 1e-6);
|
|
|
|
|
2015-09-22 04:20:55 -05:00
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
|
2013-05-30 11:44:10 -05:00
|
|
|
Scalar error;
|
|
|
|
|
|
|
|
// check the pressures
|
|
|
|
error = 1 - fsRef.pressure(phaseIdx)/fsFlash.pressure(phaseIdx);
|
2016-04-17 04:28:21 -05:00
|
|
|
if (std::abs(error) > tol) {
|
|
|
|
OPM_THROW(std::runtime_error,
|
|
|
|
"pressure error phase " << phaseIdx << " is incorrect: "
|
2013-05-30 11:44:10 -05:00
|
|
|
<< fsFlash.pressure(phaseIdx) << " flash vs "
|
|
|
|
<< fsRef.pressure(phaseIdx) << " reference"
|
2016-04-17 04:28:21 -05:00
|
|
|
<< " error=" << error);
|
2013-05-30 11:44:10 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
// check the saturations
|
|
|
|
error = fsRef.saturation(phaseIdx) - fsFlash.saturation(phaseIdx);
|
2016-04-17 04:28:21 -05:00
|
|
|
if (std::abs(error) > tol)
|
|
|
|
OPM_THROW(std::runtime_error,
|
|
|
|
"saturation error phase " << phaseIdx << " is incorrect: "
|
2013-05-30 11:44:10 -05:00
|
|
|
<< fsFlash.saturation(phaseIdx) << " flash vs "
|
|
|
|
<< fsRef.saturation(phaseIdx) << " reference"
|
2016-04-17 04:28:21 -05:00
|
|
|
<< " error=" << error);
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
// check the compositions
|
2015-09-22 04:20:55 -05:00
|
|
|
for (unsigned compIdx = 0; compIdx < numComponents; ++ compIdx) {
|
2013-05-30 11:44:10 -05:00
|
|
|
error = fsRef.moleFraction(phaseIdx, compIdx) - fsFlash.moleFraction(phaseIdx, compIdx);
|
2016-04-17 04:28:21 -05:00
|
|
|
if (std::abs(error) > tol)
|
|
|
|
OPM_THROW(std::runtime_error,
|
|
|
|
"composition error phase " << phaseIdx << ", component " << compIdx
|
|
|
|
<< " is incorrect: "
|
2013-05-30 11:44:10 -05:00
|
|
|
<< fsFlash.moleFraction(phaseIdx, compIdx) << " flash vs "
|
|
|
|
<< fsRef.moleFraction(phaseIdx, compIdx) << " reference"
|
2016-04-17 04:28:21 -05:00
|
|
|
<< " error=" << error);
|
2013-05-30 11:44:10 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
template <class Scalar, class FluidSystem, class MaterialLaw, class FluidState>
|
2017-01-18 10:04:59 -06:00
|
|
|
void checkImmiscibleFlash(const FluidState& fsRef,
|
|
|
|
typename MaterialLaw::Params& matParams)
|
2013-05-30 11:44:10 -05:00
|
|
|
{
|
|
|
|
enum { numPhases = FluidSystem::numPhases };
|
|
|
|
enum { numComponents = FluidSystem::numComponents };
|
|
|
|
typedef Dune::FieldVector<Scalar, numComponents> ComponentVector;
|
|
|
|
|
|
|
|
// calculate the total amount of stuff in the reference fluid
|
|
|
|
// phase
|
|
|
|
ComponentVector globalMolarities(0.0);
|
2015-09-22 04:20:55 -05:00
|
|
|
for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
|
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
|
2013-05-30 11:44:10 -05:00
|
|
|
globalMolarities[compIdx] +=
|
|
|
|
fsRef.saturation(phaseIdx)*fsRef.molarity(phaseIdx, compIdx);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
// initialize the fluid state for the flash calculation
|
|
|
|
typedef Opm::ImmiscibleFlash<Scalar, FluidSystem> ImmiscibleFlash;
|
|
|
|
FluidState fsFlash;
|
|
|
|
|
|
|
|
fsFlash.setTemperature(fsRef.temperature(/*phaseIdx=*/0));
|
|
|
|
|
|
|
|
// run the flash calculation
|
2016-04-17 04:28:04 -05:00
|
|
|
ImmiscibleFlash::guessInitial(fsFlash, globalMolarities);
|
|
|
|
typename FluidSystem::template ParameterCache<typename FluidState::Scalar> paramCache;
|
|
|
|
ImmiscibleFlash::template solve<MaterialLaw>(fsFlash, matParams, paramCache, globalMolarities);
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
// compare the "flashed" fluid state with the reference one
|
|
|
|
checkSame<Scalar>(fsRef, fsFlash);
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
template <class Scalar, class FluidSystem, class MaterialLaw, class FluidState>
|
2017-01-18 10:04:59 -06:00
|
|
|
void completeReferenceFluidState(FluidState& fs,
|
|
|
|
typename MaterialLaw::Params& matParams,
|
2015-09-22 04:20:55 -05:00
|
|
|
unsigned refPhaseIdx)
|
2013-05-30 11:44:10 -05:00
|
|
|
{
|
|
|
|
enum { numPhases = FluidSystem::numPhases };
|
|
|
|
typedef Dune::FieldVector<Scalar, numPhases> PhaseVector;
|
|
|
|
|
2015-09-22 04:20:55 -05:00
|
|
|
unsigned otherPhaseIdx = 1 - refPhaseIdx;
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
// calculate the other saturation
|
|
|
|
fs.setSaturation(otherPhaseIdx, 1.0 - fs.saturation(refPhaseIdx));
|
|
|
|
|
|
|
|
// calulate the capillary pressure
|
|
|
|
PhaseVector pC;
|
|
|
|
MaterialLaw::capillaryPressures(pC, matParams, fs);
|
|
|
|
fs.setPressure(otherPhaseIdx,
|
|
|
|
fs.pressure(refPhaseIdx)
|
|
|
|
+ (pC[otherPhaseIdx] - pC[refPhaseIdx]));
|
|
|
|
|
|
|
|
// set all phase densities
|
2016-04-17 04:28:04 -05:00
|
|
|
typename FluidSystem::template ParameterCache<typename FluidState::Scalar> paramCache;
|
2013-05-30 11:44:10 -05:00
|
|
|
paramCache.updateAll(fs);
|
2015-09-22 04:20:55 -05:00
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
|
2013-05-30 11:44:10 -05:00
|
|
|
Scalar rho = FluidSystem::density(fs, paramCache, phaseIdx);
|
|
|
|
fs.setDensity(phaseIdx, rho);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2016-01-14 18:28:07 -06:00
|
|
|
template <class Scalar>
|
|
|
|
inline void testAll()
|
2013-05-30 11:44:10 -05:00
|
|
|
{
|
|
|
|
typedef Opm::FluidSystems::H2ON2<Scalar> FluidSystem;
|
|
|
|
typedef Opm::ImmiscibleFluidState<Scalar, FluidSystem> ImmiscibleFluidState;
|
|
|
|
|
|
|
|
enum { numPhases = FluidSystem::numPhases };
|
|
|
|
enum { numComponents = FluidSystem::numComponents };
|
2014-04-03 09:34:23 -05:00
|
|
|
enum { liquidPhaseIdx = FluidSystem::liquidPhaseIdx };
|
|
|
|
enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
enum { H2OIdx = FluidSystem::H2OIdx };
|
|
|
|
enum { N2Idx = FluidSystem::N2Idx };
|
|
|
|
|
2014-04-03 09:34:23 -05:00
|
|
|
typedef Opm::TwoPhaseMaterialTraits<Scalar, liquidPhaseIdx, gasPhaseIdx> MaterialLawTraits;
|
2013-11-06 07:32:22 -06:00
|
|
|
typedef Opm::RegularizedBrooksCorey<MaterialLawTraits> EffMaterialLaw;
|
|
|
|
typedef Opm::EffToAbsLaw<EffMaterialLaw> MaterialLaw;
|
2016-01-14 18:28:07 -06:00
|
|
|
typedef typename MaterialLaw::Params MaterialLawParams;
|
2013-05-30 11:44:10 -05:00
|
|
|
|
2016-04-17 04:28:21 -05:00
|
|
|
std::cout << "---- using " << Dune::className<Scalar>() << " as scalar ----\n";
|
2013-05-30 11:44:10 -05:00
|
|
|
Scalar T = 273.15 + 25;
|
|
|
|
|
|
|
|
// initialize the tables of the fluid system
|
|
|
|
Scalar Tmin = T - 1.0;
|
|
|
|
Scalar Tmax = T + 1.0;
|
2015-09-22 04:20:55 -05:00
|
|
|
unsigned nT = 3;
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
Scalar pmin = 0.0;
|
|
|
|
Scalar pmax = 1.25 * 2e6;
|
2015-09-22 04:20:55 -05:00
|
|
|
unsigned np = 100;
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
FluidSystem::init(Tmin, Tmax, nT, pmin, pmax, np);
|
|
|
|
|
|
|
|
// set the parameters for the capillary pressure law
|
|
|
|
MaterialLawParams matParams;
|
2014-04-03 09:34:23 -05:00
|
|
|
matParams.setResidualSaturation(MaterialLaw::wettingPhaseIdx, 0.0);
|
|
|
|
matParams.setResidualSaturation(MaterialLaw::nonWettingPhaseIdx, 0.0);
|
2013-11-06 07:32:22 -06:00
|
|
|
matParams.setEntryPressure(0);
|
2013-05-30 11:44:10 -05:00
|
|
|
matParams.setLambda(2.0);
|
2013-11-06 07:32:22 -06:00
|
|
|
matParams.finalize();
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
ImmiscibleFluidState fsRef;
|
|
|
|
|
|
|
|
// create an fluid state which is consistent
|
|
|
|
|
|
|
|
// set the fluid temperatures
|
|
|
|
fsRef.setTemperature(T);
|
|
|
|
|
|
|
|
////////////////
|
|
|
|
// only liquid
|
|
|
|
////////////////
|
|
|
|
std::cout << "testing single-phase liquid\n";
|
|
|
|
|
|
|
|
// set liquid saturation and pressure
|
2014-04-03 09:34:23 -05:00
|
|
|
fsRef.setSaturation(liquidPhaseIdx, 1.0);
|
|
|
|
fsRef.setPressure(liquidPhaseIdx, 1e6);
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
// set the remaining parameters of the reference fluid state
|
2014-04-03 09:34:23 -05:00
|
|
|
completeReferenceFluidState<Scalar, FluidSystem, MaterialLaw>(fsRef, matParams, liquidPhaseIdx);
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
// check the flash calculation
|
|
|
|
checkImmiscibleFlash<Scalar, FluidSystem, MaterialLaw>(fsRef, matParams);
|
|
|
|
|
|
|
|
////////////////
|
|
|
|
// only gas
|
|
|
|
////////////////
|
|
|
|
std::cout << "testing single-phase gas\n";
|
|
|
|
|
|
|
|
// set gas saturation and pressure
|
2014-04-03 09:34:23 -05:00
|
|
|
fsRef.setSaturation(gasPhaseIdx, 1.0);
|
|
|
|
fsRef.setPressure(gasPhaseIdx, 1e6);
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
// set the remaining parameters of the reference fluid state
|
2014-04-03 09:34:23 -05:00
|
|
|
completeReferenceFluidState<Scalar, FluidSystem, MaterialLaw>(fsRef, matParams, gasPhaseIdx);
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
// check the flash calculation
|
|
|
|
checkImmiscibleFlash<Scalar, FluidSystem, MaterialLaw>(fsRef, matParams);
|
|
|
|
|
|
|
|
////////////////
|
|
|
|
// both phases
|
|
|
|
////////////////
|
|
|
|
std::cout << "testing two-phase\n";
|
|
|
|
|
|
|
|
// set liquid saturation and pressure
|
2014-04-03 09:34:23 -05:00
|
|
|
fsRef.setSaturation(liquidPhaseIdx, 0.5);
|
|
|
|
fsRef.setPressure(liquidPhaseIdx, 1e6);
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
// set the remaining parameters of the reference fluid state
|
2014-04-03 09:34:23 -05:00
|
|
|
completeReferenceFluidState<Scalar, FluidSystem, MaterialLaw>(fsRef, matParams, liquidPhaseIdx);
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
// check the flash calculation
|
|
|
|
checkImmiscibleFlash<Scalar, FluidSystem, MaterialLaw>(fsRef, matParams);
|
|
|
|
|
|
|
|
////////////////
|
|
|
|
// with capillary pressure
|
|
|
|
////////////////
|
|
|
|
std::cout << "testing two-phase with capillary pressure\n";
|
|
|
|
|
|
|
|
MaterialLawParams matParams2;
|
2014-04-03 09:34:23 -05:00
|
|
|
matParams2.setResidualSaturation(MaterialLaw::wettingPhaseIdx, 0.0);
|
|
|
|
matParams2.setResidualSaturation(MaterialLaw::nonWettingPhaseIdx, 0.0);
|
2013-11-06 07:32:22 -06:00
|
|
|
matParams2.setEntryPressure(1e3);
|
2013-05-30 11:44:10 -05:00
|
|
|
matParams2.setLambda(2.0);
|
2013-11-06 07:32:22 -06:00
|
|
|
matParams2.finalize();
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
// set liquid saturation
|
2014-04-03 09:34:23 -05:00
|
|
|
fsRef.setSaturation(liquidPhaseIdx, 0.5);
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
// set pressure of the liquid phase
|
2014-04-03 09:34:23 -05:00
|
|
|
fsRef.setPressure(liquidPhaseIdx, 1e6);
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
// set the remaining parameters of the reference fluid state
|
2014-04-03 09:34:23 -05:00
|
|
|
completeReferenceFluidState<Scalar, FluidSystem, MaterialLaw>(fsRef, matParams2, liquidPhaseIdx);
|
2013-05-30 11:44:10 -05:00
|
|
|
|
|
|
|
// check the flash calculation
|
|
|
|
checkImmiscibleFlash<Scalar, FluidSystem, MaterialLaw>(fsRef, matParams2);
|
2016-01-14 18:28:07 -06:00
|
|
|
}
|
2013-05-30 11:44:10 -05:00
|
|
|
|
2016-04-17 04:28:06 -05:00
|
|
|
int main(int argc, char **argv)
|
2016-01-14 18:28:07 -06:00
|
|
|
{
|
2016-04-17 04:28:06 -05:00
|
|
|
Dune::MPIHelper::instance(argc, argv);
|
|
|
|
|
|
|
|
testAll<double>();
|
|
|
|
testAll<float>();
|
|
|
|
|
2013-05-30 11:44:10 -05:00
|
|
|
return 0;
|
|
|
|
}
|