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

handle petsc 3.3

parent 47d0ef33
No related branches found
No related tags found
No related merge requests found
...@@ -651,7 +651,7 @@ void Msg::ExchangeOnelabParameter(const std::string &key, ...@@ -651,7 +651,7 @@ void Msg::ExchangeOnelabParameter(const std::string &key,
bool noRange = true, noChoices = true, noLoop = true, noGraph = true; bool noRange = true, noChoices = true, noLoop = true, noGraph = true;
if(ps.size()){ if(ps.size()){
if(ps[0].getReadOnly()) if(ps[0].getReadOnly())
ps[0].setValue(val[0]); // use value from gmsh ps[0].setValue(val[0]); // use value from gmsh (so it is updated if necessary)
else else
val[0] = ps[0].getValue(); // use value from server val[0] = ps[0].getValue(); // use value from server
// keep track of these attributes, which can be changed server-side // keep track of these attributes, which can be changed server-side
......
...@@ -105,7 +105,7 @@ bool eigenSolver::solve(int numEigenValues, std::string which) ...@@ -105,7 +105,7 @@ bool eigenSolver::solve(int numEigenValues, std::string which)
Msg::Error("SLEPc diverged after %d iterations", its); Msg::Error("SLEPc diverged after %d iterations", its);
else if(reason == EPS_DIVERGED_BREAKDOWN) else if(reason == EPS_DIVERGED_BREAKDOWN)
Msg::Error("SLEPc generic breakdown in method"); Msg::Error("SLEPc generic breakdown in method");
#if !(PETSC_VERSION_RELEASE == 0 || ((PETSC_VERSION_MAJOR == 3) && (PETSC_VERSION_MINOR == 2))) // petsc-dev #if !(PETSC_VERSION_RELEASE == 0 || ((PETSC_VERSION_MAJOR == 3) && (PETSC_VERSION_MINOR >= 2)))
else if(reason == EPS_DIVERGED_NONSYMMETRIC) else if(reason == EPS_DIVERGED_NONSYMMETRIC)
Msg::Error("The operator is nonsymmetric"); Msg::Error("The operator is nonsymmetric");
#endif #endif
...@@ -154,7 +154,7 @@ bool eigenSolver::solve(int numEigenValues, std::string which) ...@@ -154,7 +154,7 @@ bool eigenSolver::solve(int numEigenValues, std::string which)
} }
_eigenVectors.push_back(ev); _eigenVectors.push_back(ev);
} }
#if (PETSC_VERSION_RELEASE == 0 || ((PETSC_VERSION_MAJOR == 3) && (PETSC_VERSION_MINOR == 2))) // petsc-dev #if (PETSC_VERSION_RELEASE == 0 || ((PETSC_VERSION_MAJOR == 3) && (PETSC_VERSION_MINOR >= 2)))
_try(VecDestroy(&xr)); _try(VecDestroy(&xr));
_try(VecDestroy(&xi)); _try(VecDestroy(&xi));
#else #else
...@@ -163,7 +163,7 @@ bool eigenSolver::solve(int numEigenValues, std::string which) ...@@ -163,7 +163,7 @@ bool eigenSolver::solve(int numEigenValues, std::string which)
#endif #endif
} }
#if (PETSC_VERSION_RELEASE == 0 || ((PETSC_VERSION_MAJOR == 3) && (PETSC_VERSION_MINOR == 2))) // petsc-dev #if (PETSC_VERSION_RELEASE == 0 || ((PETSC_VERSION_MAJOR == 3) && (PETSC_VERSION_MINOR >= 2)))
_try(EPSDestroy(&eps)); _try(EPSDestroy(&eps));
#else #else
_try(EPSDestroy(eps)); _try(EPSDestroy(eps));
......
...@@ -164,7 +164,7 @@ bool linearSystemPETScBlockDouble::isAllocated() const ...@@ -164,7 +164,7 @@ bool linearSystemPETScBlockDouble::isAllocated() const
void linearSystemPETScBlockDouble::clear() void linearSystemPETScBlockDouble::clear()
{ {
if(_isAllocated){ if(_isAllocated){
#if (PETSC_VERSION_RELEASE == 0 || ((PETSC_VERSION_MAJOR == 3) && (PETSC_VERSION_MINOR == 2))) // petsc-dev #if (PETSC_VERSION_RELEASE == 0 || ((PETSC_VERSION_MAJOR == 3) && (PETSC_VERSION_MINOR >= 2)))
MatDestroy(&_a); MatDestroy(&_a);
VecDestroy(&_x); VecDestroy(&_x);
VecDestroy(&_b); VecDestroy(&_b);
......
...@@ -3,7 +3,7 @@ ...@@ -3,7 +3,7 @@
#include <petscksp.h> #include <petscksp.h>
#include "linearSystemPETSc.h" #include "linearSystemPETSc.h"
#if (PETSC_VERSION_RELEASE == 0 || ((PETSC_VERSION_MAJOR == 3) && (PETSC_VERSION_MINOR == 2))) // petsc-dev #if (PETSC_VERSION_RELEASE == 0 || ((PETSC_VERSION_MAJOR == 3) && (PETSC_VERSION_MINOR >= 2)))
#define PetscTruth PetscBool #define PetscTruth PetscBool
#define PetscOptionsGetTruth PetscOptionsGetBool #define PetscOptionsGetTruth PetscOptionsGetBool
#endif #endif
...@@ -47,7 +47,7 @@ linearSystemPETSc<scalar>::~linearSystemPETSc() ...@@ -47,7 +47,7 @@ linearSystemPETSc<scalar>::~linearSystemPETSc()
{ {
clear(); clear();
if(_kspAllocated) if(_kspAllocated)
#if (PETSC_VERSION_RELEASE == 0 || ((PETSC_VERSION_MAJOR == 3) && (PETSC_VERSION_MINOR == 2))) // petsc-dev #if (PETSC_VERSION_RELEASE == 0 || ((PETSC_VERSION_MAJOR == 3) && (PETSC_VERSION_MINOR >= 2)))
_try(KSPDestroy(&_ksp)); _try(KSPDestroy(&_ksp));
#else #else
_try(KSPDestroy(_ksp)); _try(KSPDestroy(_ksp));
...@@ -152,7 +152,7 @@ template <class scalar> ...@@ -152,7 +152,7 @@ template <class scalar>
void linearSystemPETSc<scalar>::clear() void linearSystemPETSc<scalar>::clear()
{ {
if(_isAllocated){ if(_isAllocated){
#if (PETSC_VERSION_RELEASE == 0 || ((PETSC_VERSION_MAJOR == 3) && (PETSC_VERSION_MINOR == 2))) // petsc-dev #if (PETSC_VERSION_RELEASE == 0 || ((PETSC_VERSION_MAJOR == 3) && (PETSC_VERSION_MINOR >= 2)))
_try(MatDestroy(&_a)); _try(MatDestroy(&_a));
_try(VecDestroy(&_x)); _try(VecDestroy(&_x));
_try(VecDestroy(&_b)); _try(VecDestroy(&_b));
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment