// -*- 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::EclProblem
*/
#ifndef EWOMS_ECL_PROBLEM_HH
#define EWOMS_ECL_PROBLEM_HH
//#define DISABLE_ALUGRID_SFC_ORDERING 1
//#define EBOS_USE_ALUGRID 1
// make sure that the EBOS_USE_ALUGRID macro. using the preprocessor for this is slightly
// hacky...
#if EBOS_USE_ALUGRID
//#define DISABLE_ALUGRID_SFC_ORDERING 1
#if !HAVE_DUNE_ALUGRID
#warning "ALUGrid was indicated to be used for the ECL black oil simulator, but this "
#warning "requires the presence of dune-alugrid >= 2.4. Falling back to Dune::CpGrid"
#undef EBOS_USE_ALUGRID
#define EBOS_USE_ALUGRID 0
#endif
#else
#define EBOS_USE_ALUGRID 0
#endif
#if EBOS_USE_ALUGRID
#include "eclalugridvanguard.hh"
#elif USE_POLYHEDRALGRID
#include "eclpolyhedralgridvanguard.hh"
#else
#include "eclcpgridvanguard.hh"
#endif
#include "eclwellmanager.hh"
#include "eclequilinitializer.hh"
#include "eclwriter.hh"
#include "ecloutputblackoilmodule.hh"
#include "ecltransmissibility.hh"
#include "eclthresholdpressure.hh"
#include "ecldummygradientcalculator.hh"
#include "eclfluxmodule.hh"
#include "eclbaseaquifermodel.hh"
#include "eclnewtonmethod.hh"
#include "ecltracermodel.hh"
#include "vtkecltracermodule.hh"
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#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 {
template
class EclProblem;
}
namespace Opm::Properties {
namespace TTag {
#if EBOS_USE_ALUGRID
struct EclBaseProblem {
using InheritstFrom = std::tuple;
};
#elif USE_POLYHEDRALGRID
struct EclBaseProblem {
using InheritsFrom = std::tuple;
};
#else
struct EclBaseProblem {
using InheritsFrom = std::tuple;
};
#endif
}
// The class which deals with ECL wells
template
struct EclWellModel {
using type = UndefinedProperty;
};
// Write all solutions for visualization, not just the ones for the
// report steps...
template
struct EnableWriteAllSolutions {
using type = UndefinedProperty;
};
// The number of time steps skipped between writing two consequtive restart files
template
struct RestartWritingInterval {
using type = UndefinedProperty;
};
// Enable partial compensation of systematic mass losses via the source term of the next time
// step
template
struct EclEnableDriftCompensation {
using type = UndefinedProperty;
};
// Enable the additional checks even if compiled in debug mode (i.e., with the NDEBUG
// macro undefined). Next to a slightly better performance, this also eliminates some
// print statements in debug mode.
template
struct EnableDebuggingChecks {
using type = UndefinedProperty;
};
// if thermal flux boundaries are enabled an effort is made to preserve the initial
// thermal gradient specified via the TEMPVD keyword
template
struct EnableThermalFluxBoundaries {
using type = UndefinedProperty;
};
// Specify whether API tracking should be enabled (replaces PVT regions).
// TODO: This is not yet implemented
template
struct EnableApiTracking {
using type = UndefinedProperty;
};
// The class which deals with ECL aquifers
template
struct EclAquiferModel {
using type = UndefinedProperty;
};
// In experimental mode, decides if the aquifer model should be enabled or not
template
struct EclEnableAquifers {
using type = UndefinedProperty;
};
// time stepping parameters
template
struct EclMaxTimeStepSizeAfterWellEvent {
using type = UndefinedProperty;
};
template
struct EclRestartShrinkFactor {
using type = UndefinedProperty;
};
template
struct EclEnableTuning {
using type = UndefinedProperty;
};
template
struct OutputMode {
using type = UndefinedProperty;
};
// Set the problem property
template
struct Problem {
using type = Opm::EclProblem;
};
// Select the element centered finite volume method as spatial discretization
template
struct SpatialDiscretizationSplice {
using type = TTag::EcfvDiscretization;
};
//! for ebos, use automatic differentiation to linearize the system of PDEs
template
struct LocalLinearizerSplice {
using type = TTag::AutoDiffLocalLinearizer;
};
// Set the material law for fluid fluxes
template
struct MaterialLaw
{
private:
using Scalar = GetPropType;
using FluidSystem = GetPropType;
typedef Opm::ThreePhaseMaterialTraits Traits;
public:
typedef Opm::EclMaterialLawManager EclMaterialLawManager;
typedef typename EclMaterialLawManager::MaterialLaw type;
};
// Set the material law for energy storage in rock
template
struct SolidEnergyLaw
{
private:
using Scalar = GetPropType;
using FluidSystem = GetPropType;
public:
typedef Opm::EclThermalLawManager EclThermalLawManager;
typedef typename EclThermalLawManager::SolidEnergyLaw type;
};
// Set the material law for thermal conduction
template
struct ThermalConductionLaw
{
private:
using Scalar = GetPropType;
using FluidSystem = GetPropType;
public:
typedef Opm::EclThermalLawManager EclThermalLawManager;
typedef typename EclThermalLawManager::ThermalConductionLaw type;
};
// ebos can use a slightly faster stencil class because it does not need the normals and
// the integration points of intersections
template
struct Stencil
{
private:
using Scalar = GetPropType;
using GridView = GetPropType;
public:
typedef Opm::EcfvStencil type;
};
// by default use the dummy aquifer "model"
template
struct EclAquiferModel {
using type = Opm::EclBaseAquiferModel;
};
// use the built-in proof of concept well model by default
template
struct EclWellModel {
using type = EclWellManager;
};
// Enable aquifers by default in experimental mode
template
struct EclEnableAquifers {
static constexpr bool value = true;
};
// Enable gravity
template
struct EnableGravity {
static constexpr bool value = true;
};
// only write the solutions for the report steps to disk
template
struct EnableWriteAllSolutions {
static constexpr bool value = false;
};
// disable API tracking
template
struct EnableApiTracking {
static constexpr bool value = false;
};
// The default for the end time of the simulation [s]
//
// By default, stop it after the universe will probably have stopped
// to exist. (the ECL problem will finish the simulation explicitly
// after it simulated the last episode specified in the deck.)
template
struct EndTime {
using type = GetPropType;
static constexpr type value = 1e100;
};
// The default for the initial time step size of the simulation [s].
//
// The chosen value means that the size of the first time step is the
// one of the initial episode (if the length of the initial episode is
// not millions of trillions of years, that is...)
template
struct InitialTimeStepSize {
using type = GetPropType;
static constexpr type value = 3600*24;
};
// the default for the allowed volumetric error for oil per second
template
struct NewtonTolerance {
using type = GetPropType;
static constexpr type value = 1e-2;
};
// the tolerated amount of "incorrect" amount of oil per time step for the complete
// reservoir. this is scaled by the pore volume of the reservoir, i.e., larger reservoirs
// will tolerate larger residuals.
template
struct EclNewtonSumTolerance {
using type = GetPropType;
static constexpr type value = 1e-4;
};
// set the exponent for the volume scaling of the sum tolerance: larger reservoirs can
// tolerate a higher amount of mass lost per time step than smaller ones! since this is
// not linear, we use the cube root of the overall pore volume by default, i.e., the
// value specified by the NewtonSumTolerance parameter is the "incorrect" mass per
// timestep for an reservoir that exhibits 1 m^3 of pore volume. A reservoir with a total
// pore volume of 10^3 m^3 will tolerate 10 times as much.
template
struct EclNewtonSumToleranceExponent {
using type = GetPropType;
static constexpr type value = 1.0/3.0;
};
// set number of Newton iterations where the volumetric residual is considered for
// convergence
template
struct EclNewtonStrictIterations {
static constexpr int value = 8;
};
// set fraction of the pore volume where the volumetric residual may be violated during
// strict Newton iterations
template
struct EclNewtonRelaxedVolumeFraction {
using type = GetPropType;
static constexpr type value = 0.03;
};
// the maximum volumetric error of a cell in the relaxed region
template
struct EclNewtonRelaxedTolerance {
using type = GetPropType;
static constexpr type value = 1e9;
};
// Ignore the maximum error mass for early termination of the newton method.
template
struct NewtonMaxError {
using type = GetPropType;
static constexpr type value = 10e9;
};
// set the maximum number of Newton iterations to 14 because the likelyhood that a time
// step succeeds at more than 14 Newton iteration is rather small
template
struct NewtonMaxIterations {
static constexpr int value = 14;
};
// also, reduce the target for the "optimum" number of Newton iterations to 6. Note that
// this is only relevant if the time step is reduced from the report step size for some
// reason. (because ebos first tries to do a report step using a single time step.)
template
struct NewtonTargetIterations {
static constexpr int value = 6;
};
// Disable the VTK output by default for this problem ...
template
struct EnableVtkOutput {
static constexpr bool value = false;
};
// ... but enable the ECL output by default
template
struct EnableEclOutput {
static constexpr bool value = true;
};
// If available, write the ECL output in a non-blocking manner
template
struct EnableAsyncEclOutput {
static constexpr bool value = true;
};
// By default, use single precision for the ECL formated results
template
struct EclOutputDoublePrecision {
static constexpr bool value = false;
};
// The default location for the ECL output files
template
struct OutputDir {
static constexpr auto value = ".";
};
// the cache for intensive quantities can be used for ECL problems and also yields a
// decent speedup...
template
struct EnableIntensiveQuantityCache {
static constexpr bool value = true;
};
// the cache for the storage term can also be used and also yields a decent speedup
template
struct EnableStorageCache {
static constexpr bool value = true;
};
// Use the "velocity module" which uses the Eclipse "NEWTRAN" transmissibilities
template
struct FluxModule {
using type = Opm::EclTransFluxModule;
};
// Use the dummy gradient calculator in order not to do unnecessary work.
template
struct GradientCalculator {
using type = Opm::EclDummyGradientCalculator;
};
// Use a custom Newton-Raphson method class for ebos in order to attain more
// sophisticated update and error computation mechanisms
template
struct NewtonMethod {
using type = Opm::EclNewtonMethod;
};
// The frequency of writing restart (*.ers) files. This is the number of time steps
// between writing restart files
template
struct RestartWritingInterval {
static constexpr int value = 0xffffff; // disable
};
// Drift compensation is an experimental feature, i.e., systematic errors in the
// conservation quantities are only compensated for
// as default if experimental mode is enabled.
template
struct EclEnableDriftCompensation {
static constexpr bool value = true;
};
// By default, we enable the debugging checks if we're compiled in debug mode
template
struct EnableDebuggingChecks {
static constexpr bool value = true;
};
// store temperature (but do not conserve energy, as long as EnableEnergy is false)
template
struct EnableTemperature {
static constexpr bool value = true;
};
// disable all extensions supported by black oil model. this should not really be
// necessary but it makes things a bit more explicit
template
struct EnablePolymer {
static constexpr bool value = false;
};
template
struct EnableSolvent {
static constexpr bool value = false;
};
template
struct EnableEnergy {
static constexpr bool value = false;
};
template
struct EnableFoam {
static constexpr bool value = false;
};
template
struct EnableExtbo {
static constexpr bool value = false;
};
// disable thermal flux boundaries by default
template
struct EnableThermalFluxBoundaries {
static constexpr bool value = false;
};
template
struct EnableTracerModel {
static constexpr bool value = false;
};
// By default, simulators derived from the EclBaseProblem are production simulators,
// i.e., experimental features must be explicitly enabled at compile time
template
struct EnableExperiments {
static constexpr bool value = false;
};
// set defaults for the time stepping parameters
template
struct EclMaxTimeStepSizeAfterWellEvent {
using type = GetPropType;
static constexpr type value = 3600*24*365.25;
};
template
struct EclRestartShrinkFactor {
using type = GetPropType;
static constexpr type value = 3;
};
template
struct EclEnableTuning {
static constexpr bool value = false;
};
template
struct OutputMode {
static constexpr auto value = "all";
};
} // namespace Opm::Properties
namespace Opm {
/*!
* \ingroup EclBlackOilSimulator
*
* \brief This problem simulates an input file given in the data format used by the
* commercial ECLiPSE simulator.
*/
template
class EclProblem : public GetPropType
{
using ParentType = GetPropType;
using Implementation = GetPropType;
using Scalar = GetPropType;
using GridView = GetPropType;
using Stencil = GetPropType;
using FluidSystem = GetPropType;
using GlobalEqVector = GetPropType;
using EqVector = GetPropType;
// Grid and world dimension
enum { dim = GridView::dimension };
enum { dimWorld = GridView::dimensionworld };
// copy some indices for convenience
enum { numEq = getPropValue() };
enum { numPhases = FluidSystem::numPhases };
enum { numComponents = FluidSystem::numComponents };
enum { enableExperiments = getPropValue() };
enum { enableSolvent = getPropValue() };
enum { enablePolymer = getPropValue() };
enum { enableBrine = getPropValue() };
enum { enablePolymerMolarWeight = getPropValue() };
enum { enableFoam = getPropValue() };
enum { enableExtbo = getPropValue() };
enum { enableTemperature = getPropValue() };
enum { enableEnergy = getPropValue() };
enum { enableThermalFluxBoundaries = getPropValue() };
enum { enableApiTracking = getPropValue() };
enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
enum { gasCompIdx = FluidSystem::gasCompIdx };
enum { oilCompIdx = FluidSystem::oilCompIdx };
enum { waterCompIdx = FluidSystem::waterCompIdx };
using PrimaryVariables = GetPropType;
using RateVector = GetPropType;
using BoundaryRateVector = GetPropType;
using Simulator = GetPropType;
using Element = typename GridView::template Codim<0>::Entity;
using ElementContext = GetPropType;
using EclMaterialLawManager = typename GetProp::EclMaterialLawManager;
using EclThermalLawManager = typename GetProp::EclThermalLawManager;
using MaterialLawParams = typename EclMaterialLawManager::MaterialLawParams;
using SolidEnergyLawParams = typename EclThermalLawManager::SolidEnergyLawParams;
using ThermalConductionLawParams = typename EclThermalLawManager::ThermalConductionLawParams;
using MaterialLaw = GetPropType;
using DofMapper = GetPropType;
using Evaluation = GetPropType;
using Indices = GetPropType;
using IntensiveQuantities = GetPropType;
using EclWellModel = GetPropType;
using EclAquiferModel = GetPropType;
typedef BlackOilSolventModule SolventModule;
typedef BlackOilPolymerModule PolymerModule;
typedef BlackOilFoamModule FoamModule;
typedef BlackOilBrineModule BrineModule;
typedef BlackOilExtboModule ExtboModule;
typedef typename EclEquilInitializer::ScalarFluidState InitialFluidState;
typedef Opm::MathToolbox Toolbox;
typedef Dune::FieldMatrix DimMatrix;
typedef EclWriter EclWriterType;
typedef EclTracerModel TracerModel;
typedef typename GridView::template Codim<0>::Iterator ElementIterator;
typedef Opm::UniformXTabulated2DFunction TabulatedTwoDFunction;
typedef Opm::Tabulated1DFunction TabulatedFunction;
struct RockParams {
Scalar referencePressure;
Scalar compressibility;
};
public:
/*!
* \copydoc FvBaseProblem::registerParameters
*/
static void registerParameters()
{
ParentType::registerParameters();
EclWriterType::registerParameters();
VtkEclTracerModule::registerParameters();
EWOMS_REGISTER_PARAM(TypeTag, bool, EnableWriteAllSolutions,
"Write all solutions to disk instead of only the ones for the "
"report steps");
EWOMS_REGISTER_PARAM(TypeTag, bool, EnableEclOutput,
"Write binary output which is compatible with the commercial "
"Eclipse simulator");
EWOMS_REGISTER_PARAM(TypeTag, bool, EclOutputDoublePrecision,
"Tell the output writer to use double precision. Useful for 'perfect' restarts");
EWOMS_REGISTER_PARAM(TypeTag, unsigned, RestartWritingInterval,
"The frequencies of which time steps are serialized to disk");
EWOMS_REGISTER_PARAM(TypeTag, bool, EnableTracerModel,
"Transport tracers found in the deck.");
EWOMS_REGISTER_PARAM(TypeTag, bool, EclEnableDriftCompensation,
"Enable partial compensation of systematic mass losses via the source term of the next time step");
if (enableExperiments)
EWOMS_REGISTER_PARAM(TypeTag, bool, EclEnableAquifers,
"Enable analytic and numeric aquifer models");
EWOMS_REGISTER_PARAM(TypeTag, Scalar, EclMaxTimeStepSizeAfterWellEvent,
"Maximum time step size after an well event");
EWOMS_REGISTER_PARAM(TypeTag, Scalar, EclRestartShrinkFactor,
"Factor by which the time step is reduced after convergence failure");
EWOMS_REGISTER_PARAM(TypeTag, bool, EclEnableTuning,
"Honor some aspects of the TUNING keyword from the ECL deck.");
EWOMS_REGISTER_PARAM(TypeTag, std::string, OutputMode,
"Specify which messages are going to be printed. Valid values are: none, log, all (default)");
}
/*!
* \copydoc FvBaseProblem::prepareOutputDir
*/
std::string prepareOutputDir() const
{ return this->simulator().vanguard().eclState().getIOConfig().getOutputDir(); }
/*!
* \copydoc FvBaseProblem::handlePositionalParameter
*/
static int handlePositionalParameter(std::set& seenParams,
std::string& errorMsg,
int argc OPM_UNUSED,
const char** argv,
int paramIdx,
int posParamIdx OPM_UNUSED)
{
using ParamsMeta = GetProp;
Dune::ParameterTree& tree = ParamsMeta::tree();
std::string param = argv[paramIdx];
size_t i = param.find('=');
if (i != std::string::npos) {
std::string oldParamName = param.substr(0, i);
std::string oldParamValue = param.substr(i+1);
std::string newParamName = "--" + oldParamName;
for (size_t j = 0; j < newParamName.size(); ++j)
if (newParamName[j] == '_')
newParamName[j] = '-';
errorMsg =
"The old syntax to specify parameters on the command line is no longer supported: "
"Try replacing '"+oldParamName+"="+oldParamValue+"' with "+
"'"+newParamName+"="+oldParamValue+"'!";
return 0;
}
if (seenParams.count("EclDeckFileName") > 0) {
errorMsg =
"Parameter 'EclDeckFileName' specified multiple times"
" as a command line parameter";
return 0;
}
tree["EclDeckFileName"] = argv[paramIdx];
seenParams.insert("EclDeckFileName");
return 1;
}
/*!
* \copydoc FvBaseProblem::helpPreamble
*/
static std::string helpPreamble(int argc OPM_UNUSED,
const char **argv)
{
std::string desc = Implementation::briefDescription();
if (!desc.empty())
desc = desc + "\n";
return
"Usage: "+std::string(argv[0]) + " [OPTIONS] [ECL_DECK_FILENAME]\n"
+ desc;
}
/*!
* \copydoc FvBaseProblem::briefDescription
*/
static std::string briefDescription()
{
if (briefDescription_.empty())
return
"The Ecl-deck Black-Oil reservoir Simulator (ebos); a hydrocarbon "
"reservoir simulation program that processes ECL-formatted input "
"files that is part of the Open Porous Media project "
"(https://opm-project.org).\n"
"\n"
"THE GOAL OF THE `ebos` SIMULATOR IS TO CATER FOR THE NEEDS OF "
"DEVELOPMENT AND RESEARCH. No guarantees are made for production use!";
else
return briefDescription_;
}
/*!
* \brief Specifies the string returned by briefDescription()
*
* This string appears in the usage message.
*/
static void setBriefDescription(const std::string& msg)
{ briefDescription_ = msg; }
/*!
* \copydoc Doxygen::defaultProblemConstructor
*/
EclProblem(Simulator& simulator)
: ParentType(simulator)
, transmissibilities_(simulator.vanguard())
, thresholdPressures_(simulator)
, wellModel_(simulator)
, aquiferModel_(simulator)
, pffDofData_(simulator.gridView(), this->elementMapper())
, tracerModel_(simulator)
{
this->model().addOutputModule(new VtkEclTracerModule(simulator));
// Tell the black-oil extensions to initialize their internal data structures
const auto& vanguard = simulator.vanguard();
SolventModule::initFromState(vanguard.eclState(), vanguard.schedule());
PolymerModule::initFromState(vanguard.eclState());
FoamModule::initFromState(vanguard.eclState());
BrineModule::initFromState(vanguard.eclState());
ExtboModule::initFromState(vanguard.eclState());
// create the ECL writer
eclWriter_.reset(new EclWriterType(simulator));
enableDriftCompensation_ = EWOMS_GET_PARAM(TypeTag, bool, EclEnableDriftCompensation);
enableEclOutput_ = EWOMS_GET_PARAM(TypeTag, bool, EnableEclOutput);
if (enableExperiments)
enableAquifers_ = EWOMS_GET_PARAM(TypeTag, bool, EclEnableAquifers);
else
enableAquifers_ = true;
enableTuning_ = EWOMS_GET_PARAM(TypeTag, bool, EclEnableTuning);
initialTimeStepSize_ = EWOMS_GET_PARAM(TypeTag, Scalar, InitialTimeStepSize);
minTimeStepSize_ = EWOMS_GET_PARAM(TypeTag, Scalar, MinTimeStepSize);
maxTimeStepSize_ = EWOMS_GET_PARAM(TypeTag, Scalar, MaxTimeStepSize);
maxTimeStepAfterWellEvent_ = EWOMS_GET_PARAM(TypeTag, Scalar, EclMaxTimeStepSizeAfterWellEvent);
restartShrinkFactor_ = EWOMS_GET_PARAM(TypeTag, Scalar, EclRestartShrinkFactor);
maxFails_ = EWOMS_GET_PARAM(TypeTag, unsigned, MaxTimeStepDivisions);
}
/*!
* \copydoc FvBaseProblem::finishInit
*/
void finishInit()
{
ParentType::finishInit();
auto& simulator = this->simulator();
const auto& eclState = simulator.vanguard().eclState();
const auto& schedule = simulator.vanguard().schedule();
const auto& timeMap = schedule.getTimeMap();
// Set the start time of the simulation
simulator.setStartTime(timeMap.getStartTime(/*reportStepIdx=*/0));
simulator.setEndTime(timeMap.getTotalTime());
// We want the episode index to be the same as the report step index to make
// things simpler, so we have to set the episode index to -1 because it is
// incremented by endEpisode(). The size of the initial time step and
// length of the initial episode is set to zero for the same reason.
simulator.setEpisodeIndex(-1);
simulator.setEpisodeLength(0.0);
// the "NOGRAV" keyword from Frontsim or setting the EnableGravity to false
// disables gravity, else the standard value of the gravity constant at sea level
// on earth is used
this->gravity_ = 0.0;
if (EWOMS_GET_PARAM(TypeTag, bool, EnableGravity))
this->gravity_[dim - 1] = 9.80665;
if (!eclState.getInitConfig().hasGravity())
this->gravity_[dim - 1] = 0.0;
if (enableTuning_) {
// if support for the TUNING keyword is enabled, we get the initial time
// steping parameters from it instead of from command line parameters
const auto& tuning = schedule.getTuning(0);
initialTimeStepSize_ = tuning.TSINIT;
maxTimeStepAfterWellEvent_ = tuning.TMAXWC;
maxTimeStepSize_ = tuning.TSMAXZ;
restartShrinkFactor_ = 1./tuning.TSFCNV;
minTimeStepSize_ = tuning.TSMINZ;
}
initFluidSystem_();
// deal with DRSDT
unsigned ntpvt = eclState.runspec().tabdims().getNumPVTTables();
size_t numDof = this->model().numGridDof();
if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
maxDRs_.resize(ntpvt, 1e30);
dRsDtOnlyFreeGas_.resize(ntpvt, false);
lastRs_.resize(numDof, 0.0);
maxDRv_.resize(ntpvt, 1e30);
lastRv_.resize(numDof, 0.0);
maxOilSaturation_.resize(numDof, 0.0);
}
readRockParameters_();
readMaterialParameters_();
readThermalParameters_();
transmissibilities_.finishInit();
const auto& initconfig = eclState.getInitConfig();
if (initconfig.restartRequested())
readEclRestartSolution_();
else
readInitialCondition_();
updatePffDofData_();
if (getPropValue()) {
const auto& vanguard = this->simulator().vanguard();
const auto& gridView = vanguard.gridView();
int numElements = gridView.size(/*codim=*/0);
maxPolymerAdsorption_.resize(numElements, 0.0);
}
tracerModel_.init();
readBoundaryConditions_();
if (enableDriftCompensation_) {
drift_.resize(numDof);
drift_ = 0.0;
}
if (enableExperiments)
{
int success = 1;
const auto& cc = simulator.vanguard().grid().comm();
try
{
checkDeckCompatibility_();
}
catch(const std::exception& e)
{
success = 0;
success = cc.min(success);
throw;
}
success = cc.min(success);
if (!success)
{
throw std::runtime_error("Checking deck compatibility failed");
}
}
// write the static output files (EGRID, INIT, SMSPEC, etc.)
if (enableEclOutput_)
eclWriter_->writeInit();
simulator.vanguard().releaseGlobalTransmissibilities();
// after finishing the initialization and writing the initial solution, we move
// to the first "real" episode/report step
// for restart the episode index and start is already set
if (!initconfig.restartRequested()) {
simulator.startNextEpisode(timeMap.getTimeStepLength(0));
simulator.setEpisodeIndex(0);
}
}
void prefetch(const Element& elem) const
{ pffDofData_.prefetch(elem); }
/*!
* \brief This method restores the complete state of the problem and its sub-objects
* from disk.
*
* The serialization format used by this method is ad-hoc. It is the inverse of the
* serialize() method.
*
* \tparam Restarter The deserializer type
*
* \param res The deserializer object
*/
template
void deserialize(Restarter& res)
{
// reload the current episode/report step from the deck
beginEpisode();
// deserialize the wells
wellModel_.deserialize(res);
if (enableAquifers_)
// deserialize the aquifer
aquiferModel_.deserialize(res);
}
/*!
* \brief This method writes the complete state of the problem and its subobjects to
* disk.
*
* The file format used here is ad-hoc.
*/
template
void serialize(Restarter& res)
{
wellModel_.serialize(res);
if (enableAquifers_)
aquiferModel_.serialize(res);
}
/*!
* \brief Called by the simulator before an episode begins.
*/
void beginEpisode()
{
// Proceed to the next report step
auto& simulator = this->simulator();
int episodeIdx = simulator.episodeIndex();
auto& eclState = simulator.vanguard().eclState();
const auto& schedule = simulator.vanguard().schedule();
const auto& events = schedule.getEvents();
const auto& timeMap = schedule.getTimeMap();
if (episodeIdx >= 0 && events.hasEvent(Opm::ScheduleEvents::GEO_MODIFIER, episodeIdx)) {
// bring the contents of the keywords to the current state of the SCHEDULE
// section.
//
// TODO (?): make grid topology changes possible (depending on what exactly
// has changed, the grid may need be re-created which has some serious
// implications on e.g., the solution of the simulation.)
const auto& miniDeck = schedule.getModifierDeck(episodeIdx);
eclState.applyModifierDeck(miniDeck);
// re-compute all quantities which may possibly be affected.
transmissibilities_.update(true);
referencePorosity_[1] = referencePorosity_[0];
updateReferencePorosity_();
updatePffDofData_();
}
if (enableExperiments && this->gridView().comm().rank() == 0 && episodeIdx >= 0) {
// print some useful information in experimental mode. (the production
// simulator does this externally.)
std::ostringstream ss;
boost::posix_time::time_facet* facet = new boost::posix_time::time_facet("%d-%b-%Y");
boost::posix_time::ptime curDateTime =
boost::posix_time::from_time_t(timeMap.getStartTime(episodeIdx));
ss.imbue(std::locale(std::locale::classic(), facet));
ss << "Report step " << episodeIdx + 1
<< "/" << timeMap.numTimesteps()
<< " at day " << timeMap.getTimePassedUntil(episodeIdx)/(24*3600)
<< "/" << timeMap.getTotalTime()/(24*3600)
<< ", date = " << curDateTime.date()
<< "\n ";
OpmLog::info(ss.str());
}
// react to TUNING changes
bool tuningEvent = false;
if (episodeIdx > 0 && enableTuning_ && events.hasEvent(Opm::ScheduleEvents::TUNING_CHANGE, episodeIdx))
{
const auto& tuning = schedule.getTuning(episodeIdx);
initialTimeStepSize_ = tuning.TSINIT;
maxTimeStepAfterWellEvent_ = tuning.TMAXWC;
maxTimeStepSize_ = tuning.TSMAXZ;
restartShrinkFactor_ = 1./tuning.TSFCNV;
minTimeStepSize_ = tuning.TSMINZ;
tuningEvent = true;
}
// set up the wells for the next episode.
wellModel_.beginEpisode();
// set up the aquifers for the next episode.
if (enableAquifers_)
// set up the aquifers for the next episode.
aquiferModel_.beginEpisode();
// set the size of the initial time step of the episode
Scalar dt = limitNextTimeStepSize_(simulator.episodeLength());
if (episodeIdx == 0 || tuningEvent)
// allow the size of the initial time step to be set via an external parameter
// if TUNING is enabled, also limit the time step size after a tuning event to TSINIT
dt = std::min(dt, initialTimeStepSize_);
simulator.setTimeStepSize(dt);
}
/*!
* \brief Called by the simulator before each time integration.
*/
void beginTimeStep()
{
const auto& simulator = this->simulator();
int episodeIdx = simulator.episodeIndex();
if (enableExperiments && this->gridView().comm().rank() == 0 && episodeIdx >= 0) {
std::ostringstream ss;
boost::posix_time::time_facet* facet = new boost::posix_time::time_facet("%d-%b-%Y");
boost::posix_time::ptime date = boost::posix_time::from_time_t((this->simulator().startTime())) + boost::posix_time::milliseconds(static_cast(this->simulator().time() / Opm::prefix::milli));
ss.imbue(std::locale(std::locale::classic(), facet));
ss <<"\nTime step " << this->simulator().timeStepIndex() << ", stepsize "
<< unit::convert::to(this->simulator().timeStepSize(), unit::day) << " days,"
<< " at day " << (double)unit::convert::to(this->simulator().time(), unit::day)
<< "/" << (double)unit::convert::to(this->simulator().endTime(), unit::day)
<< ", date = " << date;
OpmLog::info(ss.str());
}
// update explicit quantities between timesteps.
const auto& oilVaporizationControl = simulator.vanguard().schedule().getOilVaporizationProperties(episodeIdx);
if (drsdtActive_())
// DRSDT is enabled
for (size_t pvtRegionIdx = 0; pvtRegionIdx < maxDRs_.size(); ++pvtRegionIdx)
maxDRs_[pvtRegionIdx] = oilVaporizationControl.getMaxDRSDT(pvtRegionIdx)*simulator.timeStepSize();
if (drvdtActive_())
// DRVDT is enabled
for (size_t pvtRegionIdx = 0; pvtRegionIdx < maxDRv_.size(); ++pvtRegionIdx)
maxDRv_[pvtRegionIdx] = oilVaporizationControl.getMaxDRVDT(pvtRegionIdx)*this->simulator().timeStepSize();
// update maximum water saturation and minimum pressure
// used when ROCKCOMP is activated
const bool invalidateFromMaxWaterSat = updateMaxWaterSaturation_();
const bool invalidateFromMinPressure = updateMinPressure_();
// update hysteresis and max oil saturation used in vappars
const bool invalidateFromHyst = updateHysteresis_();
const bool invalidateFromMaxOilSat = updateMaxOilSaturation_();
// the derivatives may have change
bool invalidateIntensiveQuantities = invalidateFromMaxWaterSat || invalidateFromMinPressure || invalidateFromHyst || invalidateFromMaxOilSat;
if (invalidateIntensiveQuantities)
this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
if (getPropValue())
updateMaxPolymerAdsorption_();
wellModel_.beginTimeStep();
if (enableAquifers_)
aquiferModel_.beginTimeStep();
tracerModel_.beginTimeStep();
}
/*!
* \brief Return if the storage term of the first iteration is identical to the storage
* term for the solution of the previous time step.
*
* For quite technical reasons, the storage term cannot be recycled if either DRSDT
* or DRVDT are active in ebos. Nor if the porosity is changes between timesteps
* using a pore volume multiplier (i.e., poreVolumeMultiplier() != 1.0)
*/
bool recycleFirstIterationStorage() const
{ return !drsdtActive_() && !drvdtActive_() && rockCompPoroMultWc_.empty() && rockCompPoroMult_.empty(); }
/*!
* \brief Called by the simulator before each Newton-Raphson iteration.
*/
void beginIteration()
{
wellModel_.beginIteration();
if (enableAquifers_)
aquiferModel_.beginIteration();
}
/*!
* \brief Called by the simulator after each Newton-Raphson iteration.
*/
void endIteration()
{
wellModel_.endIteration();
if (enableAquifers_)
aquiferModel_.endIteration();
}
/*!
* \brief Called by the simulator after each time integration.
*/
void endTimeStep()
{
#ifndef NDEBUG
if (getPropValue()) {
// in debug mode, we don't care about performance, so we check if the model does
// the right thing (i.e., the mass change inside the whole reservoir must be
// equivalent to the fluxes over the grid's boundaries plus the source rates
// specified by the problem)
int rank = this->simulator().gridView().comm().rank();
if (rank == 0)
std::cout << "checking conservativeness of solution\n";
this->model().checkConservativeness(/*tolerance=*/-1, /*verbose=*/true);
if (rank == 0)
std::cout << "solution is sufficiently conservative\n";
}
#endif // NDEBUG
auto& simulator = this->simulator();
wellModel_.endTimeStep();
if (enableAquifers_)
aquiferModel_.endTimeStep();
tracerModel_.endTimeStep();
// deal with DRSDT and DRVDT
updateCompositionChangeLimits_();
if (enableDriftCompensation_) {
const auto& residual = this->model().linearizer().residual();
for (unsigned globalDofIdx = 0; globalDofIdx < residual.size(); globalDofIdx ++) {
drift_[globalDofIdx] = residual[globalDofIdx];
drift_[globalDofIdx] *= simulator.timeStepSize();
if (getPropValue())
drift_[globalDofIdx] *= this->model().dofTotalVolume(globalDofIdx);
}
}
bool isSubStep = !EWOMS_GET_PARAM(TypeTag, bool, EnableWriteAllSolutions) && !this->simulator().episodeWillBeOver();
eclWriter_->evalSummaryState(isSubStep);
auto& schedule = simulator.vanguard().schedule();
auto& ecl_state = simulator.vanguard().eclState();
int episodeIdx = simulator.episodeIndex();
this->applyActions(episodeIdx,
simulator.time() + simulator.timeStepSize(),
ecl_state,
schedule,
simulator.vanguard().actionState(),
simulator.vanguard().summaryState());
}
/*!
* \brief Called by the simulator after the end of an episode.
*/
void endEpisode()
{
auto& simulator = this->simulator();
auto& schedule = simulator.vanguard().schedule();
const auto& timeMap = schedule.getTimeMap();
int episodeIdx = simulator.episodeIndex();
// check if we're finished ...
if (episodeIdx + 1 >= static_cast(timeMap.numTimesteps())) {
simulator.setFinished(true);
return;
}
// .. if we're not yet done, start the next episode (report step)
simulator.startNextEpisode(timeMap.getTimeStepLength(episodeIdx + 1));
}
/*!
* \brief Always returns true. The ecl output writer takes care of the rest
*/
bool shouldWriteOutput() const
{ return true; }
/*!
* \brief Returns true if an eWoms restart file should be written to disk.
*
* The EclProblem does not write any restart files using the ad-hoc format, only ones
* using the ECL format.
*/
bool shouldWriteRestartFile() const
{ return false; }
/*!
* \brief Write the requested quantities of the current solution into the output
* files.
*/
void writeOutput(bool verbose = true)
{
// use the generic code to prepare the output fields and to
// write the desired VTK files.
ParentType::writeOutput(verbose);
bool isSubStep = !EWOMS_GET_PARAM(TypeTag, bool, EnableWriteAllSolutions) && !this->simulator().episodeWillBeOver();
if (enableEclOutput_)
eclWriter_->writeOutput(isSubStep);
}
void finalizeOutput() {
// this will write all pending output to disk
// to avoid corruption of output files
eclWriter_.reset();
}
void applyActions(int reportStep,
double sim_time,
Opm::EclipseState& ecl_state,
Opm::Schedule& schedule,
Opm::Action::State& actionState,
Opm::SummaryState& summaryState) {
const auto& actions = schedule.actions(reportStep);
if (actions.empty())
return;
Opm::Action::Context context( summaryState, schedule.getWListManager(reportStep) );
auto now = Opm::TimeStampUTC( schedule.getStartTime() ) + std::chrono::duration(sim_time);
std::string ts;
{
std::ostringstream os;
os << std::setw(4) << std::to_string(now.year()) << '/'
<< std::setw(2) << std::setfill('0') << std::to_string(now.month()) << '/'
<< std::setw(2) << std::setfill('0') << std::to_string(now.day()) << " report:" << std::to_string(reportStep);
ts = os.str();
}
for (const auto& pyaction : actions.pending_python()) {
pyaction->run(ecl_state, schedule, reportStep, summaryState);
}
auto simTime = schedule.simTime(reportStep);
for (const auto& action : actions.pending(actionState, simTime)) {
auto actionResult = action->eval(context);
if (actionResult) {
std::string wells_string;
const auto& matching_wells = actionResult.wells();
if (matching_wells.size() > 0) {
for (std::size_t iw = 0; iw < matching_wells.size() - 1; iw++)
wells_string += matching_wells[iw] + ", ";
wells_string += matching_wells.back();
}
std::string msg = "The action: " + action->name() + " evaluated to true at " + ts + " wells: " + wells_string;
Opm::OpmLog::info(msg);
schedule.applyAction(reportStep, *action, actionResult);
actionState.add_run(*action, simTime);
} else {
std::string msg = "The action: " + action->name() + " evaluated to false at " + ts;
Opm::OpmLog::info(msg);
}
}
}
/*!
* \copydoc FvBaseMultiPhaseProblem::intrinsicPermeability
*/
template
const DimMatrix& intrinsicPermeability(const Context& context,
unsigned spaceIdx,
unsigned timeIdx) const
{
unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
return transmissibilities_.permeability(globalSpaceIdx);
}
/*!
* \brief This method returns the intrinsic permeability tensor
* given a global element index.
*
* Its main (only?) usage is the ECL transmissibility calculation code...
*/
const DimMatrix& intrinsicPermeability(unsigned globalElemIdx) const
{ return transmissibilities_.permeability(globalElemIdx); }
/*!
* \copydoc EclTransmissiblity::transmissibility
*/
template
Scalar transmissibility(const Context& context,
unsigned OPM_OPTIM_UNUSED fromDofLocalIdx,
unsigned toDofLocalIdx) const
{
assert(fromDofLocalIdx == 0);
return pffDofData_.get(context.element(), toDofLocalIdx).transmissibility;
}
/*!
* \copydoc EclTransmissiblity::transmissibilityBoundary
*/
template
Scalar transmissibilityBoundary(const Context& elemCtx,
unsigned boundaryFaceIdx) const
{
unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
return transmissibilities_.transmissibilityBoundary(elemIdx, boundaryFaceIdx);
}
/*!
* \copydoc EclTransmissiblity::thermalHalfTransmissibility
*/
template
Scalar thermalHalfTransmissibility(const Context& context,
unsigned faceIdx,
unsigned timeIdx) const
{
const auto& face = context.stencil(timeIdx).interiorFace(faceIdx);
unsigned toDofLocalIdx = face.exteriorIndex();
return *pffDofData_.get(context.element(), toDofLocalIdx).thermalHalfTrans;
}
/*!
* \copydoc EclTransmissiblity::thermalHalfTransmissibility
*/
template
Scalar thermalHalfTransmissibilityBoundary(const Context& elemCtx,
unsigned boundaryFaceIdx) const
{
unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
return transmissibilities_.thermalHalfTransBoundary(elemIdx, boundaryFaceIdx);
}
/*!
* \brief Return a reference to the object that handles the "raw" transmissibilities.
*/
const EclTransmissibility& eclTransmissibilities() const
{ return transmissibilities_; }
/*!
* \copydoc BlackOilBaseProblem::thresholdPressure
*/
Scalar thresholdPressure(unsigned elem1Idx, unsigned elem2Idx) const
{ return thresholdPressures_.thresholdPressure(elem1Idx, elem2Idx); }
const EclThresholdPressure& thresholdPressure() const
{ return thresholdPressures_; }
EclThresholdPressure& thresholdPressure()
{ return thresholdPressures_; }
const EclTracerModel& tracerModel() const
{ return tracerModel_; }
/*!
* \copydoc FvBaseMultiPhaseProblem::porosity
*
* For the EclProblem, this method is identical to referencePorosity(). The intensive
* quantities object may apply various multipliers (e.g. ones which model rock
* compressibility and water induced rock compaction) to it which depend on the
* current physical conditions.
*/
template