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