let the aquifiers be managed by core ebos

also, clean them up a bit:
- do not use the intensive quantities cache directly anymore. (i.e.,
  that code should now work if the IQ cache is disabled)
- do not fiddle with the global Jacobian matrix and residual vector
  directly. Instead, implement the water fluxes to the reservoir as a
  source term like wells.

one thing that did not become fully clear to me is if each aquifer
ought to be assumed to be in contact with the whole reservoir or just
a few cells on the boundary. The current implementation goes down the
former path, while, without any deeper knowledge, I would rather
suppose that the latter applies. maybe my understanding of this is
just too limited, though.
This commit is contained in:
Andreas Lauser
2018-10-25 16:47:00 +02:00
parent 985b1b17f0
commit 7c81dbdaab
5 changed files with 69 additions and 159 deletions

View File

@@ -45,6 +45,7 @@ namespace Opm
typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
typedef typename GET_PROP_TYPE(TypeTag, Indices) BlackoilIndices;
typedef typename GET_PROP_TYPE(TypeTag, RateVector) RateVector;
typedef typename GET_PROP_TYPE(TypeTag, IntensiveQuantities) IntensiveQuantities;
enum { enableTemperature = GET_PROP_VALUE(TypeTag, EnableTemperature) };
enum { enableEnergy = GET_PROP_VALUE(TypeTag, EnableEnergy) };
@@ -61,8 +62,8 @@ namespace Opm
AquiferCarterTracy( const AquiferCT::AQUCT_data& aquct_data,
const Aquancon::AquanconOutput& connection,
Simulator& ebosSimulator )
const Aquancon::AquanconOutput& connection,
const Simulator& ebosSimulator )
: ebos_simulator_ (ebosSimulator),
aquct_data_ (aquct_data),
gravity_ (ebos_simulator_.problem().gravity()[2])
@@ -70,56 +71,45 @@ namespace Opm
initQuantities(connection);
}
inline void assembleAquiferEq(const SimulatorTimerInterface& timer)
template <class Context>
void addToSource(RateVector& rates, const Context& context, unsigned spaceIdx, unsigned timeIdx)
{
auto& ebosJac = ebos_simulator_.model().linearizer().matrix();
auto& ebosResid = ebos_simulator_.model().linearizer().residual();
unsigned idx = context.globalSpaceIndex(spaceIdx, timeIdx);
size_t cellID;
for ( size_t idx = 0; idx < cell_idx_.size(); ++idx )
{
Eval qinflow = 0.0;
cellID = cell_idx_.at(idx);
// We are dereferencing the value of IntensiveQuantities because cachedIntensiveQuantities return a const pointer to
// IntensiveQuantities of that particular cell_id
const IntensiveQuantities intQuants = *(ebos_simulator_.model().cachedIntensiveQuantities(cellID, /*timeIdx=*/ 0));
// This is the pressure at td + dt
updateCellPressure(pressure_current_,idx,intQuants);
updateCellDensity(idx,intQuants);
calculateInflowRate(idx, timer);
qinflow = Qai_.at(idx);
ebosResid[cellID][waterCompIdx] -= qinflow.value();
// We are dereferencing the value of IntensiveQuantities because cachedIntensiveQuantities return a const pointer to
// IntensiveQuantities of that particular cell_id
const IntensiveQuantities intQuants = context.intensiveQuantities(spaceIdx, timeIdx);
// This is the pressure at td + dt
updateCellPressure(pressure_current_,idx,intQuants);
updateCellDensity(idx,intQuants);
calculateInflowRate(idx, context.simulator());
for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
{
// also need to consider the efficiency factor when manipulating the jacobians.
ebosJac[cellID][cellID][waterCompIdx][pvIdx] -= qinflow.derivative(pvIdx);
}
}
rates[BlackoilIndices::conti0EqIdx + FluidSystem::waterCompIdx] -=
Qai_.at(idx)/context.dofVolume(spaceIdx, timeIdx);
}
inline void beforeTimeStep(const SimulatorTimerInterface& timer)
inline void beforeTimeStep(const Simulator& simulator)
{
auto cellID = cell_idx_.begin();
size_t idx;
for ( idx = 0; cellID != cell_idx_.end(); ++cellID, ++idx )
{
const auto& intQuants = *(ebos_simulator_.model().cachedIntensiveQuantities(*cellID, /*timeIdx=*/ 0));
const auto& intQuants = *(simulator.model().cachedIntensiveQuantities(*cellID, /*timeIdx=*/ 0));
updateCellPressure(pressure_previous_ ,idx,intQuants);
}
}
inline void afterTimeStep(const SimulatorTimerInterface& timer)
inline void afterTimeStep(const Simulator& simulator)
{
for (auto Qai = Qai_.begin(); Qai != Qai_.end(); ++Qai)
{
W_flux_ += (*Qai)*timer.currentStepLength();
W_flux_ += (*Qai)*simulator.timeStepSize();
}
}
private:
Simulator& ebos_simulator_;
const Simulator& ebos_simulator_;
// Grid variables
std::vector<size_t> cell_idx_;
@@ -133,6 +123,8 @@ namespace Opm
std::vector<Eval> rhow_;
std::vector<Scalar> alphai_;
std::vector<Eval> aquiferWaterInflux_;
// Variables constants
const AquiferCT::AQUCT_data aquct_data_;
@@ -164,6 +156,7 @@ namespace Opm
calculateAquiferConstants();
aquiferWaterInflux_.resize(cell_idx_.size());
pressure_previous_.resize(cell_idx_.size(), 0.);
pressure_current_.resize(cell_idx_.size(), 0.);
Qai_.resize(cell_idx_.size(), 0.0);
@@ -194,10 +187,10 @@ namespace Opm
}
// This function implements Eqs 5.8 and 5.9 of the EclipseTechnicalDescription
inline void calculateEqnConstants(Scalar& a, Scalar& b, const int idx, const SimulatorTimerInterface& timer)
inline void calculateEqnConstants(Scalar& a, Scalar& b, const int idx, const Simulator& simulator)
{
const Scalar td_plus_dt = (timer.currentStepLength() + timer.simulationTimeElapsed()) / Tc_;
const Scalar td = timer.simulationTimeElapsed() / Tc_;
const Scalar td_plus_dt = (simulator.timeStepSize() + simulator.time()) / Tc_;
const Scalar td = simulator.time() / Tc_;
Scalar PItdprime = 0.;
Scalar PItd = 0.;
getInfluenceTableValues(PItd, PItdprime, td_plus_dt);
@@ -206,10 +199,10 @@ namespace Opm
}
// This function implements Eq 5.7 of the EclipseTechnicalDescription
inline void calculateInflowRate(int idx, const SimulatorTimerInterface& timer)
inline void calculateInflowRate(int idx, const Simulator& simulator)
{
Scalar a, b;
calculateEqnConstants(a,b,idx,timer);
calculateEqnConstants(a,b,idx,simulator);
Qai_.at(idx) = alphai_.at(idx)*( a - b * ( pressure_current_.at(idx) - pressure_previous_.at(idx) ) );
}
@@ -365,4 +358,4 @@ namespace Opm
} // namespace Opm
#endif
#endif

