mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
ebos: add an API for aquifiers
this is similar to the mechanism for well models: the API is defined here, but can be overloaded by downstream modules. In contrast to the well model where the default class follows a simplistic approach the default class for aquifiers does nothing, i.e., the actual implemenentation must be provided by downstream.
This commit is contained in:
parent
cee74bb2b5
commit
23670c6797
128
ebos/eclbaseaquifermodel.hh
Normal file
128
ebos/eclbaseaquifermodel.hh
Normal file
@ -0,0 +1,128 @@
|
||||
// -*- 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
|
||||
* \copydoc Ewoms::EclBaseAquiferModel
|
||||
*/
|
||||
#ifndef EWOMS_ECL_BASE_AQUIFER_MODEL_HH
|
||||
#define EWOMS_ECL_BASE_AQUIFER_MODEL_HH
|
||||
|
||||
#include <ewoms/common/propertysystem.hh>
|
||||
|
||||
namespace Ewoms {
|
||||
|
||||
/*!
|
||||
* \ingroup EclBaseAquiferModel
|
||||
*
|
||||
* \brief The base class which specifies the API of aquifer models.
|
||||
*
|
||||
* This class only provides the API for the actual aquifer model, it does not do
|
||||
* anything on its own.
|
||||
*/
|
||||
template <class TypeTag>
|
||||
class EclBaseAquiferModel
|
||||
{
|
||||
typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
|
||||
typedef typename GET_PROP_TYPE(TypeTag, RateVector) RateVector;
|
||||
|
||||
public:
|
||||
EclBaseAquiferModel(Simulator& simulator)
|
||||
: simulator_(simulator)
|
||||
{}
|
||||
|
||||
/*!
|
||||
* \brief This method is called when a new episode (report step) starts.
|
||||
*/
|
||||
void beginEpisode()
|
||||
{ }
|
||||
|
||||
/*!
|
||||
* \brief This method is called when a new time step (substep) starts.
|
||||
*/
|
||||
void beginTimeStep()
|
||||
{ }
|
||||
|
||||
/*!
|
||||
* \brief This method is called before each Newton-Raphson iteration.
|
||||
*/
|
||||
void beginIteration()
|
||||
{ }
|
||||
|
||||
/*!
|
||||
* \brief Add the water which enters or leaves the reservoir due to aquifiers.
|
||||
*/
|
||||
template <class Context>
|
||||
void addToSource(RateVector& rate,
|
||||
const Context& context,
|
||||
unsigned spaceIdx,
|
||||
unsigned timeIdx) const
|
||||
{ }
|
||||
|
||||
/*!
|
||||
* \brief This method is called after each Newton-Raphson successful iteration.
|
||||
*
|
||||
* I.e., no exceptions were thrown during the linearization and linear solution
|
||||
* procedures.
|
||||
*/
|
||||
void endIteration()
|
||||
{ }
|
||||
|
||||
/*!
|
||||
* \brief This method is called after each successful time step (substep).
|
||||
*
|
||||
* I.e., all iterations of the time step were successful and the Newton-Raphson
|
||||
* algorithm converged.
|
||||
*/
|
||||
void endTimeStep()
|
||||
{ }
|
||||
|
||||
/*!
|
||||
* \brief This method is called once an episode (report step) has been finished
|
||||
* successfully.
|
||||
*/
|
||||
void endEpisode()
|
||||
{ }
|
||||
|
||||
/*!
|
||||
* \brief Write the internal state of the aquifer model to disk using an ad-hoc file
|
||||
* format.
|
||||
*/
|
||||
template <class Restarter>
|
||||
void serialize(Restarter& res)
|
||||
{ }
|
||||
|
||||
/*!
|
||||
* \brief Load the internal state of the aquifer model to disk using an ad-hoc file
|
||||
* format.
|
||||
*/
|
||||
template <class Restarter>
|
||||
void deserialize(Restarter& res)
|
||||
{ }
|
||||
|
||||
protected:
|
||||
Simulator& simulator_;
|
||||
};
|
||||
|
||||
} // namespace Ewoms
|
||||
|
||||
#endif
|
@ -59,6 +59,7 @@
|
||||
#include "eclthresholdpressure.hh"
|
||||
#include "ecldummygradientcalculator.hh"
|
||||
#include "eclfluxmodule.hh"
|
||||
#include "eclbaseaquifermodel.hh"
|
||||
|
||||
#include <ewoms/common/pffgridvector.hh>
|
||||
#include <ewoms/models/blackoil/blackoilmodel.hh>
|
||||
@ -133,6 +134,9 @@ NEW_PROP_TAG(EnableDebuggingChecks);
|
||||
// thermal gradient specified via the TEMPVD keyword
|
||||
NEW_PROP_TAG(EnableThermalFluxBoundaries);
|
||||
|
||||
// The class which deals with ECL aquifers
|
||||
NEW_PROP_TAG(EclAquiferModel);
|
||||
|
||||
// Set the problem property
|
||||
SET_TYPE_PROP(EclBaseProblem, Problem, Ewoms::EclProblem<TypeTag>);
|
||||
|
||||
@ -201,6 +205,9 @@ public:
|
||||
/*needNormal=*/false> type;
|
||||
};
|
||||
|
||||
// by default use the dummy aquifer "model"
|
||||
SET_TYPE_PROP(EclBaseProblem, EclAquiferModel, Ewoms::EclBaseAquiferModel<TypeTag>);
|
||||
|
||||
// use the built-in proof of concept well model by default
|
||||
SET_TYPE_PROP(EclBaseProblem, EclWellModel, EclWellManager<TypeTag>);
|
||||
|
||||
@ -350,6 +357,7 @@ class EclProblem : public GET_PROP_TYPE(TypeTag, BaseProblem)
|
||||
typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
|
||||
typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
|
||||
typedef typename GET_PROP_TYPE(TypeTag, EclWellModel) EclWellModel;
|
||||
typedef typename GET_PROP_TYPE(TypeTag, EclAquiferModel) EclAquiferModel;
|
||||
|
||||
typedef BlackOilSolventModule<TypeTag> SolventModule;
|
||||
typedef BlackOilPolymerModule<TypeTag> PolymerModule;
|
||||
@ -485,6 +493,7 @@ public:
|
||||
, transmissibilities_(simulator.vanguard())
|
||||
, thresholdPressures_(simulator)
|
||||
, wellModel_(simulator)
|
||||
, aquiferModel_(simulator)
|
||||
, pffDofData_(simulator.gridView(), this->elementMapper())
|
||||
{
|
||||
// Tell the black-oil extensions to initialize their internal data structures
|
||||
@ -609,10 +618,11 @@ public:
|
||||
{ pffDofData_.prefetch(elem); }
|
||||
|
||||
/*!
|
||||
* \brief This method restores the complete state of the well
|
||||
* \brief This method restores the complete state of the problem and its sub-objects
|
||||
* from disk.
|
||||
*
|
||||
* It is the inverse of the serialize() method.
|
||||
* The serialization format used by this method is ad-hoc. It is the inverse of the
|
||||
* serialize() method.
|
||||
*
|
||||
* \tparam Restarter The deserializer type
|
||||
*
|
||||
@ -626,15 +636,23 @@ public:
|
||||
|
||||
// deserialize the wells
|
||||
wellModel_.deserialize(res);
|
||||
|
||||
// deserialize the aquifer
|
||||
aquiferModel_.deserialize(res);
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief This method writes the complete state of the well
|
||||
* to the harddisk.
|
||||
* \brief This method writes the complete state of the problem and its subobjects to
|
||||
* disk.
|
||||
*
|
||||
* The file format used here is ad-hoc.
|
||||
*/
|
||||
template <class Restarter>
|
||||
void serialize(Restarter& res)
|
||||
{ wellModel_.serialize(res); }
|
||||
{
|
||||
wellModel_.serialize(res);
|
||||
aquiferModel_.serialize(res);
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief Called by the simulator before an episode begins.
|
||||
@ -702,11 +720,15 @@ public:
|
||||
updateMaxPolymerAdsorption_();
|
||||
|
||||
if (!GET_PROP_VALUE(TypeTag, DisableWells))
|
||||
// set up the wells
|
||||
// set up the wells for the next episode.
|
||||
//
|
||||
// TODO: the first two arguments seem to be unnecessary
|
||||
wellModel_.beginEpisode(this->simulator().vanguard().eclState(),
|
||||
this->simulator().vanguard().schedule(),
|
||||
isOnRestart);
|
||||
|
||||
aquiferModel_.beginEpisode();
|
||||
|
||||
if (doInvalidate)
|
||||
this->model().invalidateIntensiveQuantitiesCache(/*timeIdx=*/0);
|
||||
}
|
||||
@ -724,9 +746,10 @@ public:
|
||||
// DRVDT is enabled
|
||||
maxDRv_ = maxDRvDt_*this->simulator().timeStepSize();
|
||||
|
||||
if (!GET_PROP_VALUE(TypeTag, DisableWells)) {
|
||||
if (!GET_PROP_VALUE(TypeTag, DisableWells))
|
||||
wellModel_.beginTimeStep();
|
||||
}
|
||||
|
||||
aquiferModel_.beginTimeStep();
|
||||
}
|
||||
|
||||
/*!
|
||||
@ -736,6 +759,8 @@ public:
|
||||
{
|
||||
if (!GET_PROP_VALUE(TypeTag, DisableWells))
|
||||
wellModel_.beginIteration();
|
||||
|
||||
aquiferModel_.beginIteration();
|
||||
}
|
||||
|
||||
/*!
|
||||
@ -745,6 +770,8 @@ public:
|
||||
{
|
||||
if (!GET_PROP_VALUE(TypeTag, DisableWells))
|
||||
wellModel_.endIteration();
|
||||
|
||||
aquiferModel_.endIteration();
|
||||
}
|
||||
|
||||
/*!
|
||||
@ -762,9 +789,10 @@ public:
|
||||
}
|
||||
#endif // NDEBUG
|
||||
|
||||
if (!GET_PROP_VALUE(TypeTag, DisableWells)) {
|
||||
if (!GET_PROP_VALUE(TypeTag, DisableWells))
|
||||
wellModel_.endTimeStep();
|
||||
}
|
||||
|
||||
aquiferModel_.endTimeStep();
|
||||
|
||||
// we no longer need the initial soluiton
|
||||
if (this->simulator().episodeIndex() == 0 && !initialFluidStates_.empty()) {
|
||||
@ -782,7 +810,6 @@ public:
|
||||
initialFluidStates_.clear();
|
||||
}
|
||||
|
||||
|
||||
updateCompositionChangeLimits_();
|
||||
}
|
||||
|
||||
@ -802,6 +829,11 @@ public:
|
||||
simulator.setFinished(true);
|
||||
return;
|
||||
}
|
||||
|
||||
if (!GET_PROP_VALUE(TypeTag, DisableWells))
|
||||
wellModel_.endEpisode();
|
||||
|
||||
aquiferModel_.endEpisode();
|
||||
}
|
||||
|
||||
/*!
|
||||
@ -1314,6 +1346,8 @@ public:
|
||||
assert(Opm::isfinite(rate[eqIdx]));
|
||||
}
|
||||
}
|
||||
|
||||
aquiferModel_.addToSource(rate, context, spaceIdx, timeIdx);
|
||||
}
|
||||
|
||||
/*!
|
||||
@ -1361,8 +1395,6 @@ public:
|
||||
*
|
||||
* This a hack on top of the maxOilSaturation() hack but it is currently required to
|
||||
* do restart externally. i.e. from the flow code.
|
||||
*
|
||||
* TODO: move the restart-from-ECL-restart-files functionality to EclProblem!
|
||||
*/
|
||||
void setMaxOilSaturation(unsigned globalDofIdx, Scalar value)
|
||||
{
|
||||
@ -1913,8 +1945,6 @@ private:
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
|
||||
// update the hysteresis parameters of the material laws for the whole grid
|
||||
bool updateHysteresis_()
|
||||
{
|
||||
@ -2109,6 +2139,7 @@ private:
|
||||
std::vector<Scalar> maxOilSaturation_;
|
||||
|
||||
EclWellModel wellModel_;
|
||||
EclAquiferModel aquiferModel_;
|
||||
|
||||
std::unique_ptr<EclWriterType> eclWriter_;
|
||||
|
||||
|
Loading…
Reference in New Issue
Block a user