18 #ifndef LIBMESH_PETSC_MACRO_H
19 #define LIBMESH_PETSC_MACRO_H
22 #include "libmesh/libmesh_config.h"
24 #ifdef LIBMESH_HAVE_PETSC
33 #define PETSC_VERSION_LESS_THAN(major,minor,subminor) \
34 ((LIBMESH_DETECTED_PETSC_VERSION_MAJOR < (major) || \
35 (LIBMESH_DETECTED_PETSC_VERSION_MAJOR == (major) && (LIBMESH_DETECTED_PETSC_VERSION_MINOR < (minor) || \
36 (LIBMESH_DETECTED_PETSC_VERSION_MINOR == (minor) && \
37 LIBMESH_DETECTED_PETSC_VERSION_SUBMINOR < (subminor))))) ? 1 : 0)
41 #ifdef LIBMESH_DETECTED_PETSC_VERSION_RELEASE
43 #define PETSC_RELEASE_LESS_THAN(major,minor,subminor) \
44 (PETSC_VERSION_LESS_THAN(major,minor,subminor) && LIBMESH_DETECTED_PETSC_VERSION_RELEASE)
48 #define PETSC_RELEASE_LESS_THAN(major,minor,subminor) \
49 (PETSC_VERSION_LESS_THAN(major,minor,subminor))
57 #define EXTERN_C_FOR_PETSC_BEGIN
58 #define EXTERN_C_FOR_PETSC_END
62 #include <libmesh/ignore_warnings.h>
64 # define LIBMESH_SAW_I
68 # undef I // Avoid complex.h contamination
70 #include <libmesh/restore_warnings.h>
72 #if PETSC_RELEASE_LESS_THAN(3,1,1)
76 #if PETSC_RELEASE_LESS_THAN(3,1,1)
77 # define LibMeshVecDestroy(x) VecDestroy(*(x))
78 # define LibMeshVecScatterDestroy(x) VecScatterDestroy(*(x))
79 # define LibMeshMatDestroy(x) MatDestroy(*(x))
80 # define LibMeshISDestroy(x) ISDestroy(*(x))
81 # define LibMeshKSPDestroy(x) KSPDestroy(*(x))
82 # define LibMeshSNESDestroy(x) SNESDestroy(*(x))
83 # define LibMeshPetscViewerDestroy(x) PetscViewerDestroy(*(x))
84 # define LibMeshPCDestroy(x) PCDestroy(*(x))
86 # define LibMeshVecDestroy(x) VecDestroy(x)
87 # define LibMeshVecScatterDestroy(x) VecScatterDestroy(x)
88 # define LibMeshMatDestroy(x) MatDestroy(x)
89 # define LibMeshISDestroy(x) ISDestroy(x)
90 # define LibMeshKSPDestroy(x) KSPDestroy(x)
91 # define LibMeshSNESDestroy(x) SNESDestroy(x)
92 # define LibMeshPetscViewerDestroy(x) PetscViewerDestroy(x)
93 # define LibMeshPCDestroy(x) PCDestroy(x)
100 #define LibMeshVecScatterCreate(xin,ix,yin,iy,newctx) VecScatterCreate(xin,ix,yin,iy,newctx)
102 #if PETSC_RELEASE_LESS_THAN(3,1,1)
104 # define ISCreateLibMesh(comm,n,idx,mode,is) \
105 ((mode) == PETSC_USE_POINTER \
106 ? ISCreateGeneralWithArray((comm),(n),(idx),(is)) \
107 : ((mode) == PETSC_OWN_POINTER \
108 ? ISCreateGeneralNC((comm),(n),(idx),(is)) \
109 : ISCreateGeneral((comm),(n),(idx),(is))))
111 # define ISCreateLibMesh(comm,n,idx,mode,is) ISCreateGeneral((comm),(n),(idx),(mode),(is))
115 #if PETSC_RELEASE_LESS_THAN(3,8,0)
116 # define LibMeshCreateSubMatrix MatGetSubMatrix
118 # define LibMeshCreateSubMatrix MatCreateSubMatrix
124 #if LIBMESH_DEFAULT_QUADRUPLE_PRECISION
125 # include <boost/multiprecision/float128.hpp>
130 std::ostream &
operator<< (std::ostream & os,
const PetscScalar in)
132 os << (boost::multiprecision::float128(in));
136 #define LIBMESH_PETSCSCALAR_UNARY(funcname) \
137 inline PetscScalar funcname \
138 (const PetscScalar in) \
140 return boost::multiprecision::funcname \
141 (boost::multiprecision::float128(in)).backend().value(); \
144 LIBMESH_PETSCSCALAR_UNARY(
sqrt)
145 LIBMESH_PETSCSCALAR_UNARY(exp)
146 LIBMESH_PETSCSCALAR_UNARY(log)
147 LIBMESH_PETSCSCALAR_UNARY(log10)
148 LIBMESH_PETSCSCALAR_UNARY(sin)
149 LIBMESH_PETSCSCALAR_UNARY(cos)
150 LIBMESH_PETSCSCALAR_UNARY(tan)
151 LIBMESH_PETSCSCALAR_UNARY(asin)
152 LIBMESH_PETSCSCALAR_UNARY(acos)
153 LIBMESH_PETSCSCALAR_UNARY(atan)
154 LIBMESH_PETSCSCALAR_UNARY(sinh)
155 LIBMESH_PETSCSCALAR_UNARY(cosh)
156 LIBMESH_PETSCSCALAR_UNARY(tanh)
157 LIBMESH_PETSCSCALAR_UNARY(
abs)
158 LIBMESH_PETSCSCALAR_UNARY(fabs)
159 LIBMESH_PETSCSCALAR_UNARY(ceil)
160 LIBMESH_PETSCSCALAR_UNARY(floor)
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_RELEASE_LESS_THAN(major,minor,subminor) 1
240 #endif // LIBMESH_HAVE_PETSC
242 #endif // LIBMESH_PETSC_MACRO_H