18 #ifndef LIBMESH_PETSC_MACRO_H    19 #define LIBMESH_PETSC_MACRO_H    22 #include "libmesh/libmesh_config.h"    24 #ifdef LIBMESH_HAVE_PETSC    26 #include <petscversion.h>    30 #define PETSC_VERSION_LESS_THAN(major,minor,subminor) \    31   PETSC_VERSION_LT(major,minor,subminor)    33 #define PETSC_VERSION_EQUALS(major,minor,subminor) \    34   PETSC_VERSION_EQ(major,minor,subminor)    40 #define EXTERN_C_FOR_PETSC_BEGIN    41 #define EXTERN_C_FOR_PETSC_END    45 #include <libmesh/ignore_warnings.h>    47 # define LIBMESH_SAW_I    51 # undef I // Avoid complex.h contamination    53 #include <libmesh/restore_warnings.h>    57 #define LibMeshVecDestroy(x)         VecDestroy(x)    58 #define LibMeshVecScatterDestroy(x)  VecScatterDestroy(x)    59 #define LibMeshMatDestroy(x)         MatDestroy(x)    60 #define LibMeshISDestroy(x)          ISDestroy(x)    61 #define LibMeshKSPDestroy(x)         KSPDestroy(x)    62 #define LibMeshSNESDestroy(x)        SNESDestroy(x)    63 #define LibMeshPetscViewerDestroy(x) PetscViewerDestroy(x)    64 #define LibMeshPCDestroy(x)          PCDestroy(x)    70 #define LibMeshVecScatterCreate(xin,ix,yin,iy,newctx)  VecScatterCreate(xin,ix,yin,iy,newctx)    71 #define ISCreateLibMesh(comm,n,idx,mode,is) ISCreateGeneral((comm),(n),(idx),(mode),(is))    74 #if PETSC_VERSION_LESS_THAN(3,8,0)    75 # define LibMeshCreateSubMatrix MatGetSubMatrix    77 # define LibMeshCreateSubMatrix MatCreateSubMatrix    82 #if PETSC_VERSION_LESS_THAN(3,17,0)    83 # define LIBMESH_SETERRQ1 SETERRQ1    84 # define LIBMESH_SETERRQ2 SETERRQ2    85 # define LIBMESH_SETERRQ3 SETERRQ3    87 # define LIBMESH_SETERRQ1 SETERRQ    88 # define LIBMESH_SETERRQ2 SETERRQ    89 # define LIBMESH_SETERRQ3 SETERRQ    94 #if PETSC_VERSION_LESS_THAN(3,8,0)    95 # define LIBMESH_PETSCINT_FMT "D"    97 # define LIBMESH_PETSCINT_FMT PetscInt_FMT   101 #if PETSC_VERSION_LESS_THAN(3,19,0)   102 # define LIBMESH_PETSC_NULLPTR PETSC_NULL   103 # define LIBMESH_PETSC_SUCCESS static_cast<PetscErrorCode>(0)   105 # define LIBMESH_PETSC_NULLPTR PETSC_NULLPTR   106 # define LIBMESH_PETSC_SUCCESS PETSC_SUCCESS   112 #if LIBMESH_DEFAULT_QUADRUPLE_PRECISION   113 # include <boost/multiprecision/float128.hpp>   122 #if BOOST_VERSION < 107000   124 std::ostream & 
operator<< (std::ostream & os, 
const PetscScalar in)
   126   os << (boost::multiprecision::float128(in));
   167 template <
typename T>
   168 PetscScalar 
PS(T val)
   170   return val.backend().value();
   173 template <
typename T>
   174 PetscScalar * 
pPS(T * ptr)
   176   return &(ptr->backend().value());
   179 template <
typename T>
   180 const PetscScalar * 
pPS(
const T * ptr)
   182   return &(ptr->backend().value());
   185 template <
typename T>
   188   return &(ptr->backend().value());
   191 template <
typename T>
   192 const PetscReal * 
pPR(
const T * ptr)
   194   return &(ptr->backend().value());
   202 template <
typename T>
   203 PetscScalar 
PS(T val)
   208 template <
typename T>
   209 PetscScalar * 
pPS(T * ptr)
   214 template <
typename T>
   215 const PetscScalar * 
pPS(
const T * ptr)
   220 template <
typename T>
   221 PetscReal * 
pPR(T * ptr)
   226 template <
typename T>
   227 const PetscReal * 
pPR(
const T * ptr)
   233 #endif // LIBMESH_ENABLE_QUADRUPLE_PRECISION   235 #else // LIBMESH_HAVE_PETSC   237 #define PETSC_VERSION_LESS_THAN(major,minor,subminor) 1   238 #define PETSC_VERSION_EQUALS(major,minor,revision) 0   240 #endif // LIBMESH_HAVE_PETSC   242 #define PETSC_RELEASE_LESS_THAN(major, minor, subminor) \   243   PETSC_VERSION_LESS_THAN(major, minor, subminor)   244 #define PETSC_RELEASE_EQUALS(major, minor, subminor) \   245   PETSC_VERSION_EQUALS(major, minor, subminor)   247 #define PETSC_RELEASE_LESS_EQUALS(major, minor, subminor) \   248   (PETSC_RELEASE_LESS_THAN(major, minor, subminor) || PETSC_RELEASE_EQUALS(major, minor, subminor))   250 #define PETSC_RELEASE_GREATER_EQUALS(major, minor, subminor) \   251   (0 == PETSC_RELEASE_LESS_THAN(major, minor, subminor))   253 #define PETSC_RELEASE_GREATER_THAN(major, minor, subminor) \   254   (0 == PETSC_RELEASE_LESS_EQUALS(major, minor, subminor))   256 #endif // LIBMESH_PETSC_MACRO_H 
std::ostream & operator<<(std::ostream &os, const PetscScalar in)
PetscScalar * pPS(T *ptr)
The libMesh namespace provides an interface to certain functionality in the library.