Actual source code: petsc_p4est_package.h

  1: #ifndef 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

 12: #if defined(PETSC_HAVE_SETJMP_H) && defined(PETSC_USE_DEBUG)
 13:   #include <setjmp.h>
 14: PETSC_INTERN jmp_buf PetscScJumpBuf;

 16:   #define PetscCallP4est(func, args) \
 17:     do { \
 18:       if (setjmp(PetscScJumpBuf)) { \
 19:         return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_LIB, PETSC_ERROR_REPEAT, "Error in p4est/libsc call %s()", #func); \
 20:       } else { \
 21:         PetscStackPushExternal(#func); \
 22:         func args; \
 23:         PetscStackPop; \
 24:       } \
 25:     } while (0)
 26:   #define PetscCallP4estReturn(ret, func, args) \
 27:     do { \
 28:       if (setjmp(PetscScJumpBuf)) { \
 29:         return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_LIB, PETSC_ERROR_REPEAT, "Error in p4est/libsc call %s()", #func); \
 30:       } else { \
 31:         PetscStackPushExternal(#func); \
 32:         ret = func args; \
 33:         PetscStackPop; \
 34:       } \
 35:     } while (0)
 36: #else
 37:   #define PetscCallP4est(func, args) \
 38:     do { \
 39:       PetscStackPushExternal(#func); \
 40:       func args; \
 41:       PetscStackPop; \
 42:     } while (0)
 43:   #define PetscCallP4estReturn(ret, func, args) \
 44:     do { \
 45:       PetscStackPushExternal(#func); \
 46:       ret = func args; \
 47:       PetscStackPop; \
 48:     } while (0)
 49: #endif

 51: #if defined(P4EST_ENABLE_DEBUG)
 52:   #define PETSC_P4EST_ASSERT(x) P4EST_ASSERT(x)
 53: #else
 54:   #define PETSC_P4EST_ASSERT(x) (void)(x)
 55: #endif

 57: PETSC_EXTERN PetscErrorCode PetscP4estInitialize(void);

 59: static inline PetscErrorCode P4estLocidxCast(PetscInt a, p4est_locidx_t *b)
 60: {
 61:   *b = (p4est_locidx_t)(a);
 62: #if defined(PETSC_USE_64BIT_INDICES)
 64: #endif
 65:   return 0;
 66: }

 68: static inline PetscErrorCode P4estTopidxCast(PetscInt a, p4est_topidx_t *b)
 69: {
 70:   *b = (p4est_topidx_t)(a);
 71: #if defined(PETSC_USE_64BIT_INDICES)
 73: #endif
 74:   return 0;
 75: }

 77: #endif