Merge pull request #466 from blattms/fix-MatrixBlockError-on-some-processes

Introduces parallel exception handling for ILU0 preconditioner.
This commit is contained in:
Atgeirr Flø Rasmussen
2015-09-14 08:19:03 +02:00
2 changed files with 55 additions and 8 deletions

View File

@@ -178,7 +178,31 @@ createILU0Ptr(const M& A, double relax,
Dune::OwnerOverlapCopyCommunication<I1,I2>,
Dune::SeqILU0<M,X,X>
> PointerType;
Dune::SeqILU0<M,X,X>* ilu = new Dune::SeqILU0<M,X,X>(A, relax);
Dune::SeqILU0<M,X,X>* ilu = nullptr;
int ilu_setup_successful = 1;
std::string message;
try {
ilu = new Dune::SeqILU0<M,X,X>(A, relax);
}
catch ( Dune::MatrixBlockError error )
{
message = error.what();
std::cerr<<"Exception occured on process " <<
comm.communicator().rank() << " during " <<
"setup of ILU0 preconditioner with message: " <<
message<<std::endl;
ilu_setup_successful = 0;
}
// Check whether there was a problem on some process
if ( comm.communicator().min(ilu_setup_successful) == 0 )
{
if ( ilu ) // not null if constructor succeeded
{
// prevent memory leak
delete ilu;
}
throw Dune::MatrixBlockError();
}
return typename SelectParallelILUSharedPtr<Dune::SeqILU0<M,X,X>, I1, I2>
::type ( new PointerType(*ilu, comm), createParallelDeleter(*ilu, comm));
}

View File

@@ -132,14 +132,37 @@ namespace Opm
AdditionalObjectDeleter<SeqPreconditioner> >
constructPrecond(Operator& opA, const Comm& comm) const
{
const double relax = 1.0;
SeqPreconditioner* seq_precond= new SeqPreconditioner(opA.getmat(),
relax);
typedef AdditionalObjectDeleter<SeqPreconditioner> Deleter;
std::unique_ptr<ParPreconditioner, Deleter>
precond(new ParPreconditioner(*seq_precond, comm),
Deleter(*seq_precond));
return precond;
typedef std::unique_ptr<ParPreconditioner, Deleter> Pointer;
int ilu_setup_successful = 1;
std::string message;
const double relax = 1.0;
SeqPreconditioner* seq_precond = nullptr;
try {
seq_precond = new SeqPreconditioner(opA.getmat(),
relax);
}
catch ( Dune::MatrixBlockError error )
{
message = error.what();
std::cerr<<"Exception occured on process " <<
comm.communicator().rank() << " during " <<
"setup of ILU0 preconditioner with message: " <<
message<<std::endl;
ilu_setup_successful = 0;
}
// Check whether there was a problem on some process
if ( comm.communicator().min(ilu_setup_successful) == 0 )
{
if ( seq_precond ) // not null if constructor succeeded
{
// prevent memory leak
delete seq_precond;
}
throw Dune::MatrixBlockError();
}
return Pointer(new ParPreconditioner(*seq_precond, comm),
Deleter(*seq_precond));
}
#endif