1*a4963045SJacob Faibussowitsch #pragma once
2cd375237SToby Isaac #include <petscsys.h>
3ed587871SSatish Balay #if defined(PETSC_HAVE_MPIUNI)
4ed587871SSatish Balay #undef MPI_SUCCESS
5ed587871SSatish Balay #endif
65bc33bd4SToby Isaac #include <p4est_base.h>
7ed587871SSatish Balay #if defined(PETSC_HAVE_MPIUNI)
8ed587871SSatish Balay #define MPI_SUCCESS 0
9ed587871SSatish Balay #endif
10cd375237SToby Isaac
118f3f2bf1SBarry Smith #if defined(PETSC_HAVE_SETJMP_H) && defined(PETSC_USE_DEBUG)
12cd375237SToby Isaac #include <setjmp.h>
13cd375237SToby Isaac PETSC_INTERN jmp_buf PetscScJumpBuf;
14cd375237SToby Isaac
159371c9d4SSatish Balay #define PetscCallP4est(func, args) \
169371c9d4SSatish Balay do { \
17cd375237SToby Isaac if (setjmp(PetscScJumpBuf)) { \
18bdeda3f1SToby Isaac return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_LIB, PETSC_ERROR_REPEAT, "Error in p4est/libsc call %s()", #func); \
199371c9d4SSatish Balay } else { \
20792fecdfSBarry Smith PetscStackPushExternal(#func); \
21cd375237SToby Isaac func args; \
22cd375237SToby Isaac PetscStackPop; \
23cd375237SToby Isaac } \
24cd375237SToby Isaac } while (0)
259371c9d4SSatish Balay #define PetscCallP4estReturn(ret, func, args) \
269371c9d4SSatish Balay do { \
27cd375237SToby Isaac if (setjmp(PetscScJumpBuf)) { \
28bdeda3f1SToby Isaac return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_LIB, PETSC_ERROR_REPEAT, "Error in p4est/libsc call %s()", #func); \
299371c9d4SSatish Balay } else { \
30792fecdfSBarry Smith PetscStackPushExternal(#func); \
31cd375237SToby Isaac ret = func args; \
32cd375237SToby Isaac PetscStackPop; \
33cd375237SToby Isaac } \
34cd375237SToby Isaac } while (0)
35cd375237SToby Isaac #else
369371c9d4SSatish Balay #define PetscCallP4est(func, args) \
379371c9d4SSatish Balay do { \
38792fecdfSBarry Smith PetscStackPushExternal(#func); \
39cd375237SToby Isaac func args; \
40cd375237SToby Isaac PetscStackPop; \
41cd375237SToby Isaac } while (0)
429371c9d4SSatish Balay #define PetscCallP4estReturn(ret, func, args) \
439371c9d4SSatish Balay do { \
44792fecdfSBarry Smith PetscStackPushExternal(#func); \
45cd375237SToby Isaac ret = func args; \
46cd375237SToby Isaac PetscStackPop; \
47cd375237SToby Isaac } while (0)
48cd375237SToby Isaac #endif
49cd375237SToby Isaac
50be929b80SToby Isaac #if defined(P4EST_ENABLE_DEBUG)
51be929b80SToby Isaac #define PETSC_P4EST_ASSERT(x) P4EST_ASSERT(x)
52be929b80SToby Isaac #else
53be929b80SToby Isaac #define PETSC_P4EST_ASSERT(x) (void)(x)
54be929b80SToby Isaac #endif
55be929b80SToby Isaac
563274405dSPierre Jolivet PETSC_EXTERN PetscErrorCode PetscP4estInitialize(void);
57cd375237SToby Isaac
P4estLocidxCast(PetscInt a,p4est_locidx_t * b)58d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode P4estLocidxCast(PetscInt a, p4est_locidx_t *b)
59d71ae5a4SJacob Faibussowitsch {
605bc33bd4SToby Isaac PetscFunctionBegin;
61438c3a38SToby Isaac *b = (p4est_locidx_t)(a);
625bc33bd4SToby Isaac #if defined(PETSC_USE_64BIT_INDICES)
6308401ef6SPierre Jolivet PetscCheck((a) <= P4EST_LOCIDX_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index to large for p4est_locidx_t");
645bc33bd4SToby Isaac #endif
653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
665bc33bd4SToby Isaac }
675bc33bd4SToby Isaac
P4estTopidxCast(PetscInt a,p4est_topidx_t * b)68d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode P4estTopidxCast(PetscInt a, p4est_topidx_t *b)
69d71ae5a4SJacob Faibussowitsch {
705bc33bd4SToby Isaac PetscFunctionBegin;
71438c3a38SToby Isaac *b = (p4est_topidx_t)(a);
725bc33bd4SToby Isaac #if defined(PETSC_USE_64BIT_INDICES)
7308401ef6SPierre Jolivet PetscCheck((a) <= P4EST_TOPIDX_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index to large for p4est_topidx_t");
745bc33bd4SToby Isaac #endif
753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
765bc33bd4SToby Isaac }
77