Before 3.2 Petsc used PetscTruth; afterwards it used PetscTruth.
namespace PETScWrappers
{
+#if DEAL_II_PETSC_VERSION_LT(3,2,0)
+ typedef PetscTruth PetscBooleanType;
+#else
+ typedef PetscBool PetscBooleanType;
+#endif
+
/**
* Set an option in the global PETSc database. This function just wraps
* PetscOptionsSetValue with a version check (the signature of this function
#ifdef DEAL_II_WITH_PETSC
# include <deal.II/base/subscriptor.h>
-# include <deal.II/lac/full_matrix.h>
# include <deal.II/lac/exceptions.h>
+# include <deal.II/lac/full_matrix.h>
+# include <deal.II/lac/petsc_compatibility.h>
# include <deal.II/lac/vector.h>
# include <petscmat.h>
* Test whether a matrix is symmetric. Default tolerance is
* $1000\times32$-bit machine precision.
*/
-#if DEAL_II_PETSC_VERSION_LT(3,2,0)
- PetscTruth
-#else
- PetscBool
-#endif
+ PetscBooleanType
is_symmetric (const double tolerance = 1.e-12);
/**
* its transpose. Default tolerance is $1000\times32$-bit machine
* precision.
*/
-#if DEAL_II_PETSC_VERSION_LT(3,2,0)
- PetscTruth
-#else
- PetscBool
-#endif
+ PetscBooleanType
is_hermitian (const double tolerance = 1.e-12);
/**
#ifdef DEAL_II_WITH_PETSC
# include <deal.II/lac/exceptions.h>
+# include <deal.II/lac/petsc_compatibility.h>
# include <deal.II/lac/petsc_full_matrix.h>
# include <deal.II/lac/petsc_sparse_matrix.h>
# include <deal.II/lac/petsc_parallel_sparse_matrix.h>
AssertThrow (ierr == 0, ExcPETScError(ierr));
}
-#if DEAL_II_PETSC_VERSION_LT(3,2,0)
- PetscTruth
-#else
- PetscBool
-#endif
+ PetscBooleanType
MatrixBase::is_symmetric (const double tolerance)
{
-#if DEAL_II_PETSC_VERSION_LT(3,2,0)
- PetscTruth
-#else
- PetscBool
-#endif
- truth;
+ PetscBooleanType truth;
assert_is_compressed ();
int ierr = MatIsSymmetric (matrix, tolerance, &truth);
(void)ierr;
return truth;
}
-#if DEAL_II_PETSC_VERSION_LT(3,2,0)
- PetscTruth
-#else
- PetscBool
-#endif
+ PetscBooleanType
MatrixBase::is_hermitian (const double tolerance)
{
-#if DEAL_II_PETSC_VERSION_LT(3,2,0)
- PetscTruth
-#else
- PetscBool
-#endif
- truth;
+ PetscBooleanType truth;
assert_is_compressed ();
int ierr = MatIsHermitian (matrix, tolerance, &truth);
# include <deal.II/base/memory_consumption.h>
# include <deal.II/lac/exceptions.h>
+# include <deal.II/lac/petsc_compatibility.h>
# include <deal.II/lac/petsc_vector.h>
# include <deal.II/lac/petsc_parallel_vector.h>
# include <cmath>
Assert (size() == v.size(),
ExcDimensionMismatch(size(), v.size()));
-#if DEAL_II_PETSC_VERSION_LT(3,2,0)
- PetscTruth
-#else
- PetscBool
-#endif
- flag;
-
+ PetscBooleanType flag;
const int ierr = VecEqual (vector, v.vector, &flag);
AssertThrow (ierr == 0, ExcPETScError(ierr));
Assert (size() == v.size(),
ExcDimensionMismatch(size(), v.size()));
-#if DEAL_II_PETSC_VERSION_LT(3,2,0)
- PetscTruth
-#else
- PetscBool
-#endif
- flag;
-
+ PetscBooleanType flag;
const int ierr = VecEqual (vector, v.vector, &flag);
AssertThrow (ierr == 0, ExcPETScError(ierr));