Skip to content
Snippets Groups Projects
Commit f0ae40e8 authored by Christophe Geuzaine's avatar Christophe Geuzaine
Browse files

support for petsc-dev

parent 48310585
No related branches found
No related tags found
No related merge requests found
......@@ -105,8 +105,10 @@ bool eigenSolver::solve(int numEigenValues, std::string which)
Msg::Error("SLEPc diverged after %d iterations", its);
else if(reason == EPS_DIVERGED_BREAKDOWN)
Msg::Error("SLEPc generic breakdown in method");
#if !(PETSC_VERSION_RELEASE == 0) // petsc-dev
else if(reason == EPS_DIVERGED_NONSYMMETRIC)
Msg::Error("The operator is nonsymmetric");
#endif
// get number of converged approximate eigenpairs
PetscInt nconv;
......@@ -152,11 +154,20 @@ bool eigenSolver::solve(int numEigenValues, std::string which)
}
_eigenVectors.push_back(ev);
}
#if (PETSC_VERSION_RELEASE == 0) // petsc-dev
_try(VecDestroy(&xr));
_try(VecDestroy(&xi));
#else
_try(VecDestroy(xr));
_try(VecDestroy(xi));
#endif
}
#if (PETSC_VERSION_RELEASE == 0) // petsc-dev
_try(EPSDestroy(&eps));
#else
_try(EPSDestroy(eps));
#endif
if(reason == EPS_CONVERGED_TOL){
Msg::Debug("SLEPc done");
......
......@@ -145,9 +145,15 @@ bool linearSystemPETScBlockDouble::isAllocated() const
void linearSystemPETScBlockDouble::clear()
{
if(_isAllocated){
#if (PETSC_VERSION_RELEASE == 0) // petsc-dev
MatDestroy(&_a);
VecDestroy(&_x);
VecDestroy(&_b);
#else
MatDestroy(_a);
VecDestroy(_x);
VecDestroy(_b);
#endif
}
_isAllocated = false;
}
......
......@@ -2,6 +2,12 @@
#include <petsc.h>
#include <petscksp.h>
#include "linearSystemPETSc.h"
#if (PETSC_VERSION_RELEASE == 0) // petsc-dev
#define PetscTruth PetscBool
#define PetscOptionsGetTruth PetscOptionsGetBool
#endif
static void _try(int ierr)
{
CHKERRABORT(PETSC_COMM_WORLD, ierr);
......@@ -40,7 +46,11 @@ linearSystemPETSc<scalar>::~linearSystemPETSc()
{
clear();
if(_kspAllocated)
#if (PETSC_VERSION_RELEASE == 0) // petsc-dev
_try(KSPDestroy(&_ksp));
#else
_try(KSPDestroy(_ksp));
#endif
}
template <class scalar>
......@@ -141,9 +151,15 @@ template <class scalar>
void linearSystemPETSc<scalar>::clear()
{
if(_isAllocated){
#if (PETSC_VERSION_RELEASE == 0) // petsc-dev
_try(MatDestroy(&_a));
_try(VecDestroy(&_x));
_try(VecDestroy(&_b));
#else
_try(MatDestroy(_a));
_try(VecDestroy(_x));
_try(VecDestroy(_b));
#endif
}
_isAllocated = false;
}
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment