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 PetscStackPush(#func); \ 33 func args; \ 34 PetscStackPop; \ 35 } while (0) 36 #define PetscStackCallP4estReturn(ret,func,args) do { \ 37 PetscStackPush(#func); \ 38 ret = func args; \ 39 PetscStackPop; \ 40 } while (0) 41 #endif 42 43 PETSC_EXTERN PetscErrorCode PetscP4estInitialize(); 44 45 #undef __FUNCT__ 46 #define __FUNCT__ "P4estLocidxCast" 47 PETSC_STATIC_INLINE PetscErrorCode P4estLocidxCast(PetscInt a,p4est_locidx_t *b) 48 { 49 PetscFunctionBegin; 50 *b = (p4est_locidx_t)(a); 51 #if defined(PETSC_USE_64BIT_INDICES) 52 if ((a) > P4EST_LOCIDX_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index to large for p4est_locidx_t"); 53 #endif 54 PetscFunctionReturn(0); 55 } 56 57 #undef __FUNCT__ 58 #define __FUNCT__ "P4estTopidxCast" 59 PETSC_STATIC_INLINE PetscErrorCode P4estTopidxCast(PetscInt a,p4est_topidx_t *b) 60 { 61 PetscFunctionBegin; 62 *b = (p4est_topidx_t)(a); 63 #if defined(PETSC_USE_64BIT_INDICES) 64 if ((a) > P4EST_TOPIDX_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index to large for p4est_topidx_t"); 65 #endif 66 PetscFunctionReturn(0); 67 } 68 69 #endif 70