Simulator, Model: add everthing which is required by the polymer simulators

basically, this adds some hooks to the SimulatorBase class and the
model and introduces a few types to the SimulatorTraits.
This commit is contained in:
Andreas Lauser 2015-05-27 23:12:34 +02:00
parent 45fb80f547
commit a154c8394d
3 changed files with 82 additions and 39 deletions

View File

@ -74,13 +74,22 @@
namespace Opm
{
template<class Simulator>
struct SimulatorTraits;
/// Class collecting all necessary components for a two-phase simulation.
template<class GridT, class Implementation>
template<class Implementation>
class SimulatorBase
{
typedef SimulatorTraits<Implementation> Traits;
public:
/// \brief The type of the grid that we use.
typedef GridT Grid;
typedef typename Traits::ReservoirState ReservoirState;
typedef typename Traits::WellState WellState;
typedef typename Traits::OutputWriter OutputWriter;
typedef typename Traits::Grid Grid;
typedef typename Traits::Model Model;
/// Initialise from parameters and objects to observe.
/// \param[in] param parameters, this class accepts the following:
/// parameter (default) effect
@ -118,7 +127,7 @@ namespace Opm
const bool disgas,
const bool vapoil,
std::shared_ptr<EclipseState> eclipse_state,
BlackoilOutputWriter& output_writer,
OutputWriter& output_writer,
const std::vector<double>& threshold_pressures_by_face);
/// Run the simulation.
@ -129,17 +138,39 @@ namespace Opm
/// \param[in,out] well_state state of wells: bhp, perforation rates
/// \return simulation report, with timing data
SimulatorReport run(SimulatorTimer& timer,
BlackoilState& state);
ReservoirState& state);
protected:
Implementation& asImp_() { return *static_cast<Implementation*>(this); }
const Implementation& asImp_() const { return *static_cast<const Implementation*>(this); }
void handleAdditionalWellInflow(SimulatorTimer& timer,
WellsManager& wells_manager,
WellState& well_state,
const Wells* wells)
{};
std::shared_ptr<Model> createModel(const typename Model::ModelParameters &modelParams,
const Wells* wells)
{
return std::make_shared<Model>(modelParams,
grid_,
props_,
geo_,
rock_comp_props_,
wells,
solver_,
has_disgas_,
has_vapoil_,
terminal_output_);
}
void
computeRESV(const std::size_t step,
const Wells* wells,
const BlackoilState& x,
WellStateFullyImplicitBlackoil& xw);
WellState& xw);
// Data.
typedef RateConverter::
@ -164,7 +195,7 @@ namespace Opm
// eclipse_state
std::shared_ptr<EclipseState> eclipse_state_;
// output_writer
BlackoilOutputWriter& output_writer_;
OutputWriter& output_writer_;
RateConverterType rateConverter_;
// Threshold pressures.
std::vector<double> threshold_pressures_by_face_;

View File

@ -20,8 +20,8 @@
namespace Opm
{
template<class GridT, class Implementation>
SimulatorBase<GridT, Implementation>::SimulatorBase(const parameter::ParameterGroup& param,
template<class Implementation>
SimulatorBase<Implementation>::SimulatorBase(const parameter::ParameterGroup& param,
const Grid& grid,
const DerivedGeology& geo,
BlackoilPropsAdInterface& props,
@ -31,9 +31,8 @@ namespace Opm
const bool has_disgas,
const bool has_vapoil,
std::shared_ptr<EclipseState> eclipse_state,
BlackoilOutputWriter& output_writer,
OutputWriter& output_writer,
const std::vector<double>& threshold_pressures_by_face)
: param_(param),
grid_(grid),
props_(props),
@ -69,11 +68,11 @@ namespace Opm
#endif
}
template<class GridT, class Implementation>
SimulatorReport SimulatorBase<GridT, Implementation>::run(SimulatorTimer& timer,
BlackoilState& state)
template<class Implementation>
SimulatorReport SimulatorBase<Implementation>::run(SimulatorTimer& timer,
ReservoirState& state)
{
WellStateFullyImplicitBlackoil prev_well_state;
WellState prev_well_state;
// Create timers and file for writing timing info.
Opm::time::StopWatch solver_timer;
@ -84,8 +83,6 @@ namespace Opm
std::string tstep_filename = output_writer_.outputDirectory() + "/step_timing.txt";
std::ofstream tstep_os(tstep_filename.c_str());
typedef GridT Grid;
typedef BlackoilModel<Grid> Model;
typedef typename Model::ModelParameters ModelParams;
ModelParams modelParams( param_ );
typedef NewtonSolver<Model> Solver;
@ -134,9 +131,12 @@ namespace Opm
props_.permeability(),
is_parallel_run_);
const Wells* wells = wells_manager.c_wells();
WellStateFullyImplicitBlackoil well_state;
WellState well_state;
well_state.init(wells, state, prev_well_state);
// give the polymer and surfactant simulators the chance to do their stuff
asImp_().handleAdditionalWellInflow(timer, wells_manager, well_state, wells);
// write simulation state at the report stage
output_writer_.writeTimeStep( timer, state, well_state );
@ -145,16 +145,16 @@ namespace Opm
props_.updateSatHyst(state.saturation(), allcells_);
// Compute reservoir volumes for RESV controls.
computeRESV(timer.currentStepNum(), wells, state, well_state);
asImp_().computeRESV(timer.currentStepNum(), wells, state, well_state);
// Run a multiple steps of the solver depending on the time step control.
solver_timer.start();
Model model(modelParams, grid_, props_, geo_, rock_comp_props_, wells, solver_, has_disgas_, has_vapoil_, terminal_output_);
auto model = asImp_().createModel(modelParams, wells);
if (!threshold_pressures_by_face_.empty()) {
model.setThresholdPressures(threshold_pressures_by_face_);
model->setThresholdPressures(threshold_pressures_by_face_);
}
Solver solver(solverParams, model);
Solver solver(solverParams, *model);
// If sub stepping is enabled allow the solver to sub cycle
// in case the report steps are to large for the solver to converge
@ -327,11 +327,11 @@ namespace Opm
}
} // namespace SimFIBODetails
template <class GridT, class Implementation>
void SimulatorBase<GridT, Implementation>::computeRESV(const std::size_t step,
template <class Implementation>
void SimulatorBase<Implementation>::computeRESV(const std::size_t step,
const Wells* wells,
const BlackoilState& x,
WellStateFullyImplicitBlackoil& xw)
WellState& xw)
{
typedef SimFIBODetails::WellMap WellMap;

View File

@ -23,13 +23,25 @@
#include "SimulatorBase.hpp"
namespace Opm {
template<class GridT>
class SimulatorFullyImplicitBlackoil;
template<class GridT>
struct SimulatorTraits<SimulatorFullyImplicitBlackoil<GridT> >
{
typedef WellStateFullyImplicitBlackoil WellState;
typedef BlackoilState ReservoirState;
typedef BlackoilOutputWriter OutputWriter;
typedef GridT Grid;
typedef BlackoilModel<Grid> Model;
};
/// a simulator for the blackoil model
template<class GridT>
class SimulatorFullyImplicitBlackoil
: public SimulatorBase<GridT, SimulatorFullyImplicitBlackoil<GridT> >
: public SimulatorBase<SimulatorFullyImplicitBlackoil<GridT> >
{
typedef SimulatorBase<GridT, SimulatorFullyImplicitBlackoil<GridT> > Base;
typedef SimulatorBase<SimulatorFullyImplicitBlackoil<GridT> > Base;
public:
// forward the constructor to the base class
SimulatorFullyImplicitBlackoil(const parameter::ParameterGroup& param,