mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Merge pull request #2785 from blattms/improve-parallel-deck-error-handling
Improves error handling and reporting in parallel.
This commit is contained in:
commit
453dec3e26
@ -42,6 +42,8 @@
|
||||
#include <opm/simulators/utils/ParallelEclipseState.hpp>
|
||||
#include <opm/simulators/utils/ParallelSerialization.hpp>
|
||||
|
||||
#include <cstdlib>
|
||||
|
||||
namespace Opm
|
||||
{
|
||||
|
||||
@ -170,9 +172,7 @@ void readDeck(int rank, std::string& deckFilename, std::unique_ptr<Opm::Deck>& d
|
||||
errorGuard = std::make_unique<ErrorGuard>();
|
||||
}
|
||||
|
||||
#if HAVE_MPI
|
||||
int parseSuccess = 0;
|
||||
#endif
|
||||
int parseSuccess = 1; // > 0 is success
|
||||
std::string failureMessage;
|
||||
|
||||
if (rank==0) {
|
||||
@ -223,13 +223,13 @@ void readDeck(int rank, std::string& deckFilename, std::unique_ptr<Opm::Deck>& d
|
||||
}
|
||||
if (!summaryConfig)
|
||||
summaryConfig = std::make_unique<Opm::SummaryConfig>(*deck, *schedule, eclipseState->getTableManager(), *parseContext, *errorGuard);
|
||||
#if HAVE_MPI
|
||||
parseSuccess = 1;
|
||||
#endif
|
||||
|
||||
Opm::checkConsistentArrayDimensions(*eclipseState, *schedule, *parseContext, *errorGuard);
|
||||
}
|
||||
catch(const std::exception& e)
|
||||
{
|
||||
failureMessage = e.what();
|
||||
parseSuccess = 0;
|
||||
}
|
||||
}
|
||||
#if HAVE_MPI
|
||||
@ -242,27 +242,39 @@ void readDeck(int rank, std::string& deckFilename, std::unique_ptr<Opm::Deck>& d
|
||||
eclipseState = std::make_unique<Opm::ParallelEclipseState>();
|
||||
}
|
||||
|
||||
auto comm = Dune::MPIHelper::getCollectiveCommunication();
|
||||
parseSuccess = comm.max(parseSuccess);
|
||||
if (!parseSuccess)
|
||||
try
|
||||
{
|
||||
if (*errorGuard) {
|
||||
errorGuard->dump();
|
||||
errorGuard->clear();
|
||||
}
|
||||
MPI_Abort(MPI_COMM_WORLD, EXIT_FAILURE);
|
||||
Opm::eclStateBroadcast(*eclipseState, *schedule, *summaryConfig);
|
||||
}
|
||||
catch(const std::exception& e)
|
||||
{
|
||||
failureMessage = e.what();
|
||||
parseSuccess = 0;
|
||||
}
|
||||
|
||||
Opm::eclStateBroadcast(*eclipseState, *schedule, *summaryConfig);
|
||||
#endif
|
||||
|
||||
Opm::checkConsistentArrayDimensions(*eclipseState, *schedule, *parseContext, *errorGuard);
|
||||
if (*errorGuard) { // errors encountered
|
||||
parseSuccess = 0;
|
||||
}
|
||||
|
||||
if (*errorGuard) {
|
||||
errorGuard->dump();
|
||||
errorGuard->clear();
|
||||
// print errors and warnings!
|
||||
errorGuard->dump();
|
||||
errorGuard->clear();
|
||||
|
||||
throw std::runtime_error("Unrecoverable errors were encountered while loading input.");
|
||||
auto comm = Dune::MPIHelper::getCollectiveCommunication();
|
||||
parseSuccess = comm.min(parseSuccess);
|
||||
|
||||
if (!parseSuccess)
|
||||
{
|
||||
if (rank == 0)
|
||||
{
|
||||
OpmLog::error(std::string("Unrecoverable errors were encountered while loading input: ")+failureMessage);
|
||||
}
|
||||
#if HAVE_MPI
|
||||
MPI_Finalize();
|
||||
#endif
|
||||
std::exit(EXIT_FAILURE);
|
||||
}
|
||||
}
|
||||
} // end namespace Opm
|
||||
|
Loading…
Reference in New Issue
Block a user