diff options
author | schnetter <schnetter@1d96b42b-98df-4a6a-9d84-1b24288d4588> | 2006-03-06 19:43:02 +0000 |
---|---|---|
committer | schnetter <schnetter@1d96b42b-98df-4a6a-9d84-1b24288d4588> | 2006-03-06 19:43:02 +0000 |
commit | c14d14eca019ab32134593325466c734b2d87b3c (patch) | |
tree | 6d18ad5d7a504d703942b529b1ea865de88cb98d /src/petsc_defines.h | |
parent | 24ad3e9442606dfe2e3aeaf3778a52d1673e2dc3 (diff) |
Clean up the backwards compatibility macros for PETSc.
Incorporate Steve's changes for PETSc 2.2.x.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusElliptic/EllPETSc/trunk@94 1d96b42b-98df-4a6a-9d84-1b24288d4588
Diffstat (limited to 'src/petsc_defines.h')
-rw-r--r-- | src/petsc_defines.h | 63 |
1 files changed, 63 insertions, 0 deletions
diff --git a/src/petsc_defines.h b/src/petsc_defines.h new file mode 100644 index 0000000..19c321e --- /dev/null +++ b/src/petsc_defines.h @@ -0,0 +1,63 @@ +#ifndef PETSC_DEFINES_H +#define PETSC_DEFINES_H + +#include "cctk.h" + +#include <petsc.h> +#include <petscversion.h> + + + +/* The PETSc API changes frequently, sometimes even in minor releases. + The #defines below map calls to the PETSc linear solver to the + corresponding PETSc call (or sequence of calls), depending on the + version of PETSc. */ + + + +#if PETSC_VERSION_MAJOR < 2 || \ + (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 0) + +/* Up to and including version 2.0.x */ +# include <sles.h> + +#elif PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR <= 1 + +/* Up to and including version 2.1.x */ +# include <petscsles.h> + +#elif PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2 + +/* Version 2.2.x */ +# include <petscmat.h> +# include <petscksp.h> +# define SLES KSP +# define SLESCreate KSPCreate +# define SLESSetOperators KSPSetOperators +# define SLESGetKSP(a,b) ((*b) = (a), 0) +# define SLESGetPC KSPGetPC +# define SLESSetFromOptions KSPSetFromOptions +# define SLESSolve(a,b,c,d) (KSPSetRhs(a,b), \ + KSPSetSolution(a,c), \ + KSPSolve(a), \ + KSPGetIterationNumber(a,d)) +# define SLESDestroy KSPDestroy + +#else + +/* Later versions */ +# include <petscmat.h> +# include <petscksp.h> +# define SLES KSP +# define SLESCreate KSPCreate +# define SLESSetOperators KSPSetOperators +# define SLESGetKSP(a,b) ((*b) = (a), 0) +# define SLESGetPC KSPGetPC +# define SLESSetFromOptions KSPSetFromOptions +# define SLESSolve(a,b,c,d) (KSPSolve(a,b,c), \ + KSPGetIterationNumber(a,d)) +# define SLESDestroy KSPDestroy + +#endif + +#endif /* #ifndef PETSC_DEFINES_H */ |