Refactor flow_ebos_blackoil.cpp (2)

NOTE: this pull request depends on #2555 which should be merged first.

A rewrite of the outdated PR #2543.

Refactors flow_ebos_blackoil.cpp such that we can choose not to execute
the whole simulation using the flowEbosBlackoilMain() function but
instead only initialize by calling flowEbosBlackoilMainInit(). This is
necessary to implement a Python step() method that can advance the
simulator one report step at a time.

Also adds a method initFlowEbosBlackoil() to Main.hpp that can be used
directly from the Python interface's BlackOilSimulator object to gain
access to the FlowMainEbos object before it has initialized the
simulation main loop.
This commit is contained in:
Håkon Hægland 2020-05-11 17:57:43 +02:00
parent df58d44f70
commit 946b5f5806
3 changed files with 40 additions and 4 deletions

View File

@ -46,8 +46,8 @@ void flowEbosBlackoilSetDeck(double setupTime, Deck *deck, EclipseState& eclStat
Vanguard::setExternalSummaryConfig(&summaryConfig); Vanguard::setExternalSummaryConfig(&summaryConfig);
} }
// ----------------- Main program ----------------- std::unique_ptr<Opm::FlowMainEbos<TTAG(EclFlowProblem)>>
int flowEbosBlackoilMain(int argc, char** argv, bool outputCout, bool outputFiles) flowEbosBlackoilMainInit(int argc, char** argv, bool outputCout, bool outputFiles)
{ {
// we always want to use the default locale, and thus spare us the trouble // we always want to use the default locale, and thus spare us the trouble
// with incorrect locale settings. // with incorrect locale settings.
@ -59,8 +59,14 @@ int flowEbosBlackoilMain(int argc, char** argv, bool outputCout, bool outputFile
Dune::MPIHelper::instance(argc, argv); Dune::MPIHelper::instance(argc, argv);
#endif #endif
Opm::FlowMainEbos<TTAG(EclFlowProblem)> mainfunc; return std::make_unique<Opm::FlowMainEbos<TTAG(EclFlowProblem)>>();
return mainfunc.execute(argc, argv, outputCout, outputFiles); }
// ----------------- Main program -----------------
int flowEbosBlackoilMain(int argc, char** argv, bool outputCout, bool outputFiles)
{
auto mainfunc = flowEbosBlackoilMainInit(argc, argv, outputCout, outputFiles);
return mainfunc->execute(argc, argv, outputCout, outputFiles);
} }
} }

View File

@ -21,10 +21,15 @@
#include <opm/parser/eclipse/EclipseState/EclipseState.hpp> #include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Schedule.hpp> #include <opm/parser/eclipse/EclipseState/Schedule/Schedule.hpp>
#include <opm/parser/eclipse/EclipseState/SummaryConfig/SummaryConfig.hpp> #include <opm/parser/eclipse/EclipseState/SummaryConfig/SummaryConfig.hpp>
#include <opm/simulators/flow/FlowMainEbos.hpp>
namespace Opm { namespace Opm {
void flowEbosBlackoilSetDeck(double setupTime, Deck *deck, EclipseState& eclState, Schedule& schedule, SummaryConfig& summaryConfig); void flowEbosBlackoilSetDeck(double setupTime, Deck *deck, EclipseState& eclState, Schedule& schedule, SummaryConfig& summaryConfig);
int flowEbosBlackoilMain(int argc, char** argv, bool outputCout, bool outputFiles); int flowEbosBlackoilMain(int argc, char** argv, bool outputCout, bool outputFiles);
std::unique_ptr<Opm::FlowMainEbos<TTAG(EclFlowProblem)>>
flowEbosBlackoilMainInit(int argc, char** argv, bool outputCout, bool outputFiles);
} }
#endif // FLOW_EBOS_BLACKOIL_HPP #endif // FLOW_EBOS_BLACKOIL_HPP

View File

@ -117,6 +117,7 @@ namespace Opm
class Main class Main
{ {
private: private:
using FlowMainEbosType = Opm::FlowMainEbos<TTAG(EclFlowProblem)>;
enum class FileOutputMode { enum class FileOutputMode {
//! \brief No output to files. //! \brief No output to files.
OUTPUT_NONE = 0, OUTPUT_NONE = 0,
@ -174,6 +175,30 @@ namespace Opm
} }
} }
// To be called from the Python interface code. Only do the
// initialization and then return a pointer to the FlowEbosMain
// object that can later be accessed directly from the Python interface
// to e.g. advance the simulator one report step
std::unique_ptr<FlowMainEbosType> initFlowEbosBlackoil(int& exitCode)
{
exitCode = EXIT_SUCCESS;
if (initialize_<TTAG(FlowEarlyBird)>(exitCode)) {
// TODO: check that this deck really represents a blackoil
// case. E.g. check that number of phases == 3
Opm::flowEbosBlackoilSetDeck(
setupTime_,
deck_.get(),
*eclipseState_,
*schedule_,
*summaryConfig_);
return Opm::flowEbosBlackoilMainInit(
argc_, argv_, outputCout_, outputFiles_);
} else {
exitCode = EXIT_FAILURE;
return std::unique_ptr<FlowMainEbosType>(); // nullptr
}
}
private: private:
int dispatchDynamic_() int dispatchDynamic_()
{ {