View File

@@ -24,6 +24,8 @@
#ifndef OPM_BLACKOILAQUIFERMODEL_HEADER_INCLUDED
#define OPM_BLACKOILAQUIFERMODEL_HEADER_INCLUDED
#include <ebos/eclbaseaquifermodel.hh>
#include <opm/parser/eclipse/EclipseState/AquiferCT.hpp>
#include <opm/parser/eclipse/EclipseState/Aquancon.hpp>
#include <opm/simulators/timestepping/SimulatorTimer.hpp>
@@ -34,45 +36,46 @@ namespace Opm {
/// Class for handling the blackoil well model.
template<typename TypeTag>
class BlackoilAquiferModel {
class BlackoilAquiferModel : public Ewoms::EclBaseAquiferModel<TypeTag>
{
typedef Ewoms::EclBaseAquiferModel<TypeTag> ParentType;
typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
typedef typename GET_PROP_TYPE(TypeTag, RateVector) RateVector;
public:
explicit BlackoilAquiferModel(Simulator& ebosSimulator);
// at the beginning of each time step (Not report step)
void beginTimeStep();
// add the water rate due to aquifers to the source term.
template <class Context>
void addToSource(RateVector& rates,
const Context& context,
unsigned spaceIdx,
unsigned timeIdx) const
{
for (auto aquifer = aquifers_.begin(); aquifer != aquifers_.end(); ++aquifer)
aquifer->addToSource(rates, context, spaceIdx, timeIdx);
}
protected:
// --------- Types ---------
typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
typedef AquiferCarterTracy<TypeTag> Aquifer_object;
explicit BlackoilAquiferModel(Simulator& ebosSimulator);
// compute the well fluxes and assemble them in to the reservoir equations as source terms
// and in the well equations.
void assemble( const SimulatorTimerInterface& timer,
const int iterationIdx );
// called at the end of a time step
void timeStepSucceeded(const SimulatorTimerInterface& timer);
protected:
Simulator& ebosSimulator_;
std::vector<Aquifer_object> aquifers_;
// TODO: declaring this to be mutable is a hack which should be fixed in the
// long term
mutable std::vector<Aquifer_object> aquifers_;
// This initialization function is used to connect the parser objects with the ones needed by AquiferCarterTracy
void init();
void updateConnectionIntensiveQuantities() const;
void assembleAquiferEq(const SimulatorTimerInterface& timer);
// at the beginning of each time step (Not report step)
void prepareTimeStep(const SimulatorTimerInterface& timer);
bool aquiferActive() const;
};

View File

@@ -4,87 +4,22 @@ namespace Opm {
template<typename TypeTag>
BlackoilAquiferModel<TypeTag>::
BlackoilAquiferModel(Simulator& ebosSimulator)
: ebosSimulator_(ebosSimulator)
: ParentType(ebosSimulator)
{
init();
}
// called at the end of a time step
template<typename TypeTag>
void
BlackoilAquiferModel<TypeTag>:: timeStepSucceeded(const SimulatorTimerInterface& timer)
{
if ( !aquiferActive() ) {
return;
}
for (auto aquifer = aquifers_.begin(); aquifer != aquifers_.end(); ++aquifer)
{
aquifer->afterTimeStep(timer);
}
}
template<typename TypeTag>
void
BlackoilAquiferModel<TypeTag>::
assemble( const SimulatorTimerInterface& timer,
const int iterationIdx )
{
if ( !aquiferActive() ) {
return;
}
// We need to update the reservoir pressures connected to the aquifer
updateConnectionIntensiveQuantities();
if (iterationIdx == 0) {
// We can do the Table check and coefficients update in this function
// For now, it does nothing!
prepareTimeStep(timer);
}
assembleAquiferEq(timer);
}
template<typename TypeTag>
void
BlackoilAquiferModel<TypeTag>:: updateConnectionIntensiveQuantities() const
{
ElementContext elemCtx(ebosSimulator_);
const auto& gridView = ebosSimulator_.gridView();
const auto& elemEndIt = gridView.template end</*codim=*/0, Dune::Interior_Partition>();
for (auto elemIt = gridView.template begin</*codim=*/0, Dune::Interior_Partition>();
elemIt != elemEndIt;
++elemIt)
{
elemCtx.updatePrimaryStencil(*elemIt);
elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
}
}
// Protected function which calls the individual aquifer models
template<typename TypeTag>
void
BlackoilAquiferModel<TypeTag>:: assembleAquiferEq(const SimulatorTimerInterface& timer)
{
for (auto aquifer = aquifers_.begin(); aquifer != aquifers_.end(); ++aquifer)
{
aquifer->assembleAquiferEq(timer);
}
}
// Protected function
// some preparation work, mostly related to group control and RESV,
// at the beginning of each time step (Not report step)
template<typename TypeTag>
void BlackoilAquiferModel<TypeTag>:: prepareTimeStep(const SimulatorTimerInterface& timer)
void
BlackoilAquiferModel<TypeTag>::beginTimeStep()
{
// Here we can ask each carter tracy aquifers to get the current previous time step's pressure
for (auto aquifer = aquifers_.begin(); aquifer != aquifers_.end(); ++aquifer)
{
aquifer->beforeTimeStep(timer);
aquifer->beforeTimeStep(this->simulator_);
}
}
@@ -93,14 +28,14 @@ namespace Opm {
void
BlackoilAquiferModel<TypeTag>:: init()
{
const auto& deck = ebosSimulator_.vanguard().deck();
const auto& deck = this->simulator_.vanguard().deck();
if ( !deck.hasKeyword("AQUCT") ) {
return ;
}
updateConnectionIntensiveQuantities();
const auto& eclState = ebosSimulator_.vanguard().eclState();
//updateConnectionIntensiveQuantities();
const auto& eclState = this->simulator_.vanguard().eclState();
// Get all the carter tracy aquifer properties data and put it in aquifers vector
const AquiferCT aquiferct = AquiferCT(eclState,deck);
@@ -115,7 +50,7 @@ namespace Opm {
for (size_t i = 0; i < aquifersData.size(); ++i)
{
aquifers_.push_back(
AquiferCarterTracy<TypeTag> (aquifersData.at(i), aquifer_connection.at(i), ebosSimulator_)
AquiferCarterTracy<TypeTag> (aquifersData.at(i), aquifer_connection.at(i), this->simulator_)
);
}
}

View File

@@ -78,6 +78,8 @@ SET_BOOL_PROP(EclFlowProblem, EnableDebuggingChecks, false);
SET_BOOL_PROP(EclFlowProblem, BlackoilConserveSurfaceVolume, true);
SET_BOOL_PROP(EclFlowProblem, UseVolumetricResidual, false);
SET_TYPE_PROP(EclFlowProblem, EclAquiferModel, Opm::BlackoilAquiferModel<TypeTag>);
// disable all extensions supported by black oil model. this should not really be
// necessary but it makes things a bit more explicit
SET_BOOL_PROP(EclFlowProblem, EnablePolymer, false);
@@ -145,7 +147,6 @@ namespace Opm {
BlackoilModelEbos(Simulator& ebosSimulator,
const ModelParameters& param,
BlackoilWellModel<TypeTag>& well_model,
BlackoilAquiferModel<TypeTag>& aquifer_model,
const NewtonIterationBlackoilInterface& linsolver,
const bool terminal_output)
: ebosSimulator_(ebosSimulator)
@@ -159,7 +160,6 @@ namespace Opm {
, has_energy_(GET_PROP_VALUE(TypeTag, EnableEnergy))
, param_( param )
, well_model_ (well_model)
, aquifer_model_(aquifer_model)
, terminal_output_ (terminal_output)
, current_relaxation_(1.0)
, dx_old_(UgGridHelpers::numCells(grid_))
@@ -342,7 +342,6 @@ namespace Opm {
void afterStep(const SimulatorTimerInterface& OPM_UNUSED timer)
{
wellModel().timeStepSucceeded(timer.simulationTimeElapsed());
aquiferModel().timeStepSucceeded(timer);
ebosSimulator_.problem().endTimeStep();
}
@@ -360,17 +359,6 @@ namespace Opm {
ebosSimulator_.model().linearizer().linearize();
ebosSimulator_.problem().endIteration();
// -------- Aquifer models ----------
try
{
// Modify the Jacobian and residuals according to the aquifer models
aquiferModel().assemble(timer, iterationIdx);
}
catch( ... )
{
OPM_THROW(Opm::NumericalIssue,"Error when assembling aquifer models");
}
// -------- Current time step length ----------
const double dt = timer.currentStepLength();
@@ -959,9 +947,6 @@ namespace Opm {
// Well Model
BlackoilWellModel<TypeTag>& well_model_;
// Aquifer Model
BlackoilAquiferModel<TypeTag>& aquifer_model_;
/// \brief Whether we print something to std::cout
bool terminal_output_;
/// \brief The number of cells of the global grid.
@@ -981,9 +966,6 @@ namespace Opm {
const BlackoilWellModel<TypeTag>&
wellModel() const { return well_model_; }
BlackoilAquiferModel<TypeTag>&
aquiferModel() { return aquifer_model_; }
void beginReportStep()
{
ebosSimulator_.problem().beginEpisode();

View File

@@ -201,8 +201,6 @@ public:
ebosSimulator_.model().addAuxiliaryModule(wellAuxMod_.get());
}
AquiferModel aquifer_model(ebosSimulator_);
// Main simulation loop.
while (!timer.done()) {
// Report timestep.
@@ -217,7 +215,7 @@ public:
wellModel.beginReportStep(timer.currentStepNum());
auto solver = createSolver(wellModel, aquifer_model);
auto solver = createSolver(wellModel);
// write the inital state at the report stage
if (timer.initialStep()) {
@@ -343,12 +341,11 @@ public:
protected:
std::unique_ptr<Solver> createSolver(WellModel& wellModel, AquiferModel& aquifer_model)
std::unique_ptr<Solver> createSolver(WellModel& wellModel)
{
auto model = std::unique_ptr<Model>(new Model(ebosSimulator_,
modelParam_,
wellModel,
aquifer_model,
linearSolver_,
terminalOutput_));