1 #if !defined(PETSC_P4EST_PACKAGE_H) 2 #define PETSC_P4EST_PACKAGE_H 3 #include <petscsys.h> 4 #if defined(PETSC_HAVE_MPIUNI) 5 #undef MPI_SUCCESS 6 #endif 7 #include <p4est_base.h> 8 #if defined(PETSC_HAVE_MPIUNI) 9 #define MPI_SUCCESS 0 10 #endif 11 12 #if defined(PETSC_HAVE_SETJMP_H) && defined(PETSC_USE_DEBUG) 13 #include <setjmp.h> 14 PETSC_INTERN jmp_buf PetscScJumpBuf; 15 16 #define PetscStackCallP4est(func,args) do { \ 17 if (setjmp(PetscScJumpBuf)) { \ 18 return PetscError(PETSC_COMM_SELF,__LINE__,PETSC_FUNCTION_NAME,__FILE__,PETSC_ERR_LIB,PETSC_ERROR_REPEAT,"Error in p4est/libsc call %s()",#func); \ 19 } \ 20 else { \ 21 PetscStackPush(#func); \ 22 func args; \ 23 PetscStackPop; \ 24 } \ 25 } while (0) 26 #define PetscStackCallP4estReturn(ret,func,args) do { \ 27 if (setjmp(PetscScJumpBuf)) { \ 28 return PetscError(PETSC_COMM_SELF,__LINE__,PETSC_FUNCTION_NAME,__FILE__,PETSC_ERR_LIB,PETSC_ERROR_REPEAT,"Error in p4est/libsc call %s()",#func); \ 29 } \ 30 else { \ 31 PetscStackPush(#func); \ 32 ret = func args; \ 33 PetscStackPop; \ 34 } \ 35 } while (0) 36 #else 37 #define PetscStackCallP4est(func,args) do { \ 38 PetscStackPush(#func); \ 39 func args; \ 40 PetscStackPop; \ 41 } while (0) 42 #define PetscStackCallP4estReturn(ret,func,args) do { \ 43 PetscStackPush(#func); \ 44 ret = func args; \ 45 PetscStackPop; \ 46 } while (0) 47 #endif 48 49 #if defined(P4EST_ENABLE_DEBUG) 50 #define PETSC_P4EST_ASSERT(x) P4EST_ASSERT(x) 51 #else 52 #define PETSC_P4EST_ASSERT(x) (void)(x) 53 #endif 54 55 PETSC_EXTERN PetscErrorCode PetscP4estInitialize(); 56 57 PETSC_STATIC_INLINE PetscErrorCode P4estLocidxCast(PetscInt a,p4est_locidx_t *b) 58 { 59 PetscFunctionBegin; 60 *b = (p4est_locidx_t)(a); 61 #if defined(PETSC_USE_64BIT_INDICES) 62 if ((a) > P4EST_LOCIDX_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index to large for p4est_locidx_t"); 63 #endif 64 PetscFunctionReturn(0); 65 } 66 67 PETSC_STATIC_INLINE PetscErrorCode P4estTopidxCast(PetscInt a,p4est_topidx_t *b) 68 { 69 PetscFunctionBegin; 70 *b = (p4est_topidx_t)(a); 71 #if defined(PETSC_USE_64BIT_INDICES) 72 if ((a) > P4EST_TOPIDX_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index to large for p4est_topidx_t"); 73 #endif 74 PetscFunctionReturn(0); 75 } 76 77 #endif 78