Actual source code: petsc_p4est_package.h
1: #pragma once
2: #include <petscsys.h>
3: #if defined(PETSC_HAVE_MPIUNI)
4: #undef MPI_SUCCESS
5: #endif
6: #include <p4est_base.h>
7: #if defined(PETSC_HAVE_MPIUNI)
8: #define MPI_SUCCESS 0
9: #endif
11: #if defined(PETSC_HAVE_SETJMP_H) && defined(PETSC_USE_DEBUG)
12: #include <setjmp.h>
13: PETSC_INTERN jmp_buf PetscScJumpBuf;
15: #define PetscCallP4est(func, args) \
16: 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: } else { \
20: PetscStackPushExternal(#func); \
21: func args; \
22: PetscStackPop; \
23: } \
24: } while (0)
25: #define PetscCallP4estReturn(ret, func, args) \
26: 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: } else { \
30: PetscStackPushExternal(#func); \
31: ret = func args; \
32: PetscStackPop; \
33: } \
34: } while (0)
35: #else
36: #define PetscCallP4est(func, args) \
37: do { \
38: PetscStackPushExternal(#func); \
39: func args; \
40: PetscStackPop; \
41: } while (0)
42: #define PetscCallP4estReturn(ret, func, args) \
43: do { \
44: PetscStackPushExternal(#func); \
45: ret = func args; \
46: PetscStackPop; \
47: } while (0)
48: #endif
50: #if defined(P4EST_ENABLE_DEBUG)
51: #define PETSC_P4EST_ASSERT(x) P4EST_ASSERT(x)
52: #else
53: #define PETSC_P4EST_ASSERT(x) (void)(x)
54: #endif
56: PETSC_EXTERN PetscErrorCode PetscP4estInitialize(void);
58: static inline PetscErrorCode P4estLocidxCast(PetscInt a, p4est_locidx_t *b)
59: {
60: PetscFunctionBegin;
61: *b = (p4est_locidx_t)(a);
62: #if defined(PETSC_USE_64BIT_INDICES)
63: PetscCheck((a) <= P4EST_LOCIDX_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index to large for p4est_locidx_t");
64: #endif
65: PetscFunctionReturn(PETSC_SUCCESS);
66: }
68: static inline PetscErrorCode P4estTopidxCast(PetscInt a, p4est_topidx_t *b)
69: {
70: PetscFunctionBegin;
71: *b = (p4est_topidx_t)(a);
72: #if defined(PETSC_USE_64BIT_INDICES)
73: PetscCheck((a) <= P4EST_TOPIDX_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index to large for p4est_topidx_t");
74: #endif
75: PetscFunctionReturn(PETSC_SUCCESS);
76: }