Merge pull request #585 from totto82/trinemykk-smartwater

Trinemykk smartwater
This commit is contained in:
Joakim Hove 2020-01-11 06:50:44 +01:00 committed by GitHub
commit 98f41467c5
11 changed files with 517 additions and 13 deletions

View File

@ -0,0 +1,425 @@
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*
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/>.
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.
*/
/*!
* \file
*
* \brief Contains the classes required to extend the black-oil model by brine.
*/
#ifndef EWOMS_BLACK_OIL_BRINE_MODULE_HH
#define EWOMS_BLACK_OIL_BRINE_MODULE_HH
#include "blackoilproperties.hh"
#include <opm/models/common/quantitycallbacks.hh>
#include <opm/material/common/Tabulated1DFunction.hpp>
#include <opm/material/common/IntervalTabulated2DFunction.hpp>
#if HAVE_ECL_INPUT
#include <opm/parser/eclipse/Deck/Deck.hpp>
#include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/PvtwsaltTable.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/BrineDensityTable.hpp>
#endif
#include <opm/material/common/Valgrind.hpp>
#include <opm/material/common/Unused.hpp>
#include <opm/material/common/Exceptions.hpp>
#include <dune/common/fvector.hh>
#include <string>
#include <math.h>
namespace Opm {
/*!
* \ingroup BlackOil
* \brief Contains the high level supplements required to extend the black oil
* model by brine.
*/
template <class TypeTag, bool enableBrineV = GET_PROP_VALUE(TypeTag, EnableBrine)>
class BlackOilBrineModule
{
typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) PrimaryVariables;
typedef typename GET_PROP_TYPE(TypeTag, IntensiveQuantities) IntensiveQuantities;
typedef typename GET_PROP_TYPE(TypeTag, ExtensiveQuantities) ExtensiveQuantities;
typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
typedef typename GET_PROP_TYPE(TypeTag, Model) Model;
typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
typedef typename GET_PROP_TYPE(TypeTag, EqVector) EqVector;
typedef typename GET_PROP_TYPE(TypeTag, RateVector) RateVector;
typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
typedef Opm::MathToolbox<Evaluation> Toolbox;
typedef typename Opm::Tabulated1DFunction<Scalar> TabulatedFunction;
typedef typename Opm::IntervalTabulated2DFunction<Scalar> TabulatedTwoDFunction;
static constexpr unsigned saltConcentrationIdx = Indices::saltConcentrationIdx;
static constexpr unsigned contiBrineEqIdx = Indices::contiBrineEqIdx;
static constexpr unsigned waterPhaseIdx = FluidSystem::waterPhaseIdx;
static constexpr unsigned enableBrine = enableBrineV;
static constexpr unsigned numEq = GET_PROP_VALUE(TypeTag, NumEq);
static constexpr unsigned numPhases = FluidSystem::numPhases;
public:
#if HAVE_ECL_INPUT
/*!
* \brief Initialize all internal data structures needed by the brine module
*/
static void initFromDeck(const Opm::Deck& deck, const Opm::EclipseState& eclState)
{
// some sanity checks: if brine are enabled, the BRINE keyword must be
// present, if brine are disabled the keyword must not be present.
if (enableBrine && !deck.hasKeyword("BRINE")) {
throw std::runtime_error("Non-trivial brine treatment requested at compile time, but "
"the deck does not contain the BRINE keyword");
}
else if (!enableBrine && deck.hasKeyword("BRINE")) {
throw std::runtime_error("Brine treatment disabled at compile time, but the deck "
"contains the BRINE keyword");
}
if (!deck.hasKeyword("BRINE"))
return; // brine treatment is supposed to be disabled
const auto& tableManager = eclState.getTableManager();
unsigned numPvtRegions = tableManager.getTabdims().getNumPVTTables();
referencePressure_.resize(numPvtRegions);
const auto& pvtwsaltTables = tableManager.getPvtwSaltTables();
// initialize the objects which deal with the BDENSITY keyword
const auto& bdensityTables = tableManager.getBrineDensityTables();
if (!bdensityTables.empty()) {
bdensityTable_.resize(numPvtRegions);
assert(numPvtRegions == bdensityTables.size());
for (unsigned pvtRegionIdx = 0; pvtRegionIdx < numPvtRegions; ++ pvtRegionIdx) {
const auto& bdensityTable = bdensityTables[pvtRegionIdx];
const auto& pvtwsaltTable = pvtwsaltTables[pvtRegionIdx];
const auto& c = pvtwsaltTable.getSaltConcentrationColumn();
const auto& density = bdensityTable.getBrineDensityColumn();
bdensityTable_[pvtRegionIdx].setXYContainers(c, density);
}
}
}
#endif
/*!
* \brief Register all run-time parameters for the black-oil brine module.
*/
static void registerParameters()
{
if (!enableBrine)
// brine have been disabled at compile time
return;
}
static bool primaryVarApplies(unsigned pvIdx)
{
if (!enableBrine)
// brine have been disabled at compile time
return false;
return pvIdx == saltConcentrationIdx;
}
/*!
* \brief Assign the brine specific primary variables to a PrimaryVariables object
*/
template <class FluidState>
static void assignPrimaryVars(PrimaryVariables& priVars,
const FluidState& fluidState)
{
if (!enableBrine)
return;
priVars[saltConcentrationIdx] = fluidState.saltConcentration();
}
static std::string primaryVarName(unsigned pvIdx)
{
assert(primaryVarApplies(pvIdx));
return "saltConcentration";
}
static Scalar primaryVarWeight(unsigned pvIdx OPM_OPTIM_UNUSED)
{
assert(primaryVarApplies(pvIdx));
// TODO: it may be beneficial to chose this differently.
return static_cast<Scalar>(1.0);
}
static bool eqApplies(unsigned eqIdx)
{
if (!enableBrine)
return false;
return eqIdx == contiBrineEqIdx;
}
static std::string eqName(unsigned eqIdx)
{
assert(eqApplies(eqIdx));
return "conti^brine";
}
static Scalar eqWeight(unsigned eqIdx OPM_OPTIM_UNUSED)
{
assert(eqApplies(eqIdx));
// TODO: it may be beneficial to chose this differently.
return static_cast<Scalar>(1.0);
}
// must be called after water storage is computed
template <class LhsEval>
static void addStorage(Dune::FieldVector<LhsEval, numEq>& storage,
const IntensiveQuantities& intQuants)
{
if (!enableBrine)
return;
const auto& fs = intQuants.fluidState();
LhsEval surfaceVolumeWater =
Toolbox::template decay<LhsEval>(fs.saturation(waterPhaseIdx))
* Toolbox::template decay<LhsEval>(fs.invB(waterPhaseIdx))
* Toolbox::template decay<LhsEval>(intQuants.porosity());
// avoid singular matrix if no water is present.
surfaceVolumeWater = Opm::max(surfaceVolumeWater, 1e-10);
// Brine in water phase
const LhsEval massBrine = surfaceVolumeWater
* Toolbox::template decay<LhsEval>(fs.saltConcentration());
storage[contiBrineEqIdx] += massBrine;
}
static void computeFlux(RateVector& flux,
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
{
if (!enableBrine)
return;
const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
const unsigned upIdx = extQuants.upstreamIndex(FluidSystem::waterPhaseIdx);
const unsigned inIdx = extQuants.interiorIndex();
const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
if (upIdx == inIdx) {
flux[contiBrineEqIdx] =
extQuants.volumeFlux(waterPhaseIdx)
*up.fluidState().invB(waterPhaseIdx)
*up.fluidState().saltConcentration();
}
else {
flux[contiBrineEqIdx] =
extQuants.volumeFlux(waterPhaseIdx)
*Opm::decay<Scalar>(up.fluidState().invB(waterPhaseIdx))
*Opm::decay<Scalar>(up.fluidState().saltConcentration());
}
}
/*!
* \brief Return how much a Newton-Raphson update is considered an error
*/
static Scalar computeUpdateError(const PrimaryVariables& oldPv OPM_OPTIM_UNUSED,
const EqVector& delta OPM_OPTIM_UNUSED)
{
// do not consider consider the cange of Brine primary variables for
// convergence
// TODO: maybe this should be changed
return static_cast<Scalar>(0.0);
}
template <class DofEntity>
static void serializeEntity(const Model& model, std::ostream& outstream, const DofEntity& dof)
{
if (!enableBrine)
return;
unsigned dofIdx = model.dofMapper().index(dof);
const PrimaryVariables& priVars = model.solution(/*timeIdx=*/0)[dofIdx];
outstream << priVars[saltConcentrationIdx];
}
template <class DofEntity>
static void deserializeEntity(Model& model, std::istream& instream, const DofEntity& dof)
{
if (!enableBrine)
return;
unsigned dofIdx = model.dofMapper().index(dof);
PrimaryVariables& priVars0 = model.solution(/*timeIdx=*/0)[dofIdx];
PrimaryVariables& priVars1 = model.solution(/*timeIdx=*/1)[dofIdx];
instream >> priVars0[saltConcentrationIdx];
// set the primary variables for the beginning of the current time step.
priVars1[saltConcentrationIdx] = priVars0[saltConcentrationIdx];
}
static const Scalar& referencePressure(const ElementContext& elemCtx,
unsigned scvIdx,
unsigned timeIdx)
{
unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
return referencePressure_[pvtnumRegionIdx];
}
static const TabulatedFunction& bdensityTable(const ElementContext& elemCtx,
unsigned scvIdx,
unsigned timeIdx)
{
unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
return bdensityTable_[pvtnumRegionIdx];
}
static bool hasBDensityTables()
{
return !bdensityTable_.empty();
}
private:
static std::vector<TabulatedFunction> bdensityTable_;
static std::vector<Scalar> referencePressure_;
};
template <class TypeTag, bool enableBrineV>
std::vector<typename BlackOilBrineModule<TypeTag, enableBrineV>::TabulatedFunction>
BlackOilBrineModule<TypeTag, enableBrineV>::bdensityTable_;
template <class TypeTag, bool enableBrineV>
std::vector<typename BlackOilBrineModule<TypeTag, enableBrineV>::Scalar>
BlackOilBrineModule<TypeTag, enableBrineV>::referencePressure_;
/*!
* \ingroup BlackOil
* \class Ewoms::BlackOilBrineIntensiveQuantities
*
* \brief Provides the volumetric quantities required for the equations needed by the
* brine extension of the black-oil model.
*/
template <class TypeTag, bool enableBrineV = GET_PROP_VALUE(TypeTag, EnableBrine)>
class BlackOilBrineIntensiveQuantities
{
typedef typename GET_PROP_TYPE(TypeTag, IntensiveQuantities) Implementation;
typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) PrimaryVariables;
typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw;
typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
typedef BlackOilBrineModule<TypeTag> BrineModule;
enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
static constexpr int saltConcentrationIdx = Indices::saltConcentrationIdx;
static constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
static constexpr int oilPhaseIdx = FluidSystem::oilPhaseIdx;
static constexpr unsigned enableBrine = enableBrineV;
static constexpr int contiBrineEqIdx = Indices::contiBrineEqIdx;
public:
/*!
* \brief Update the intensive properties needed to handle brine from the
* primary variables
*
*/
void updateSaltConcentration_(const ElementContext& elemCtx,
unsigned dofIdx,
unsigned timeIdx)
{
const PrimaryVariables& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
auto& fs = asImp_().fluidState_;
// set saltconcentration
fs.setSaltConcentration(priVars.makeEvaluation(saltConcentrationIdx, timeIdx));
}
const Evaluation& saltConcentration() const
{ return saltConcentration_; }
const Evaluation& brineRefDensity() const
{ return refDensity_; }
protected:
Implementation& asImp_()
{ return *static_cast<Implementation*>(this); }
Evaluation saltConcentration_;
Evaluation refDensity_;
};
template <class TypeTag>
class BlackOilBrineIntensiveQuantities<TypeTag, false>
{
typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
public:
void updateSaltConcentration_(const ElementContext& elemCtx OPM_UNUSED,
unsigned dofIdx OPM_UNUSED,
unsigned timeIdx OPM_UNUSED)
{ }
const Evaluation& saltConcentration() const
{ throw std::runtime_error("saltConcentration() called but brine are disabled"); }
const Evaluation& brineRefDensity() const
{ throw std::runtime_error("brineRefDensity() called but brine are disabled"); }
};
} // namespace Ewoms
#endif

