// -*- 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 .
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
*
* \copydoc Opm::Co2InjectionProblem
*/
#ifndef EWOMS_CO2_INJECTION_PROBLEM_HH
#define EWOMS_CO2_INJECTION_PROBLEM_HH
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
namespace Opm {
//! \cond SKIP_THIS
template
class Co2InjectionProblem;
//! \endcond
}
namespace Opm::Properties {
namespace TTag {
struct Co2InjectionBaseProblem {};
}
// Set the grid type
template
struct Grid { using type = Dune::YaspGrid<2>; };
// Set the problem property
template
struct Problem
{ using type = Opm::Co2InjectionProblem; };
// Set fluid configuration
template
struct FluidSystem
{
private:
using Scalar = GetPropType;
public:
using type = Opm::BrineCO2FluidSystem;
//using type = Opm::H2ON2FluidSystem;
};
// Set the material Law
template
struct MaterialLaw
{
private:
using FluidSystem = GetPropType;
enum { liquidPhaseIdx = FluidSystem::liquidPhaseIdx };
enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
using Scalar = GetPropType;
using Traits = Opm::TwoPhaseMaterialTraits;
// define the material law which is parameterized by effective
// saturations
using EffMaterialLaw = Opm::RegularizedBrooksCorey;
public:
// define the material law parameterized by absolute saturations
using type = Opm::EffToAbsLaw;
};
// Set the thermal conduction law
template
struct ThermalConductionLaw
{
private:
using Scalar = GetPropType;
using FluidSystem = GetPropType;
public:
// define the material law parameterized by absolute saturations
using type = Opm::SomertonThermalConductionLaw;
};
// set the energy storage law for the solid phase
template
struct SolidEnergyLaw
{ using type = Opm::ConstantSolidHeatCapLaw>; };
// Use the algebraic multi-grid linear solver for this problem
template
struct LinearSolverSplice { using type = TTag::ParallelAmgLinearSolver; };
} // namespace Opm::Properties
namespace Opm::Parameters {
struct FluidSystemNumPressure { static constexpr unsigned value = 100; };
struct FluidSystemNumTemperature { static constexpr unsigned value = 100; };
template
struct FluidSystemPressureHigh { static constexpr Scalar value = 4e7; };
template
struct FluidSystemPressureLow { static constexpr Scalar value = 3e7; };
template
struct FluidSystemTemperatureHigh { static constexpr Scalar value = 500.0; };
template
struct FluidSystemTemperatureLow { static constexpr Scalar value = 290.0; };
template
struct MaxDepth { static constexpr Scalar value = 2500.0; };
struct SimulationName { static constexpr auto value = "co2injection"; };
template
struct Temperature { static constexpr Scalar value = 293.15; };
} // namespace Opm::Parameters
namespace Opm {
double Co2InjectionTolerance = 1e-8;
/*!
* \ingroup TestProblems
*
* \brief Problem where \f$CO_2\f$ is injected under a low permeable
* layer at a depth of 2700m.
*
* The domain is sized 60m times 40m and consists of two layers, one
* which is moderately permeable (\f$K = 10^{-12}\;m^2\f$) for \f$ y >
* 22\; m\f$ and one with a lower intrinsic permeablility (\f$
* K=10^{-13}\;m^2\f$) in the rest of the domain.
*
* \f$CO_2\f$ gets injected by means of a forced-flow boundary
* condition into water-filled aquifer, which is situated 2700m below
* sea level, at the lower-right boundary (\f$5m
class Co2InjectionProblem : public GetPropType
{
using ParentType = GetPropType;
using Scalar = GetPropType;
using Evaluation = GetPropType;
using GridView = GetPropType;
using FluidSystem = GetPropType;
enum { dim = GridView::dimension };
enum { dimWorld = GridView::dimensionworld };
// copy some indices for convenience
using Indices = GetPropType;
enum { numPhases = FluidSystem::numPhases };
enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
enum { liquidPhaseIdx = FluidSystem::liquidPhaseIdx };
enum { CO2Idx = FluidSystem::CO2Idx };
enum { BrineIdx = FluidSystem::BrineIdx };
static constexpr int conti0EqIdx = Indices::conti0EqIdx;
static constexpr int contiCO2EqIdx = conti0EqIdx + CO2Idx;
using PrimaryVariables = GetPropType;
using RateVector = GetPropType;
using BoundaryRateVector = GetPropType;
using MaterialLaw = GetPropType;
using Simulator = GetPropType;
using Model = GetPropType;
using MaterialLawParams = GetPropType;
using ThermalConductionLaw = GetPropType;
using SolidEnergyLawParams = GetPropType;
using ThermalConductionLawParams = typename ThermalConductionLaw::Params;
using Toolbox = Opm::MathToolbox;
using CoordScalar = typename GridView::ctype;
using GlobalPosition = Dune::FieldVector;
using DimMatrix = Dune::FieldMatrix;
public:
/*!
* \copydoc Doxygen::defaultProblemConstructor
*/
Co2InjectionProblem(Simulator& simulator)
: ParentType(simulator)
{ }
/*!
* \copydoc FvBaseProblem::finishInit
*/
void finishInit()
{
ParentType::finishInit();
eps_ = 1e-6;
temperatureLow_ = Parameters::Get>();
temperatureHigh_ = Parameters::Get>();
nTemperature_ = Parameters::Get();
pressureLow_ = Parameters::Get>();
pressureHigh_ = Parameters::Get>();
nPressure_ = Parameters::Get();
maxDepth_ = Parameters::Get>();
temperature_ = Parameters::Get>();
// initialize the tables of the fluid system
// FluidSystem::init();
FluidSystem::init(/*Tmin=*/temperatureLow_,
/*Tmax=*/temperatureHigh_,
/*nT=*/nTemperature_,
/*pmin=*/pressureLow_,
/*pmax=*/pressureHigh_,
/*np=*/nPressure_);
fineLayerBottom_ = 22.0;
// intrinsic permeabilities
fineK_ = this->toDimMatrix_(1e-13);
coarseK_ = this->toDimMatrix_(1e-12);
// porosities
finePorosity_ = 0.3;
coarsePorosity_ = 0.3;
// residual saturations
fineMaterialParams_.setResidualSaturation(liquidPhaseIdx, 0.2);
fineMaterialParams_.setResidualSaturation(gasPhaseIdx, 0.0);
coarseMaterialParams_.setResidualSaturation(liquidPhaseIdx, 0.2);
coarseMaterialParams_.setResidualSaturation(gasPhaseIdx, 0.0);
// parameters for the Brooks-Corey law
fineMaterialParams_.setEntryPressure(1e4);
coarseMaterialParams_.setEntryPressure(5e3);
fineMaterialParams_.setLambda(2.0);
coarseMaterialParams_.setLambda(2.0);
fineMaterialParams_.finalize();
coarseMaterialParams_.finalize();
// parameters for the somerton law thermal conduction
computeThermalCondParams_(fineThermalCondParams_, finePorosity_);
computeThermalCondParams_(coarseThermalCondParams_, coarsePorosity_);
// assume constant heat capacity and granite
solidEnergyLawParams_.setSolidHeatCapacity(790.0 // specific heat capacity of granite [J / (kg K)]
* 2700.0); // density of granite [kg/m^3]
solidEnergyLawParams_.finalize();
}
/*!
* \copydoc FvBaseMultiPhaseProblem::registerParameters
*/
static void registerParameters()
{
ParentType::registerParameters();
Parameters::Register>
("The lower temperature [K] for tabulation of the fluid system");
Parameters::Register>
("The upper temperature [K] for tabulation of the fluid system");
Parameters::Register
("The number of intervals between the lower and upper temperature");
Parameters::Register>
("The lower pressure [Pa] for tabulation of the fluid system");
Parameters::Register>
("The upper pressure [Pa] for tabulation of the fluid system");
Parameters::Register
("The number of intervals between the lower and upper pressure");
Parameters::Register>
("The temperature [K] in the reservoir");
Parameters::Register>
("The maximum depth [m] of the reservoir");
Parameters::Register
("The name of the simulation used for the output files");
Parameters::SetDefault("data/co2injection.dgf");
Parameters::SetDefault>(1e4);
Parameters::SetDefault>(250);
Parameters::SetDefault>(Scalar{Co2InjectionTolerance});
Parameters::SetDefault(true);
}
/*!
* \name Problem parameters
*/
//! \{
/*!
* \copydoc FvBaseProblem::name
*/
std::string name() const
{
std::ostringstream oss;
oss << Parameters::Get()
<< "_" << Model::name();
if (getPropValue())
oss << "_ni";
oss << "_" << Model::discretizationName();
return oss.str();
}
/*!
* \copydoc FvBaseProblem::endTimeStep
*/
void endTimeStep()
{
#ifndef NDEBUG
Scalar tol = this->model().newtonMethod().tolerance()*1e5;
this->model().checkConservativeness(tol);
// Calculate storage terms
PrimaryVariables storageL, storageG;
this->model().globalPhaseStorage(storageL, /*phaseIdx=*/0);
this->model().globalPhaseStorage(storageG, /*phaseIdx=*/1);
// Write mass balance information for rank 0
if (this->gridView().comm().rank() == 0) {
std::cout << "Storage: liquid=[" << storageL << "]"
<< " gas=[" << storageG << "]\n" << std::flush;
}
#endif // NDEBUG
}
/*!
* \copydoc FvBaseMultiPhaseProblem::temperature
*/
template
Scalar temperature(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
{
const auto& pos = context.pos(spaceIdx, timeIdx);
if (inHighTemperatureRegion_(pos))
return temperature_ + 100;
return temperature_;
}
/*!
* \copydoc FvBaseMultiPhaseProblem::intrinsicPermeability
*/
template
const DimMatrix& intrinsicPermeability(const Context& context, unsigned spaceIdx,
unsigned timeIdx) const
{
const GlobalPosition& pos = context.pos(spaceIdx, timeIdx);
if (isFineMaterial_(pos))
return fineK_;
return coarseK_;
}
/*!
* \copydoc FvBaseMultiPhaseProblem::porosity
*/
template
Scalar porosity(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
{
const GlobalPosition& pos = context.pos(spaceIdx, timeIdx);
if (isFineMaterial_(pos))
return finePorosity_;
return coarsePorosity_;
}
/*!
* \copydoc FvBaseMultiPhaseProblem::materialLawParams
*/
template
const MaterialLawParams& materialLawParams(const Context& context,
unsigned spaceIdx, unsigned timeIdx) const
{
const GlobalPosition& pos = context.pos(spaceIdx, timeIdx);
if (isFineMaterial_(pos))
return fineMaterialParams_;
return coarseMaterialParams_;
}
/*!
* \brief Return the parameters for the heat storage law of the rock
*
* In this case, we assume the rock-matrix to be granite.
*/
template
const SolidEnergyLawParams&
solidEnergyLawParams(const Context& /*context*/,
unsigned /*spaceIdx*/,
unsigned /*timeIdx*/) const
{ return solidEnergyLawParams_; }
/*!
* \copydoc FvBaseMultiPhaseProblem::thermalConductionParams
*/
template
const ThermalConductionLawParams &
thermalConductionLawParams(const Context& context,
unsigned spaceIdx,
unsigned timeIdx) const
{
const GlobalPosition& pos = context.pos(spaceIdx, timeIdx);
if (isFineMaterial_(pos))
return fineThermalCondParams_;
return coarseThermalCondParams_;
}
//! \}
/*!
* \name Boundary conditions
*/
//! \{
/*!
* \copydoc FvBaseProblem::boundary
*/
template
void boundary(BoundaryRateVector& values, const Context& context,
unsigned spaceIdx, unsigned timeIdx) const
{
const auto& pos = context.pos(spaceIdx, timeIdx);
if (onLeftBoundary_(pos)) {
Opm::CompositionalFluidState fs;
initialFluidState_(fs, context, spaceIdx, timeIdx);
fs.checkDefined();
// impose an freeflow boundary condition
values.setFreeFlow(context, spaceIdx, timeIdx, fs);
}
else if (onInlet_(pos)) {
RateVector massRate(0.0);
massRate[contiCO2EqIdx] = -1e-3; // [kg/(m^3 s)]
using FluidState = Opm::ImmiscibleFluidState;
FluidState fs;
fs.setSaturation(gasPhaseIdx, 1.0);
const auto& pg =
context.intensiveQuantities(spaceIdx, timeIdx).fluidState().pressure(gasPhaseIdx);
fs.setPressure(gasPhaseIdx, Toolbox::value(pg));
fs.setTemperature(temperature(context, spaceIdx, timeIdx));
typename FluidSystem::template ParameterCache paramCache;
paramCache.updatePhase(fs, gasPhaseIdx);
Scalar h = FluidSystem::template enthalpy(fs, paramCache, gasPhaseIdx);
// impose an forced inflow boundary condition for pure CO2
values.setMassRate(massRate);
values.setEnthalpyRate(massRate[contiCO2EqIdx] * h);
}
else
// no flow on top and bottom
values.setNoFlow();
}
// \}
/*!
* \name Volumetric terms
*/
//! \{
/*!
* \copydoc FvBaseProblem::initial
*/
template
void initial(PrimaryVariables& values, const Context& context, unsigned spaceIdx,
unsigned timeIdx) const
{
Opm::CompositionalFluidState fs;
initialFluidState_(fs, context, spaceIdx, timeIdx);
// const auto& matParams = this->materialLawParams(context, spaceIdx,
// timeIdx);
// values.assignMassConservative(fs, matParams, /*inEquilibrium=*/true);
values.assignNaive(fs);
}
/*!
* \copydoc FvBaseProblem::source
*
* For this problem, the source term of all components is 0
* everywhere.
*/
template
void source(RateVector& rate,
const Context& /*context*/,
unsigned /*spaceIdx*/,
unsigned /*timeIdx*/) const
{ rate = Scalar(0.0); }
//! \}
private:
template
void initialFluidState_(FluidState& fs,
const Context& context,
unsigned spaceIdx,
unsigned timeIdx) const
{
const GlobalPosition& pos = context.pos(spaceIdx, timeIdx);
//////
// set temperature
//////
fs.setTemperature(temperature(context, spaceIdx, timeIdx));
//////
// set saturations
//////
fs.setSaturation(FluidSystem::liquidPhaseIdx, 1.0);
fs.setSaturation(FluidSystem::gasPhaseIdx, 0.0);
//////
// set pressures
//////
Scalar densityL = FluidSystem::Brine::liquidDensity(temperature_, Scalar(1e5));
Scalar depth = maxDepth_ - pos[dim - 1];
Scalar pl = 1e5 - densityL * this->gravity()[dim - 1] * depth;
Scalar pC[numPhases];
const auto& matParams = this->materialLawParams(context, spaceIdx, timeIdx);
MaterialLaw::capillaryPressures(pC, matParams, fs);
fs.setPressure(liquidPhaseIdx, pl + (pC[liquidPhaseIdx] - pC[liquidPhaseIdx]));
fs.setPressure(gasPhaseIdx, pl + (pC[gasPhaseIdx] - pC[liquidPhaseIdx]));
//////
// set composition of the liquid phase
//////
fs.setMoleFraction(liquidPhaseIdx, CO2Idx, 0.005);
fs.setMoleFraction(liquidPhaseIdx, BrineIdx,
1.0 - fs.moleFraction(liquidPhaseIdx, CO2Idx));
typename FluidSystem::template ParameterCache paramCache;
using CFRP = Opm::ComputeFromReferencePhase;
CFRP::solve(fs, paramCache,
/*refPhaseIdx=*/liquidPhaseIdx,
/*setViscosity=*/true,
/*setEnthalpy=*/true);
}
bool onLeftBoundary_(const GlobalPosition& pos) const
{ return pos[0] < eps_; }
bool onRightBoundary_(const GlobalPosition& pos) const
{ return pos[0] > this->boundingBoxMax()[0] - eps_; }
bool onInlet_(const GlobalPosition& pos) const
{ return onRightBoundary_(pos) && (5 < pos[1]) && (pos[1] < 15); }
bool inHighTemperatureRegion_(const GlobalPosition& pos) const
{ return (pos[0] > 20) && (pos[0] < 30) && (pos[1] > 5) && (pos[1] < 35); }
void computeThermalCondParams_(ThermalConductionLawParams& params, Scalar poro)
{
Scalar lambdaWater = 0.6;
Scalar lambdaGranite = 2.8;
Scalar lambdaWet = std::pow(lambdaGranite, (1 - poro))
* std::pow(lambdaWater, poro);
Scalar lambdaDry = std::pow(lambdaGranite, (1 - poro));
params.setFullySaturatedLambda(gasPhaseIdx, lambdaDry);
params.setFullySaturatedLambda(liquidPhaseIdx, lambdaWet);
params.setVacuumLambda(lambdaDry);
}
bool isFineMaterial_(const GlobalPosition& pos) const
{ return pos[dim - 1] > fineLayerBottom_; }
DimMatrix fineK_;
DimMatrix coarseK_;
Scalar fineLayerBottom_;
Scalar finePorosity_;
Scalar coarsePorosity_;
MaterialLawParams fineMaterialParams_;
MaterialLawParams coarseMaterialParams_;
ThermalConductionLawParams fineThermalCondParams_;
ThermalConductionLawParams coarseThermalCondParams_;
SolidEnergyLawParams solidEnergyLawParams_;
Scalar temperature_;
Scalar maxDepth_;
Scalar eps_;
unsigned nTemperature_;
unsigned nPressure_;
Scalar pressureLow_, pressureHigh_;
Scalar temperatureLow_, temperatureHigh_;
};
} // namespace Opm
#endif