// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
* Copyright (C) 2013 by Andreas Lauser *
* *
* This program 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. *
* *
* This program 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 this program. If not, see . *
*****************************************************************************/
/*!
* \file
*
* \brief This test makes sure that the API for fluid-matrix
* interactions is observed by all capillary pressure / relperm
* laws.
*/
#include "config.h"
// include all capillary pressure laws
#include
#include
#include
#include
#include
#include
#include
#include
#include
// include the helper classes to construct traits
#include
// include some fluid states
#include
#include
// include some fluid systems
#include
#include
// include some components
#include
#include
#include
// include the MPI header if available
#if HAVE_MPI
#include
#endif // HAVE_MPI
// class to call MPI_Init() on construction and MPI_Finalize() in the
// destructor
class MyMpiHelper
{
public:
MyMpiHelper(int &argc, char **&argv)
{
#if HAVE_MPI
MPI_Init(&argc, &argv);
#endif // HAVE_MPI
};
~MyMpiHelper()
{
#if HAVE_MPI
MPI_Finalize();
#endif // HAVE_MPI
};
};
// this function makes sure that a capillary pressure law adheres to
// the generic programming interface for such laws. This API _must_ be
// implemented by all capillary pressure laws. If there are no _very_
// strong reasons to do otherwise, numerical models should only use on
// this API.
template
void testGenericApi()
{
while (0) {
// ensure the presence of the required values
static const int numPhases = MaterialLaw::numPhases;
// check for the presence of the is*Dependent values
OPM_UNUSED static const bool isSaturationDependent = MaterialLaw::isSaturationDependent;
OPM_UNUSED static const bool isPressureDependent = MaterialLaw::isPressureDependent;
OPM_UNUSED static const bool isTemperatureDependent = MaterialLaw::isTemperatureDependent;
OPM_UNUSED static const bool isCompositionDependent = MaterialLaw::isCompositionDependent;
// Make sure that the Traits, Params and Scalar typedefs are
// exported by the material law
typedef typename MaterialLaw::Params Params;
typedef typename MaterialLaw::Traits Traits;
typedef typename MaterialLaw::Scalar Scalar;
typedef typename MaterialLaw::Traits::Scalar TraitsScalar;
static_assert(std::is_same::value,
"The traits and the material law must use the same type as scalar value");
static_assert(numPhases == Traits::numPhases,
"The traits and the material law must use the number of fluid phases");
// check the API of the parameter class. setting the actual
// parameter values is implementation specific. But all
// parameters must be default and copy constructible as well
// as exhibit the finalize() method!
Params params;
params.finalize();
const Params paramsConst(params);
// test the generic methods which need to be implemented by
// all material laws
const FluidState fs;
double destValues[numPhases];
MaterialLaw::capillaryPressures(destValues, paramsConst, fs);
MaterialLaw::saturations(destValues, paramsConst, fs);
MaterialLaw::relativePermeabilities(destValues, paramsConst, fs);
std::array dpc;
MaterialLaw::dCapillaryPressures_dSaturation(dpc, paramsConst, fs, /*phaseIdx=*/0);
MaterialLaw::dCapillaryPressures_dPressure(dpc, paramsConst, fs, /*phaseIdx=*/0);
MaterialLaw::dCapillaryPressures_dTemperature(dpc, paramsConst, fs);
MaterialLaw::dCapillaryPressures_dMoleFraction(dpc, paramsConst, fs, /*phaseIdx=*/0, /*compIdx=*/0);
std::array dkr;
MaterialLaw::dRelativePermeabilities_dSaturation(dkr, paramsConst, fs, /*phaseIdx=*/0);
MaterialLaw::dRelativePermeabilities_dPressure(dkr, paramsConst, fs, /*phaseIdx=*/0);
MaterialLaw::dRelativePermeabilities_dTemperature(dkr, paramsConst, fs);
MaterialLaw::dRelativePermeabilities_dMoleFraction(dkr, paramsConst, fs, /*phaseIdx=*/0, /*compIdx=*/0);
}
}
// this function makes ensures that a pressure law adheres to the
// covenience programming interface for two-phase material laws. The
// main purpose of this interface is to simplify the implementation of
// nested material laws.
template
void testTwoPhaseApi()
{
typedef typename MaterialLaw::Scalar Scalar;
while (0) {
static const int numPhases = MaterialLaw::numPhases;
static_assert(numPhases == 2,
"The number of fluid phases for a twophase "
"capillary pressure law must be 2");
static_assert(MaterialLaw::implementsTwoPhaseApi,
"This material law is expected to implement "
"the two-phase API!");
OPM_UNUSED static const int wPhaseIdx = MaterialLaw::wPhaseIdx;
OPM_UNUSED static const int nPhaseIdx = MaterialLaw::nPhaseIdx;
// make sure the two-phase specific methods are present
const FluidState fs;
const typename MaterialLaw::Params params;
OPM_UNUSED Scalar v;
v = MaterialLaw::pcnw(params, fs);
v = MaterialLaw::Sw(params, fs);
v = MaterialLaw::Sn(params, fs);
v = MaterialLaw::krw(params, fs);
v = MaterialLaw::krn(params, fs);
}
}
template
void testTwoPhaseSatApi()
{
typedef typename MaterialLaw::Scalar Scalar;
while (0) {
static_assert(MaterialLaw::implementsTwoPhaseSatApi,
"This material law is expected to implement "
"the two-phase saturation only API!");
static_assert(!MaterialLaw::isPressureDependent,
"Capillary pressure laws which implement the twophase saturation only "
"API cannot be dependent on the absolute phase pressures!");
static_assert(!MaterialLaw::isTemperatureDependent,
"Capillary pressure laws which implement the twophase saturation only "
"API cannot be dependent on temperature!");
static_assert(!MaterialLaw::isCompositionDependent,
"Capillary pressure laws which implement the twophase saturation only "
"API cannot be dependent on the phase compositions!");
OPM_UNUSED static const int numPhases = MaterialLaw::numPhases;
// make sure the two-phase specific methods are present
const typename MaterialLaw::Params params;
Scalar Sw = 0;
OPM_UNUSED Scalar v;
v = MaterialLaw::twoPhaseSatPcnw(params, Sw);
v = MaterialLaw::twoPhaseSatSw(params, Sw);
v = MaterialLaw::twoPhaseSatSn(params, Sw);
v = MaterialLaw::twoPhaseSatKrw(params, Sw);
v = MaterialLaw::twoPhaseSatKrn(params, Sw);
}
}
template
void testThreePhaseApi()
{
typedef typename MaterialLaw::Scalar Scalar;
while (0) {
static const int numPhases = MaterialLaw::numPhases;
static_assert(numPhases == 3,
"The number of fluid phases for a threephase "
"capillary pressure law must be 3");
OPM_UNUSED static const int wPhaseIdx = MaterialLaw::wPhaseIdx;
OPM_UNUSED static const int nPhaseIdx = MaterialLaw::nPhaseIdx;
OPM_UNUSED static const int gPhaseIdx = MaterialLaw::gPhaseIdx;
// make sure the two-phase specific methods are present
const FluidState fs;
const typename MaterialLaw::Params params;
OPM_UNUSED Scalar v;
v = MaterialLaw::pcnw(params, fs);
v = MaterialLaw::Sw(params, fs);
v = MaterialLaw::Sn(params, fs);
v = MaterialLaw::Sg(params, fs);
v = MaterialLaw::krw(params, fs);
v = MaterialLaw::krn(params, fs);
v = MaterialLaw::krg(params, fs);
}
}
template
void testThreePhaseSatApi()
{
}
int main(int argc, char **argv)
{
typedef double Scalar;
typedef Opm::SimpleH2O H2O;
typedef Opm::N2 N2;
typedef Opm::LiquidPhase Liquid;
typedef Opm::GasPhase Gas;
typedef Opm::FluidSystems::TwoPImmiscible TwoPFluidSystem;
typedef Opm::FluidSystems::BlackOil ThreePFluidSystem;
typedef Opm::TwoPhaseMaterialTraits TwoPhaseTraits;
typedef Opm::ThreePhaseMaterialTraits ThreePhaseTraits;
typedef Opm::ImmiscibleFluidState TwoPhaseFluidState;
typedef Opm::ImmiscibleFluidState ThreePhaseFluidState;
MyMpiHelper mpiHelper(argc, argv);
// test conformance to the capillary pressure APIs
{
typedef Opm::BrooksCorey MaterialLaw;
testGenericApi();
testTwoPhaseApi();
testTwoPhaseSatApi();
}
{
typedef Opm::LinearMaterial MaterialLaw;
testGenericApi();
testTwoPhaseApi();
testTwoPhaseSatApi();
typedef Opm::EffToAbsLaw TwoPAbsLaw;
testGenericApi();
testTwoPhaseApi();
testTwoPhaseSatApi();
typedef Opm::LinearMaterial ThreePMaterialLaw;
testGenericApi();
testThreePhaseApi();
//testThreePhaseSatApi();
typedef Opm::EffToAbsLaw ThreePAbsLaw;
testGenericApi();
testThreePhaseApi();
//testThreePhaseSatApi();
}
{
typedef Opm::BrooksCorey TwoPhaseMaterial;
typedef Opm::EclDefaultMaterial MaterialLaw;
testGenericApi();
testThreePhaseApi();
//testThreePhaseSatApi();
}
{
typedef Opm::NullMaterial MaterialLaw;
testGenericApi();
testTwoPhaseApi();
testTwoPhaseSatApi();
typedef Opm::NullMaterial ThreePMaterialLaw;
testGenericApi();
testThreePhaseApi();
//testThreePhaseSatApi();
}
{
typedef Opm::ParkerLenhard MaterialLaw;
testGenericApi();
testTwoPhaseApi();
testTwoPhaseSatApi();
}
{
typedef Opm::VanGenuchten MaterialLaw;
testGenericApi();
testTwoPhaseApi();
testTwoPhaseSatApi();
}
{
typedef Opm::RegularizedBrooksCorey MaterialLaw;
testGenericApi();
testTwoPhaseApi();
testTwoPhaseSatApi();
}
{
typedef Opm::RegularizedVanGenuchten MaterialLaw;
testGenericApi();
testTwoPhaseApi();
testTwoPhaseSatApi();
}
return 0;
}