View File

@ -35,7 +35,7 @@ namespace Opm {
*
* \brief The primary variable and equation indices for the black-oil model.
*/
template <unsigned numSolventsV, unsigned numPolymersV, unsigned numEnergyV, bool enableFoam, unsigned PVOffset>
template <unsigned numSolventsV, unsigned numPolymersV, unsigned numEnergyV, bool enableFoam, bool enableBrine, unsigned PVOffset>
struct BlackOilIndices
{
//! Number of phases active at all times
@ -67,8 +67,11 @@ struct BlackOilIndices
//! Number of foam equations to be considered
static const int numFoam = enableFoam? 1 : 0;
//! Number of salt equations to be considered
static const int numBrine = enableBrine? 1 : 0;
//! The number of equations
static const int numEq = numPhases + numSolvents + numPolymers + numEnergy + numFoam;
static const int numEq = numPhases + numSolvents + numPolymers + numEnergy + numFoam + numBrine;
//! \brief returns the index of "active" component
static constexpr unsigned canonicalToActiveComponentIndex(unsigned compIdx)
@ -113,9 +116,13 @@ struct BlackOilIndices
static const int foamConcentrationIdx =
enableFoam ? PVOffset + numPhases + numSolvents + numPolymers : -1000;
//! Index of the primary variable for the brine
static const int saltConcentrationIdx =
enableBrine ? PVOffset + numPhases + numSolvents + numPolymers + numFoam : -1000;
//! Index of the primary variable for temperature
static const int temperatureIdx =
enableEnergy ? PVOffset + numPhases + numSolvents + numPolymers + numFoam : - 1000;
enableEnergy ? PVOffset + numPhases + numSolvents + numPolymers + numFoam + numBrine : - 1000;
////////
@ -142,10 +149,14 @@ struct BlackOilIndices
static const int contiFoamEqIdx =
enableFoam ? PVOffset + numPhases + numSolvents + numPolymers : -1000;
//! Index of the continuity equation for the salt water component
static const int contiBrineEqIdx =
enableBrine ? PVOffset + numPhases + numSolvents + numPolymers + numFoam : -1000;
//! Index of the continuity equation for energy
static const int contiEnergyEqIdx =
enableEnergy ? PVOffset + numPhases + numSolvents + numPolymers + numFoam : -1000;
enableEnergy ? PVOffset + numPhases + numSolvents + numPolymers + numFoam + numBrine: -1000;
};
} // namespace Opm

View File

@ -32,6 +32,7 @@
#include "blackoilsolventmodules.hh"
#include "blackoilpolymermodules.hh"
#include "blackoilfoammodules.hh"
#include "blackoilbrinemodules.hh"
#include "blackoilenergymodules.hh"
#include <opm/material/fluidstates/BlackOilFluidState.hpp>
#include <opm/material/common/Valgrind.hpp>
@ -56,6 +57,7 @@ class BlackOilIntensiveQuantities
, public BlackOilSolventIntensiveQuantities<TypeTag>
, public BlackOilPolymerIntensiveQuantities<TypeTag>
, public BlackOilFoamIntensiveQuantities<TypeTag>
, public BlackOilBrineIntensiveQuantities<TypeTag>
, public BlackOilEnergyIntensiveQuantities<TypeTag>
{
typedef typename GET_PROP_TYPE(TypeTag, DiscIntensiveQuantities) ParentType;
@ -75,6 +77,7 @@ class BlackOilIntensiveQuantities
enum { enableSolvent = GET_PROP_VALUE(TypeTag, EnableSolvent) };
enum { enablePolymer = GET_PROP_VALUE(TypeTag, EnablePolymer) };
enum { enableFoam = GET_PROP_VALUE(TypeTag, EnableFoam) };
enum { enableBrine = GET_PROP_VALUE(TypeTag, EnableBrine) };
enum { enableTemperature = GET_PROP_VALUE(TypeTag, EnableTemperature) };
enum { enableEnergy = GET_PROP_VALUE(TypeTag, EnableEnergy) };
enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
@ -94,7 +97,7 @@ class BlackOilIntensiveQuantities
typedef Opm::MathToolbox<Evaluation> Toolbox;
typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimMatrix;
typedef typename FluxModule::FluxIntensiveQuantities FluxIntensiveQuantities;
typedef Opm::BlackOilFluidState<Evaluation, FluidSystem, enableTemperature, enableEnergy, compositionSwitchEnabled, Indices::numPhases > FluidState;
typedef Opm::BlackOilFluidState<Evaluation, FluidSystem, enableTemperature, enableEnergy, compositionSwitchEnabled, enableBrine, Indices::numPhases > FluidState;
public:
BlackOilIntensiveQuantities()
@ -125,6 +128,8 @@ public:
unsigned pvtRegionIdx = priVars.pvtRegionIndex();
fluidState_.setPvtRegionIndex(pvtRegionIdx);
asImp_().updateSaltConcentration_(elemCtx, dofIdx, timeIdx);
// extract the water and the gas saturations for convenience
Evaluation Sw = 0.0;
if (waterEnabled) {
@ -444,6 +449,8 @@ private:
friend BlackOilPolymerIntensiveQuantities<TypeTag>;
friend BlackOilEnergyIntensiveQuantities<TypeTag>;
friend BlackOilFoamIntensiveQuantities<TypeTag>;
friend BlackOilBrineIntensiveQuantities<TypeTag>;
Implementation& asImp_()
{ return *static_cast<Implementation*>(this); }

View File

@ -33,6 +33,7 @@
#include "blackoilpolymermodules.hh"
#include "blackoilenergymodules.hh"
#include "blackoilfoammodules.hh"
#include "blackoilbrinemodules.hh"
#include <opm/material/fluidstates/BlackOilFluidState.hpp>
@ -82,6 +83,7 @@ class BlackOilLocalResidual : public GET_PROP_TYPE(TypeTag, DiscLocalResidual)
typedef BlackOilPolymerModule<TypeTag> PolymerModule;
typedef BlackOilEnergyModule<TypeTag> EnergyModule;
typedef BlackOilFoamModule<TypeTag> FoamModule;
typedef BlackOilBrineModule<TypeTag> BrineModule;
public:
/*!
@ -149,6 +151,9 @@ public:
// deal with foam (if present)
FoamModule::addStorage(storage, intQuants);
// deal with salt (if present)
BrineModule::addStorage(storage, intQuants);
}
/*!
@ -189,6 +194,9 @@ public:
// deal with foam (if present)
FoamModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
// deal with salt (if present)
BrineModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
}
/*!

View File

@ -44,6 +44,7 @@
#include "blackoilsolventmodules.hh"
#include "blackoilpolymermodules.hh"
#include "blackoilfoammodules.hh"
#include "blackoilbrinemodules.hh"
#include "blackoildarcyfluxmodule.hh"
#include <opm/models/common/multiphasebasemodel.hh>
@ -113,6 +114,7 @@ SET_TYPE_PROP(BlackOilModel, Indices,
GET_PROP_VALUE(TypeTag, EnablePolymer),
GET_PROP_VALUE(TypeTag, EnableEnergy),
GET_PROP_VALUE(TypeTag, EnableFoam),
GET_PROP_VALUE(TypeTag, EnableBrine),
/*PVOffset=*/0>);
//! Set the fluid system to the black-oil fluid system by default
@ -131,6 +133,7 @@ SET_BOOL_PROP(BlackOilModel, EnableSolvent, false);
SET_BOOL_PROP(BlackOilModel, EnablePolymer, false);
SET_BOOL_PROP(BlackOilModel, EnablePolymerMW, false);
SET_BOOL_PROP(BlackOilModel, EnableFoam, false);
SET_BOOL_PROP(BlackOilModel, EnableBrine, false);
//! By default, the blackoil model is isothermal and does not conserve energy
SET_BOOL_PROP(BlackOilModel, EnableTemperature, false);

View File

@ -188,6 +188,7 @@ protected:
static constexpr bool enablePolymerWeight = Indices::polymerMoleWeightIdx >= 0;
static constexpr bool enableEnergy = Indices::temperatureIdx >= 0;
static constexpr bool enableFoam = Indices::foamConcentrationIdx >= 0;
static constexpr bool enableBrine = Indices::saltConcentrationIdx >= 0;
currentValue.checkDefined();
Opm::Valgrind::CheckDefined(update);
@ -289,6 +290,10 @@ protected:
if (enableFoam && pvIdx == Indices::foamConcentrationIdx)
nextValue[pvIdx] = std::max(nextValue[pvIdx], 0.0);
// keep the salt concentration above 0
if (enableBrine && pvIdx == Indices::saltConcentrationIdx)
nextValue[pvIdx] = std::max(nextValue[pvIdx], 0.0);
// keep the temperature above 100 and below 1000 Kelvin
if (enableEnergy && pvIdx == Indices::temperatureIdx)
nextValue[pvIdx] = std::max(std::min(nextValue[pvIdx], 1000.0), 100.0);

View File

@ -37,7 +37,7 @@ namespace Opm {
*
* \brief The primary variable and equation indices for the black-oil model.
*/
template <unsigned numSolventsV, unsigned numPolymersV, unsigned numEnergyV, bool enableFoam, unsigned PVOffset, unsigned canonicalCompIdx>
template <unsigned numSolventsV, unsigned numPolymersV, unsigned numEnergyV, bool enableFoam, bool enableBrine, unsigned PVOffset, unsigned canonicalCompIdx>
struct BlackOilOnePhaseIndices
{
//! Is phase enabled or not
@ -66,11 +66,14 @@ struct BlackOilOnePhaseIndices
//! Number of foam equations to be considered
static const int numFoam = enableFoam? 1 : 0;
//! Number of salt equations to be considered
static const int numBrine = enableBrine? 1 : 0;
//! The number of fluid phases
static const int numPhases = 1;
//! The number of equations
static const int numEq = numPhases + numSolvents + numPolymers + numEnergy + numFoam;
static const int numEq = numPhases + numSolvents + numPolymers + numEnergy + numFoam + numBrine;
//////////////////////////////
// Primary variable indices
@ -106,9 +109,13 @@ struct BlackOilOnePhaseIndices
static const int foamConcentrationIdx =
enableFoam ? PVOffset + numPhases + numSolvents + numPolymers : -1000;
//! Index of the primary variable for the salt
static const int saltConcentrationIdx =
enableBrine ? PVOffset + numPhases + numSolvents + numPolymers + numFoam : -1000;
//! Index of the primary variable for temperature
static const int temperatureIdx =
enableEnergy ? PVOffset + numPhases + numSolvents + numPolymers + numFoam: - 1000;
enableEnergy ? PVOffset + numPhases + numSolvents + numPolymers + numFoam + numBrine: - 1000;
//////////////////////
// Equation indices
@ -153,9 +160,13 @@ struct BlackOilOnePhaseIndices
static const int contiFoamEqIdx =
enableFoam ? PVOffset + numPhases + numSolvents + numPolymers : -1000;
//! Index of the continuity equation for the salt component
static const int contiBrineEqIdx =
enableBrine ? PVOffset + numPhases + numSolvents + numPolymers + numFoam : -1000;
//! Index of the continuity equation for energy
static const int contiEnergyEqIdx =
enableEnergy ? PVOffset + numPhases + numSolvents + numPolymers + numFoam : -1000;
enableEnergy ? PVOffset + numPhases + numSolvents + numPolymers + numFoam + numBrine: -1000;
};
} // namespace Opm

View File

@ -33,6 +33,7 @@
#include "blackoilpolymermodules.hh"
#include "blackoilenergymodules.hh"
#include "blackoilfoammodules.hh"
#include "blackoilbrinemodules.hh"
#include <opm/models/discretization/common/fvbaseprimaryvariables.hh>
@ -51,6 +52,9 @@ class BlackOilSolventModule;
template <class TypeTag, bool enablePolymer>
class BlackOilPolymerModule;
template <class TypeTag, bool enableBrine>
class BlackOilBrineModule;
/*!
* \ingroup BlackOilModel
*
@ -92,6 +96,7 @@ class BlackOilPrimaryVariables : public FvBasePrimaryVariables<TypeTag>
enum { enableSolvent = GET_PROP_VALUE(TypeTag, EnableSolvent) };
enum { enablePolymer = GET_PROP_VALUE(TypeTag, EnablePolymer) };
enum { enableFoam = GET_PROP_VALUE(TypeTag, EnableFoam) };
enum { enableBrine = GET_PROP_VALUE(TypeTag, EnableBrine) };
enum { enableEnergy = GET_PROP_VALUE(TypeTag, EnableEnergy) };
enum { gasCompIdx = FluidSystem::gasCompIdx };
enum { waterCompIdx = FluidSystem::waterCompIdx };
@ -103,6 +108,7 @@ class BlackOilPrimaryVariables : public FvBasePrimaryVariables<TypeTag>
typedef BlackOilPolymerModule<TypeTag, enablePolymer> PolymerModule;
typedef BlackOilEnergyModule<TypeTag, enableEnergy> EnergyModule;
typedef BlackOilFoamModule<TypeTag, enableFoam> FoamModule;
typedef BlackOilBrineModule<TypeTag, enableBrine> BrineModule;
static_assert(numPhases == 3, "The black-oil model assumes three phases!");
static_assert(numComponents == 3, "The black-oil model assumes three components!");
@ -619,6 +625,14 @@ private:
return (*this)[Indices::foamConcentrationIdx];
}
Scalar saltConcentration_() const
{
if (!enableBrine)
return 0.0;
return (*this)[Indices::saltConcentrationIdx];
}
Scalar temperature_() const
{
if (!enableEnergy)

View File

@ -54,6 +54,9 @@ NEW_PROP_TAG(EnablePolymerMW);
NEW_PROP_TAG(BlackoilConserveSurfaceVolume);
//! Enable the ECL-blackoil extension for foam
NEW_PROP_TAG(EnableFoam);
//! Enable the ECL-blackoil extension for salt
NEW_PROP_TAG(EnableBrine);
//! Allow the spatial and temporal domains to exhibit non-constant temperature
//! in the black-oil model

View File

@ -59,6 +59,7 @@ class BlackOilRateVector
typedef BlackOilSolventModule<TypeTag> SolventModule;
typedef BlackOilPolymerModule<TypeTag> PolymerModule;
typedef BlackOilFoamModule<TypeTag> FoamModule;
typedef BlackOilBrineModule<TypeTag> BrineModule;
enum { numEq = GET_PROP_VALUE(TypeTag, NumEq) };
enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
@ -69,6 +70,7 @@ class BlackOilRateVector
enum { enablePolymer = GET_PROP_VALUE(TypeTag, EnablePolymer) };
enum { enablePolymerMolarWeight = GET_PROP_VALUE(TypeTag, EnablePolymerMW) };
enum { enableFoam = GET_PROP_VALUE(TypeTag, EnableFoam) };
enum { enableBrine = GET_PROP_VALUE(TypeTag, EnableBrine) };
typedef Opm::MathToolbox<Evaluation> Toolbox;
typedef Dune::FieldVector<Evaluation, numEq> ParentType;
@ -139,6 +141,10 @@ public:
throw std::logic_error("setMolarRate() not implemented for foam");
}
if ( enableBrine ) {
throw std::logic_error("setMolarRate() not implemented for salt water");
}
// convert to "surface volume" if requested
if (GET_PROP_VALUE(TypeTag, BlackoilConserveSurfaceVolume)) {
if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {

View File

@ -37,7 +37,7 @@ namespace Opm {
*
* \brief The primary variable and equation indices for the black-oil model.
*/
template <unsigned numSolventsV, unsigned numPolymersV, unsigned numEnergyV, bool enableFoam, unsigned PVOffset, unsigned disabledCanonicalCompIdx>
template <unsigned numSolventsV, unsigned numPolymersV, unsigned numEnergyV, bool enableFoam, bool enableBrine, unsigned PVOffset, unsigned disabledCanonicalCompIdx>
struct BlackOilTwoPhaseIndices
{
//! Is phase enabled or not
@ -66,11 +66,14 @@ struct BlackOilTwoPhaseIndices
//! Number of foam equations to be considered
static const int numFoam = enableFoam? 1 : 0;
//! Number of salt equations to be considered
static const int numBrine = enableBrine? 1 : 0;
//! The number of fluid phases
static const int numPhases = 2;
//! The number of equations
static const int numEq = numPhases + numSolvents + numPolymers + numEnergy + numFoam;
static const int numEq = numPhases + numSolvents + numPolymers + numEnergy + numFoam + numBrine;
//////////////////////////////
// Primary variable indices
@ -106,9 +109,13 @@ struct BlackOilTwoPhaseIndices
static const int foamConcentrationIdx =
enableFoam ? PVOffset + numPhases + numSolvents + numPolymers : -1000;
//! Index of the primary variable for the salt
static const int saltConcentrationIdx =
enableBrine ? PVOffset + numPhases + numSolvents + numPolymers + numFoam : -1000;
//! Index of the primary variable for temperature
static const int temperatureIdx =
enableEnergy ? PVOffset + numPhases + numSolvents + numPolymers + numFoam : - 1000;
enableEnergy ? PVOffset + numPhases + numSolvents + numPolymers + numFoam + numBrine : - 1000;
//////////////////////
// Equation indices
@ -171,9 +178,13 @@ struct BlackOilTwoPhaseIndices
static const int contiFoamEqIdx =
enableFoam ? PVOffset + numPhases + numSolvents + numPolymers : -1000;
//! Index of the continuity equation for the salt component
static const int contiBrineEqIdx =
enableBrine ? PVOffset + numPhases + numSolvents + numPolymers + numFoam : -1000;
//! Index of the continuity equation for energy
static const int contiEnergyEqIdx =
enableEnergy ? PVOffset + numPhases + numSolvents + numPolymers + numFoam : -1000;
enableEnergy ? PVOffset + numPhases + numSolvents + numPolymers + numFoam + numBrine : -1000;
};
} // namespace Opm