154a8ef01SBarry Smith /* 2f621e05eSBarry Smith Contains all error handling interfaces for PETSc. 354a8ef01SBarry Smith */ 46524c165SJacob Faibussowitsch #ifndef PETSCERROR_H 526bd1501SBarry Smith #define PETSCERROR_H 66c7e564aSBarry Smith 7e1bf4ed2SJacob Faibussowitsch #include <petscmacros.h> 8e1bf4ed2SJacob Faibussowitsch #include <petscsystypes.h> 9e1bf4ed2SJacob Faibussowitsch 10ac09b921SBarry Smith /* SUBMANSEC = Sys */ 11ac09b921SBarry Smith 1254a8ef01SBarry Smith /* 13329ffe3dSLois Curfman McInnes These are the generic error codes. These error codes are used 14e2d1d2b7SBarry Smith many different places in the PETSc source code. The string versions are 150e5e90baSSatish Balay at src/sys/error/err.c any changes here must also be made there 16e951c290SBarry Smith These are also define in src/sys/f90-mod/petscerror.h any CHANGES here 170f9cf654SBarry Smith must be also made there. 1845d48df9SBarry Smith 1954a8ef01SBarry Smith */ 202a6744ebSBarry Smith #define PETSC_ERR_MIN_VALUE 54 /* should always be one less then the smallest value */ 212a6744ebSBarry Smith 2245d48df9SBarry Smith #define PETSC_ERR_MEM 55 /* unable to allocate requested memory */ 2347794344SBarry Smith #define PETSC_ERR_SUP 56 /* no support for requested operation */ 24e2d1d2b7SBarry Smith #define PETSC_ERR_SUP_SYS 57 /* no support for requested operation on this computer system */ 25e2d1d2b7SBarry Smith #define PETSC_ERR_ORDER 58 /* operation done in wrong order */ 2645d48df9SBarry Smith #define PETSC_ERR_SIG 59 /* signal received */ 27f1caa5a4SBarry Smith #define PETSC_ERR_FP 72 /* floating point exception */ 28a8c6a408SBarry Smith #define PETSC_ERR_COR 74 /* corrupted PETSc object */ 29a8c6a408SBarry Smith #define PETSC_ERR_LIB 76 /* error in library called by PETSc */ 30329ffe3dSLois Curfman McInnes #define PETSC_ERR_PLIB 77 /* PETSc library generated inconsistent data */ 31329ffe3dSLois Curfman McInnes #define PETSC_ERR_MEMC 78 /* memory corruption */ 32b3cc6726SBarry Smith #define PETSC_ERR_CONV_FAILED 82 /* iterative method (KSP or SNES) failed */ 331302d50aSBarry Smith #define PETSC_ERR_USER 83 /* user has not provided needed function */ 344e2ffeddSBarry Smith #define PETSC_ERR_SYS 88 /* error in system call */ 35a8b45ee7SBarry Smith #define PETSC_ERR_POINTER 70 /* pointer does not point to valid address */ 363d96e996SBarry Smith #define PETSC_ERR_MPI_LIB_INCOMP 87 /* MPI library at runtime is not compatible with MPI user compiled with */ 3745d48df9SBarry Smith 3845d48df9SBarry Smith #define PETSC_ERR_ARG_SIZ 60 /* nonconforming object sizes used in operation */ 3945d48df9SBarry Smith #define PETSC_ERR_ARG_IDN 61 /* two arguments not allowed to be the same */ 40a8c6a408SBarry Smith #define PETSC_ERR_ARG_WRONG 62 /* wrong argument (but object probably ok) */ 4145d48df9SBarry Smith #define PETSC_ERR_ARG_CORRUPT 64 /* null or corrupted PETSc object as argument */ 4245d48df9SBarry Smith #define PETSC_ERR_ARG_OUTOFRANGE 63 /* input argument, out of range */ 434f227f7cSBarry Smith #define PETSC_ERR_ARG_BADPTR 68 /* invalid pointer argument */ 444f227f7cSBarry Smith #define PETSC_ERR_ARG_NOTSAMETYPE 69 /* two args must be same object type */ 456831982aSBarry Smith #define PETSC_ERR_ARG_NOTSAMECOMM 80 /* two args must be same communicators */ 46d252947aSBarry Smith #define PETSC_ERR_ARG_WRONGSTATE 73 /* object in argument is in wrong state, e.g. unassembled mat */ 478cda6cd7SBarry Smith #define PETSC_ERR_ARG_TYPENOTSET 89 /* the type of the object has not yet been set */ 48a8c6a408SBarry Smith #define PETSC_ERR_ARG_INCOMP 75 /* two arguments are incompatible */ 494482741eSBarry Smith #define PETSC_ERR_ARG_NULL 85 /* argument is null that should not be */ 50958c9bccSBarry Smith #define PETSC_ERR_ARG_UNKNOWN_TYPE 86 /* type name doesn't match any registered type */ 514f227f7cSBarry Smith 524f227f7cSBarry Smith #define PETSC_ERR_FILE_OPEN 65 /* unable to open file */ 534f227f7cSBarry Smith #define PETSC_ERR_FILE_READ 66 /* unable to read from file */ 544f227f7cSBarry Smith #define PETSC_ERR_FILE_WRITE 67 /* unable to write to file */ 55a8c6a408SBarry Smith #define PETSC_ERR_FILE_UNEXPECTED 79 /* unexpected data in file */ 5645d48df9SBarry Smith 57329ffe3dSLois Curfman McInnes #define PETSC_ERR_MAT_LU_ZRPVT 71 /* detected a zero pivot during LU factorization */ 589e3b2f23SBarry Smith #define PETSC_ERR_MAT_CH_ZRPVT 81 /* detected a zero pivot during Cholesky factorization */ 5954a8ef01SBarry Smith 6068e69593SBarry Smith #define PETSC_ERR_INT_OVERFLOW 84 613855c12bSBarry Smith 62bf3909cdSBarry Smith #define PETSC_ERR_FLOP_COUNT 90 63e113a28aSBarry Smith #define PETSC_ERR_NOT_CONVERGED 91 /* solver did not converge */ 6492e8f287SBarry Smith #define PETSC_ERR_MISSING_FACTOR 92 /* MatGetFactor() failed */ 65691b26d3SBarry Smith #define PETSC_ERR_OPT_OVERWRITE 93 /* attempted to over write options which should not be changed */ 66691b26d3SBarry Smith #define PETSC_ERR_WRONG_MPI_SIZE 94 /* example/application run with number of MPI ranks it does not support */ 67691b26d3SBarry Smith #define PETSC_ERR_USER_INPUT 95 /* missing or incorrect user input */ 68589f383fSBarry Smith #define PETSC_ERR_GPU_RESOURCE 96 /* unable to load a GPU resource, for example cuBLAS */ 69589f383fSBarry Smith #define PETSC_ERR_GPU 97 /* An error from a GPU call, this may be due to lack of resources on the GPU or a true error in the call */ 70ffc4695bSBarry Smith #define PETSC_ERR_MPI 98 /* general MPI error */ 71e809461dSJacob Faibussowitsch #define PETSC_ERR_RETURN 99 /* PetscError() incorrectly returned an error code of 0 */ 72e809461dSJacob Faibussowitsch #define PETSC_ERR_MAX_VALUE 100 /* this is always the one more than the largest error code */ 73b9eb5ee8SHong Zhang 7498921bdaSJacob Faibussowitsch #define SETERRQ1(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__) 7598921bdaSJacob Faibussowitsch #define SETERRQ2(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__) 7698921bdaSJacob Faibussowitsch #define SETERRQ3(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__) 7798921bdaSJacob Faibussowitsch #define SETERRQ4(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__) 7898921bdaSJacob Faibussowitsch #define SETERRQ5(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__) 7998921bdaSJacob Faibussowitsch #define SETERRQ6(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__) 8098921bdaSJacob Faibussowitsch #define SETERRQ7(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__) 8198921bdaSJacob Faibussowitsch #define SETERRQ8(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__) 8298921bdaSJacob Faibussowitsch #define SETERRQ9(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__) 8398921bdaSJacob Faibussowitsch 8430de9b25SBarry Smith /*MC 851957e957SBarry Smith SETERRQ - Macro to be called when an error has been detected, 8630de9b25SBarry Smith 8730de9b25SBarry Smith Synopsis: 88aaa7dc30SBarry Smith #include <petscsys.h> 8998921bdaSJacob Faibussowitsch PetscErrorCode SETERRQ(MPI_Comm comm,PetscErrorCode ierr,char *message,...) 9030de9b25SBarry Smith 91d083f849SBarry Smith Collective 9230de9b25SBarry Smith 9330de9b25SBarry Smith Input Parameters: 9487497f52SBarry Smith + comm - A communicator, use `PETSC_COMM_SELF` unless you know all ranks of another communicator will detect the error 953af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 9630de9b25SBarry Smith - message - error message 9730de9b25SBarry Smith 9830de9b25SBarry Smith Level: beginner 9930de9b25SBarry Smith 10030de9b25SBarry Smith Notes: 10187497f52SBarry Smith This is rarely needed, one should use `PetscCheck()` and `PetscCall()` and friends to automatically handle error conditions. 10230de9b25SBarry Smith Once the error handler is called the calling function is then returned from with the given error code. 10330de9b25SBarry Smith 10487497f52SBarry Smith Experienced users can set the error handler with `PetscPushErrorHandler()`. 10530de9b25SBarry Smith 1063b1008b8SBarry Smith Fortran Notes: 1073b1008b8SBarry Smith SETERRQ() may be called from Fortran subroutines but SETERRA() must be called from the 1083b1008b8SBarry Smith Fortran main program. 1093b1008b8SBarry Smith 110db781477SPatrick Sanan .seealso: `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 111db781477SPatrick Sanan `PetscError()`, `PetscCall()`, `CHKMEMQ`, `CHKERRA()`, `PetscCallMPI()` 11230de9b25SBarry Smith M*/ 113e809461dSJacob Faibussowitsch #define SETERRQ(comm, ierr, ...) \ 114e809461dSJacob Faibussowitsch do { \ 115e809461dSJacob Faibussowitsch PetscErrorCode ierr_seterrq_petsc_ = PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \ 116e809461dSJacob Faibussowitsch return ierr_seterrq_petsc_ ? ierr_seterrq_petsc_ : PETSC_ERR_RETURN; \ 117e809461dSJacob Faibussowitsch } while (0) 118986eef2eSBarry Smith 119ffc4695bSBarry Smith /* 120ffc4695bSBarry Smith Returned from PETSc functions that are called from MPI, such as related to attributes 121ffc4695bSBarry Smith Do not confuse PETSC_MPI_ERROR_CODE and PETSC_ERR_MPI, the first is registered with MPI and returned to MPI as 122888a1fb3SBarry Smith an error code, the latter is a regular PETSc error code passed within PETSc code indicating an error was detected in an MPI call. 123ffc4695bSBarry Smith */ 124ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CLASS; 125ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CODE; 126ffc4695bSBarry Smith 127986eef2eSBarry Smith /*MC 128986eef2eSBarry Smith SETERRMPI - Macro to be called when an error has been detected within an MPI callback function 129986eef2eSBarry Smith 130986eef2eSBarry Smith Synopsis: 131986eef2eSBarry Smith #include <petscsys.h> 13298921bdaSJacob Faibussowitsch PetscErrorCode SETERRMPI(MPI_Comm comm,PetscErrorCode ierr,char *message,...) 133986eef2eSBarry Smith 134d083f849SBarry Smith Collective 135986eef2eSBarry Smith 136986eef2eSBarry Smith Input Parameters: 13787497f52SBarry Smith + comm - A communicator, use `PETSC_COMM_SELF` unless you know all ranks of another communicator will detect the error 138986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 139986eef2eSBarry Smith - message - error message 140986eef2eSBarry Smith 141986eef2eSBarry Smith Level: developer 142986eef2eSBarry Smith 143986eef2eSBarry Smith Notes: 14487497f52SBarry Smith This macro is FOR USE IN MPI CALLBACK FUNCTIONS ONLY, such as those passed to `MPI_Comm_create_keyval()`. It always returns the error code `PETSC_MPI_ERROR_CODE` 14587497f52SBarry Smith which is registered with `MPI_Add_error_code()` when PETSc is initialized. 146986eef2eSBarry Smith 147db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscCall()`, `PetscCallMPI()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ` 148986eef2eSBarry Smith M*/ 14998921bdaSJacob Faibussowitsch #define SETERRMPI(comm, ierr, ...) return (PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__), PETSC_MPI_ERROR_CODE) 150ee8199e6SMatthew G. Knepley 151ee8199e6SMatthew G. Knepley /*MC 152f388eb8bSPatrick Sanan SETERRA - Fortran-only macro that can be called when an error has been detected from the main program 153f388eb8bSPatrick Sanan 154f388eb8bSPatrick Sanan Synopsis: 155f388eb8bSPatrick Sanan #include <petscsys.h> 156f388eb8bSPatrick Sanan PetscErrorCode SETERRA(MPI_Comm comm,PetscErrorCode ierr,char *message) 157f388eb8bSPatrick Sanan 158f388eb8bSPatrick Sanan Collective 159f388eb8bSPatrick Sanan 160f388eb8bSPatrick Sanan Input Parameters: 161f388eb8bSPatrick Sanan + comm - A communicator, so that the error can be collective 162f388eb8bSPatrick Sanan . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 163f388eb8bSPatrick Sanan - message - error message in the printf format 164f388eb8bSPatrick Sanan 165f388eb8bSPatrick Sanan Level: beginner 166f388eb8bSPatrick Sanan 167f388eb8bSPatrick Sanan Notes: 16887497f52SBarry Smith This should only be used with Fortran. With C/C++, use `SETERRQ()`. 169f388eb8bSPatrick Sanan 170f388eb8bSPatrick Sanan Fortran Notes: 17187497f52SBarry Smith `SETERRQ()` may be called from Fortran subroutines but `SETERRA()` must be called from the 172f388eb8bSPatrick Sanan Fortran main program. 173f388eb8bSPatrick Sanan 174db781477SPatrick Sanan .seealso: `SETERRQ()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()` 175f388eb8bSPatrick Sanan M*/ 176f388eb8bSPatrick Sanan 177f388eb8bSPatrick Sanan /*MC 178fa190f98SMatthew G. Knepley SETERRABORT - Macro that can be called when an error has been detected, 179fa190f98SMatthew G. Knepley 180fa190f98SMatthew G. Knepley Synopsis: 181fa190f98SMatthew G. Knepley #include <petscsys.h> 18298921bdaSJacob Faibussowitsch PetscErrorCode SETERRABORT(MPI_Comm comm,PetscErrorCode ierr,char *message,...) 183fa190f98SMatthew G. Knepley 184d083f849SBarry Smith Collective 185fa190f98SMatthew G. Knepley 186fa190f98SMatthew G. Knepley Input Parameters: 187fa190f98SMatthew G. Knepley + comm - A communicator, so that the error can be collective 1883af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 189fa190f98SMatthew G. Knepley - message - error message in the printf format 190fa190f98SMatthew G. Knepley 191fa190f98SMatthew G. Knepley Level: beginner 192fa190f98SMatthew G. Knepley 193fa190f98SMatthew G. Knepley Notes: 19487497f52SBarry Smith This function just calls `MPI_Abort()`. 19587497f52SBarry Smith 19687497f52SBarry Smith This should only be called in routines that cannot return an error code, such as in C++ constructors. 197fa190f98SMatthew G. Knepley 198db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `PetscCall()`, `CHKMEMQ` 199fa190f98SMatthew G. Knepley M*/ 2009371c9d4SSatish Balay #define SETERRABORT(comm, ierr, ...) \ 2019371c9d4SSatish Balay do { \ 20298921bdaSJacob Faibussowitsch PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \ 20398921bdaSJacob Faibussowitsch MPI_Abort(comm, ierr); \ 20498921bdaSJacob Faibussowitsch } while (0) 2059a00fa46SSatish Balay 20630de9b25SBarry Smith /*MC 2072c71b3e2SJacob Faibussowitsch PetscCheck - Check that a particular condition is true 2082c71b3e2SJacob Faibussowitsch 2092c71b3e2SJacob Faibussowitsch Synopsis: 2102c71b3e2SJacob Faibussowitsch #include <petscerror.h> 2112c71b3e2SJacob Faibussowitsch void PetscCheck(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 2122c71b3e2SJacob Faibussowitsch 2132c71b3e2SJacob Faibussowitsch Collective 2142c71b3e2SJacob Faibussowitsch 2152c71b3e2SJacob Faibussowitsch Input Parameters: 2162c71b3e2SJacob Faibussowitsch + cond - The boolean condition 2172c71b3e2SJacob Faibussowitsch . comm - The communicator on which the check can be collective on 2182c71b3e2SJacob Faibussowitsch . ierr - A nonzero error code, see include/petscerror.h for the complete list 2192c71b3e2SJacob Faibussowitsch - message - Error message in printf format 2202c71b3e2SJacob Faibussowitsch 2212c71b3e2SJacob Faibussowitsch Notes: 2222c71b3e2SJacob Faibussowitsch Enabled in both optimized and debug builds. 2232c71b3e2SJacob Faibussowitsch 22487497f52SBarry Smith Calls `SETERRQ()` if the assertion fails, so can only be called from functions returning a 22587497f52SBarry Smith `PetscErrorCode` (or equivalent type after conversion). 2262c71b3e2SJacob Faibussowitsch 2272c71b3e2SJacob Faibussowitsch Level: beginner 2282c71b3e2SJacob Faibussowitsch 2294be741a6SBarry Smith .seealso: `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheckAbort()` 2309566063dSJacob Faibussowitsch M*/ 2319371c9d4SSatish Balay #define PetscCheck(cond, comm, ierr, ...) \ 2329371c9d4SSatish Balay do { \ 2339371c9d4SSatish Balay if (PetscUnlikely(!(cond))) SETERRQ(comm, ierr, __VA_ARGS__); \ 2349371c9d4SSatish Balay } while (0) 2352c71b3e2SJacob Faibussowitsch 2362c71b3e2SJacob Faibussowitsch /*MC 2374be741a6SBarry Smith PetscCheckAbort - Check that a particular condition is true, otherwise prints error and aborts 2384be741a6SBarry Smith 2394be741a6SBarry Smith Synopsis: 2404be741a6SBarry Smith #include <petscerror.h> 2414be741a6SBarry Smith void PetscCheckAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 2424be741a6SBarry Smith 2434be741a6SBarry Smith Collective 2444be741a6SBarry Smith 2454be741a6SBarry Smith Input Parameters: 2464be741a6SBarry Smith + cond - The boolean condition 2474be741a6SBarry Smith . comm - The communicator on which the check can be collective on 2484be741a6SBarry Smith . ierr - A nonzero error code, see include/petscerror.h for the complete list 2494be741a6SBarry Smith - message - Error message in printf format 2504be741a6SBarry Smith 2514be741a6SBarry Smith Notes: 2524be741a6SBarry Smith Enabled in both optimized and debug builds. 2534be741a6SBarry Smith 25487497f52SBarry Smith Calls `SETERRABORT()` if the assertion fails, can be called from a function that does not return an 25587497f52SBarry Smith error code, such as a C++ constructor. usually `PetscCheck()` should be used. 2564be741a6SBarry Smith 2574be741a6SBarry Smith Level: developer 2584be741a6SBarry Smith 2590e6b6b59SJacob Faibussowitsch .seealso: `PetscAssertAbort()`, `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheck()`, `SETTERRABORT()` 2604be741a6SBarry Smith M*/ 2619371c9d4SSatish Balay #define PetscCheckAbort(cond, comm, ierr, ...) \ 2620e6b6b59SJacob Faibussowitsch do { \ 2630e6b6b59SJacob Faibussowitsch if (PetscUnlikely(!(cond))) SETERRABORT(comm, ierr, __VA_ARGS__); \ 264c7481402SJacob Faibussowitsch } while (0) 2654be741a6SBarry Smith 2664be741a6SBarry Smith /*MC 2679ace16cdSJacob Faibussowitsch PetscAssert - Assert that a particular condition is true 2689ace16cdSJacob Faibussowitsch 2699ace16cdSJacob Faibussowitsch Synopsis: 2709ace16cdSJacob Faibussowitsch #include <petscerror.h> 2719ace16cdSJacob Faibussowitsch void PetscAssert(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 2729ace16cdSJacob Faibussowitsch 2739ace16cdSJacob Faibussowitsch Collective 2749ace16cdSJacob Faibussowitsch 2759ace16cdSJacob Faibussowitsch Input Parameters: 2769ace16cdSJacob Faibussowitsch + cond - The boolean condition 2779ace16cdSJacob Faibussowitsch . comm - The communicator on which the check can be collective on 2789ace16cdSJacob Faibussowitsch . ierr - A nonzero error code, see include/petscerror.h for the complete list 2799ace16cdSJacob Faibussowitsch - message - Error message in printf format 2809ace16cdSJacob Faibussowitsch 2819ace16cdSJacob Faibussowitsch Notes: 282c7481402SJacob Faibussowitsch Equivalent to `PetscCheck()` if debugging is enabled, and `PetscAssume(cond)` otherwise. 2832c71b3e2SJacob Faibussowitsch 28487497f52SBarry Smith See `PetscCheck()` for usage and behaviour. 28587497f52SBarry Smith 28687497f52SBarry Smith This is needed instead of simply using `assert()` because this correctly handles the collective nature of errors under MPI 2879ace16cdSJacob Faibussowitsch 2889ace16cdSJacob Faibussowitsch Level: beginner 2899ace16cdSJacob Faibussowitsch 2900e6b6b59SJacob Faibussowitsch .seealso: `PetscCheck()`, `SETERRQ()`, `PetscError()`, `PetscAssertAbort()` 2919566063dSJacob Faibussowitsch M*/ 292c7481402SJacob Faibussowitsch #if PetscDefined(USE_DEBUG) 293c7481402SJacob Faibussowitsch #define PetscAssert(cond, comm, ierr, ...) PetscCheck(cond, comm, ierr, __VA_ARGS__) 294c7481402SJacob Faibussowitsch #else 295c7481402SJacob Faibussowitsch #define PetscAssert(cond, ...) PetscAssume(cond) 296c7481402SJacob Faibussowitsch #endif 2979ace16cdSJacob Faibussowitsch 2989ace16cdSJacob Faibussowitsch /*MC 2990e6b6b59SJacob Faibussowitsch PetscAssertAbort - Assert that a particular condition is true, otherwise prints error and aborts 3000e6b6b59SJacob Faibussowitsch 3010e6b6b59SJacob Faibussowitsch Synopsis: 3020e6b6b59SJacob Faibussowitsch #include <petscerror.h> 3030e6b6b59SJacob Faibussowitsch void PetscAssertAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 3040e6b6b59SJacob Faibussowitsch 3050e6b6b59SJacob Faibussowitsch Collective 3060e6b6b59SJacob Faibussowitsch 3070e6b6b59SJacob Faibussowitsch Input Parameters: 3080e6b6b59SJacob Faibussowitsch + cond - The boolean condition 3090e6b6b59SJacob Faibussowitsch . comm - The communicator on which the check can be collective on 3100e6b6b59SJacob Faibussowitsch . ierr - A nonzero error code, see include/petscerror.h for the complete list 3110e6b6b59SJacob Faibussowitsch - message - Error message in printf format 3120e6b6b59SJacob Faibussowitsch 3130e6b6b59SJacob Faibussowitsch Notes: 3140e6b6b59SJacob Faibussowitsch Enabled only in debug builds. See `PetscCheckAbort()` for usage. 3150e6b6b59SJacob Faibussowitsch 3160e6b6b59SJacob Faibussowitsch Level: beginner 3170e6b6b59SJacob Faibussowitsch 3180e6b6b59SJacob Faibussowitsch .seealso: `PetscCheckAbort()`, `PetscAssert()`, `PetscCheck()`, `SETERRABORT()`, `PetscError()` 3190e6b6b59SJacob Faibussowitsch M*/ 3200e6b6b59SJacob Faibussowitsch #define PetscAssertAbort(cond, comm, ierr, ...) \ 3210e6b6b59SJacob Faibussowitsch do { \ 3220e6b6b59SJacob Faibussowitsch if (PetscUnlikelyDebug(!(cond))) SETERRABORT(comm, ierr, __VA_ARGS__); \ 3230e6b6b59SJacob Faibussowitsch } while (0) 3240e6b6b59SJacob Faibussowitsch 3250e6b6b59SJacob Faibussowitsch /*MC 32649c86fc7SBarry Smith PetscCall - Calls a PETSc function and then checks the resulting error code, if it is non-zero it calls the error 32749c86fc7SBarry Smith handler and returns from the current function with the error code. 3289566063dSJacob Faibussowitsch 3299566063dSJacob Faibussowitsch Synopsis: 3309566063dSJacob Faibussowitsch #include <petscerror.h> 33149c86fc7SBarry Smith void PetscCall(PetscFunction(args)) 3329566063dSJacob Faibussowitsch 3339566063dSJacob Faibussowitsch Not Collective 3349566063dSJacob Faibussowitsch 3359566063dSJacob Faibussowitsch Input Parameter: 33649c86fc7SBarry Smith . PetscFunction - any PETSc function that returns an error code 3379566063dSJacob Faibussowitsch 3389566063dSJacob Faibussowitsch Notes: 3399566063dSJacob Faibussowitsch Once the error handler is called the calling function is then returned from with the given 34087497f52SBarry Smith error code. Experienced users can set the error handler with `PetscPushErrorHandler()`. 3419566063dSJacob Faibussowitsch 34287497f52SBarry Smith `PetscCall()` cannot be used in functions returning a datatype not convertible to 34387497f52SBarry Smith `PetscErrorCode`. For example, `PetscCall()` may not be used in functions returning void, use 34487497f52SBarry Smith `PetscCallVoid()` in this case. 3459566063dSJacob Faibussowitsch 3469566063dSJacob Faibussowitsch Example Usage: 3479566063dSJacob Faibussowitsch .vb 3489566063dSJacob Faibussowitsch PetscCall(PetscInitiailize(...)); // OK to call even when PETSc is not yet initialized! 3499566063dSJacob Faibussowitsch 3509566063dSJacob Faibussowitsch struct my_struct 3519566063dSJacob Faibussowitsch { 3529566063dSJacob Faibussowitsch void *data; 3539566063dSJacob Faibussowitsch } my_complex_type; 3549566063dSJacob Faibussowitsch 3559566063dSJacob Faibussowitsch struct my_struct bar(void) 3569566063dSJacob Faibussowitsch { 3576aad120cSJose E. Roman PetscCall(foo(15)); // ERROR PetscErrorCode not convertible to struct my_struct! 3589566063dSJacob Faibussowitsch } 3599566063dSJacob Faibussowitsch 3606aad120cSJose E. Roman PetscCall(bar()) // ERROR input not convertible to PetscErrorCode 3619566063dSJacob Faibussowitsch .ve 3629566063dSJacob Faibussowitsch 36387497f52SBarry Smith It is also possible to call this directly on a `PetscErrorCode` variable 36449c86fc7SBarry Smith .vb 36549c86fc7SBarry Smith PetscCall(ierr); // check if ierr is nonzero 36649c86fc7SBarry Smith .ve 36749c86fc7SBarry Smith 368792fecdfSBarry Smith Should not be used to call callback functions provided by users, `PetscCallBack()` should be used in that situation. 369ef1023bdSBarry Smith 3706a8be23eSBarry Smith `PetscUseTypeMethod()` or `PetscTryTypeMethod()` should be used when calling functions pointers contained in a PETSc object's `ops` array 3716a8be23eSBarry Smith 37249c86fc7SBarry Smith Fortran Notes: 37349c86fc7SBarry Smith The Fortran function from which this is used must declare a variable PetscErrorCode ierr and ierr must be 37487497f52SBarry Smith the final argument to the PETSc function being called. 37549c86fc7SBarry Smith 37649c86fc7SBarry Smith In the main program and in Fortran subroutines that do not have ierr as the final return parameter one 37787497f52SBarry Smith should use `PetscCallA()` 37849c86fc7SBarry Smith 37949c86fc7SBarry Smith Example Fortran Usage: 38049c86fc7SBarry Smith .vb 38149c86fc7SBarry Smith PetscErrorCode ierr 38249c86fc7SBarry Smith Vec v 38349c86fc7SBarry Smith 38449c86fc7SBarry Smith ... 38549c86fc7SBarry Smith PetscCall(VecShift(v,1.0,ierr)) 38649c86fc7SBarry Smith PetscCallA(VecShift(v,1.0,ierr)) 38749c86fc7SBarry Smith .ve 38849c86fc7SBarry Smith 3899566063dSJacob Faibussowitsch Level: beginner 3909566063dSJacob Faibussowitsch 39149c86fc7SBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()` 392792fecdfSBarry Smith `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, `CHKERRMPI()`, `PetscCallBack()` 3939566063dSJacob Faibussowitsch M*/ 394ef1023bdSBarry Smith 395ef1023bdSBarry Smith /*MC 396792fecdfSBarry Smith PetscCallBack - Calls a user provided PETSc callback function and then checks the resulting error code, if it is non-zero it calls the error 397ef1023bdSBarry Smith handler and returns from the current function with the error code. 398ef1023bdSBarry Smith 399ef1023bdSBarry Smith Synopsis: 400ef1023bdSBarry Smith #include <petscerror.h> 401792fecdfSBarry Smith void PetscCallBack(const char *functionname,PetscFunction(args)) 402ef1023bdSBarry Smith 403ef1023bdSBarry Smith Not Collective 404ef1023bdSBarry Smith 405ef1023bdSBarry Smith Input Parameters: 406ef1023bdSBarry Smith + functionname - the name of the function being called, this can be a string with spaces that describes the meaning of the callback 40787497f52SBarry Smith - PetscFunction - user provided callback function that returns an error code 408ef1023bdSBarry Smith 409ef1023bdSBarry Smith Notes: 410ef1023bdSBarry Smith Once the error handler is called the calling function is then returned from with the given 41187497f52SBarry Smith error code. Experienced users can set the error handler with `PetscPushErrorHandler()`. 412ef1023bdSBarry Smith 41387497f52SBarry Smith `PetscCallBack()` should only be called in PETSc when a call is being made to a user provided call-back routine. 414ef1023bdSBarry Smith 415ef1023bdSBarry Smith Example Usage: 416ef1023bdSBarry Smith .vb 417792fecdfSBarry Smith PetscCallBack("XXX callback to do something",a->callback(...)); 418ef1023bdSBarry Smith .ve 419ef1023bdSBarry Smith 420ef1023bdSBarry Smith Level: developer 421ef1023bdSBarry Smith 42287497f52SBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscCall()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()` 423ef1023bdSBarry Smith `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, `CHKERRMPI()`, `PetscCall()` 424ef1023bdSBarry Smith M*/ 425ef1023bdSBarry Smith 4269566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 4279566063dSJacob Faibussowitsch void PetscCall(PetscErrorCode); 428792fecdfSBarry Smith void PetscCallBack(const char *, PetscErrorCode); 4299566063dSJacob Faibussowitsch void PetscCallVoid(PetscErrorCode); 4309566063dSJacob Faibussowitsch #else 4319371c9d4SSatish Balay #define PetscCall(...) \ 4329371c9d4SSatish Balay do { \ 433e33ced7fSLisandro Dalcin PetscErrorCode ierr_q_; \ 43437154d10SBarry Smith PetscStackUpdateLine; \ 435e33ced7fSLisandro Dalcin ierr_q_ = __VA_ARGS__; \ 4369566063dSJacob Faibussowitsch if (PetscUnlikely(ierr_q_)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_q_, PETSC_ERROR_REPEAT, " "); \ 4379566063dSJacob Faibussowitsch } while (0) 4389371c9d4SSatish Balay #define PetscCallBack(function, ...) \ 4399371c9d4SSatish Balay do { \ 440e33ced7fSLisandro Dalcin PetscErrorCode ierr_q_; \ 441e33ced7fSLisandro Dalcin PetscStackUpdateLine; \ 442e33ced7fSLisandro Dalcin PetscStackPushExternal(function); \ 443e33ced7fSLisandro Dalcin ierr_q_ = __VA_ARGS__; \ 444ef1023bdSBarry Smith PetscStackPop; \ 445ef1023bdSBarry Smith if (PetscUnlikely(ierr_q_)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_q_, PETSC_ERROR_REPEAT, " "); \ 446ef1023bdSBarry Smith } while (0) 4479371c9d4SSatish Balay #define PetscCallVoid(...) \ 4489371c9d4SSatish Balay do { \ 449e33ced7fSLisandro Dalcin PetscErrorCode ierr_void_; \ 450e33ced7fSLisandro Dalcin PetscStackUpdateLine; \ 451e33ced7fSLisandro Dalcin ierr_void_ = __VA_ARGS__; \ 4529371c9d4SSatish Balay if (PetscUnlikely(ierr_void_)) { \ 4539371c9d4SSatish Balay (void)PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_void_, PETSC_ERROR_REPEAT, " "); \ 4549371c9d4SSatish Balay return; \ 4559371c9d4SSatish Balay } \ 4569566063dSJacob Faibussowitsch } while (0) 4579566063dSJacob Faibussowitsch #endif 4589566063dSJacob Faibussowitsch 4599566063dSJacob Faibussowitsch /*MC 4609566063dSJacob Faibussowitsch CHKERRQ - Checks error code returned from PETSc function 46130de9b25SBarry Smith 46230de9b25SBarry Smith Synopsis: 463aaa7dc30SBarry Smith #include <petscsys.h> 4649566063dSJacob Faibussowitsch void CHKERRQ(PetscErrorCode ierr) 4659566063dSJacob Faibussowitsch 4669566063dSJacob Faibussowitsch Not Collective 4679566063dSJacob Faibussowitsch 4689566063dSJacob Faibussowitsch Input Parameters: 4699566063dSJacob Faibussowitsch . ierr - nonzero error code 4709566063dSJacob Faibussowitsch 4719566063dSJacob Faibussowitsch Notes: 47287497f52SBarry Smith Deprecated in favor of `PetscCall()`. This routine behaves identically to it. 4739566063dSJacob Faibussowitsch 4749566063dSJacob Faibussowitsch Level: deprecated 4759566063dSJacob Faibussowitsch 476db781477SPatrick Sanan .seealso: `PetscCall()` 4779566063dSJacob Faibussowitsch M*/ 4789566063dSJacob Faibussowitsch #define CHKERRQ(...) PetscCall(__VA_ARGS__) 4799566063dSJacob Faibussowitsch #define CHKERRV(...) PetscCallVoid(__VA_ARGS__) 4809566063dSJacob Faibussowitsch 481db9cea48SBarry Smith PETSC_EXTERN void PetscMPIErrorString(PetscMPIInt, char *); 482db9cea48SBarry Smith 4839566063dSJacob Faibussowitsch /*MC 4849566063dSJacob Faibussowitsch PetscCallMPI - Checks error code returned from MPI calls, if non-zero it calls the error 4859566063dSJacob Faibussowitsch handler and then returns 4869566063dSJacob Faibussowitsch 4879566063dSJacob Faibussowitsch Synopsis: 4889566063dSJacob Faibussowitsch #include <petscerror.h> 48949c86fc7SBarry Smith void PetscCallMPI(MPI_Function(args)) 49030de9b25SBarry Smith 491eca87e8dSBarry Smith Not Collective 49230de9b25SBarry Smith 49330de9b25SBarry Smith Input Parameters: 49449c86fc7SBarry Smith . MPI_Function - an MPI function that returns an MPI error code 49530de9b25SBarry Smith 4969566063dSJacob Faibussowitsch Notes: 49787497f52SBarry Smith Always returns the error code `PETSC_ERR_MPI`; the MPI error code and string are embedded in 4989566063dSJacob Faibussowitsch the string error message. Do not use this to call any other routines (for example PETSc 4999566063dSJacob Faibussowitsch routines), it should only be used for direct MPI calls. Due to limitations of the 5009566063dSJacob Faibussowitsch preprocessor this can unfortunately not easily be enforced, so the user should take care to 5019566063dSJacob Faibussowitsch check this themselves. 5029566063dSJacob Faibussowitsch 5039566063dSJacob Faibussowitsch Example Usage: 5049566063dSJacob Faibussowitsch .vb 5059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(...)); // OK, calling MPI function 5069566063dSJacob Faibussowitsch 5079566063dSJacob Faibussowitsch PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead! 5089566063dSJacob Faibussowitsch .ve 5099566063dSJacob Faibussowitsch 51049c86fc7SBarry Smith Fortran Notes: 51187497f52SBarry Smith The Fortran function from which this is used must declare a variable `PetscErrorCode` ierr and ierr must be 51249c86fc7SBarry Smith the final argument to the MPI function being called. 51349c86fc7SBarry Smith 51449c86fc7SBarry Smith In the main program and in Fortran subroutines that do not have ierr as the final return parameter one 51587497f52SBarry Smith should use `PetscCallMPIA()` 51649c86fc7SBarry Smith 51749c86fc7SBarry Smith Fortran Usage: 51849c86fc7SBarry Smith .vb 51949c86fc7SBarry Smith PetscErrorCode ierr or integer ierr 52049c86fc7SBarry Smith ... 52149c86fc7SBarry Smith PetscCallMPI(MPI_Comm_size(...,ierr)) 52249c86fc7SBarry Smith PetscCallMPIA(MPI_Comm_size(...,ierr)) ! Will abort after calling error handler 52349c86fc7SBarry Smith 52449c86fc7SBarry Smith PetscCallMPI(MPI_Comm_size(...,eflag)) ! ERROR, final argument must be ierr 52549c86fc7SBarry Smith .ve 52649c86fc7SBarry Smith 52730de9b25SBarry Smith Level: beginner 52830de9b25SBarry Smith 529db781477SPatrick Sanan .seealso: `SETERRMPI()`, `PetscCall()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`, 530db781477SPatrick Sanan `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ` 53130de9b25SBarry Smith M*/ 5323fcd9f07SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 53363a3b9bcSJacob Faibussowitsch void PetscCallMPI(PetscMPIInt); 534064a246eSJacob Faibussowitsch #else 5359371c9d4SSatish Balay #define PetscCallMPI(...) \ 5369371c9d4SSatish Balay do { \ 537e33ced7fSLisandro Dalcin PetscMPIInt _7_errorcode; \ 538e33ced7fSLisandro Dalcin char _7_errorstring[2 * MPI_MAX_ERROR_STRING]; \ 539ef1023bdSBarry Smith PetscStackUpdateLine; \ 540792fecdfSBarry Smith PetscStackPushExternal("MPI function"); \ 541d71ae5a4SJacob Faibussowitsch { \ 542d71ae5a4SJacob Faibussowitsch _7_errorcode = __VA_ARGS__; \ 543d71ae5a4SJacob Faibussowitsch } \ 544ef1023bdSBarry Smith PetscStackPop; \ 5459566063dSJacob Faibussowitsch if (PetscUnlikely(_7_errorcode)) { \ 546db9cea48SBarry Smith PetscMPIErrorString(_7_errorcode, (char *)_7_errorstring); \ 547db9cea48SBarry Smith SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MPI, "MPI error %d %s", (int)_7_errorcode, _7_errorstring); \ 5483fcd9f07SJacob Faibussowitsch } \ 5493fcd9f07SJacob Faibussowitsch } while (0) 5509566063dSJacob Faibussowitsch #endif 551064a246eSJacob Faibussowitsch 5527037b0edSPatrick Sanan /*MC 5539566063dSJacob Faibussowitsch CHKERRMPI - Checks error code returned from MPI calls, if non-zero it calls the error 5549566063dSJacob Faibussowitsch handler and then returns 5559566063dSJacob Faibussowitsch 5569566063dSJacob Faibussowitsch Synopsis: 5579566063dSJacob Faibussowitsch #include <petscerror.h> 5589566063dSJacob Faibussowitsch void CHKERRMPI(PetscErrorCode ierr) 5599566063dSJacob Faibussowitsch 5609566063dSJacob Faibussowitsch Not Collective 5619566063dSJacob Faibussowitsch 5629566063dSJacob Faibussowitsch Input Parameter: 5639566063dSJacob Faibussowitsch . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 5649566063dSJacob Faibussowitsch 5659566063dSJacob Faibussowitsch Notes: 56687497f52SBarry Smith Deprecated in favor of `PetscCallMPI()`. This routine behaves identically to it. 5679566063dSJacob Faibussowitsch 5689566063dSJacob Faibussowitsch Level: deprecated 5699566063dSJacob Faibussowitsch 570db781477SPatrick Sanan .seealso: `PetscCallMPI()` 5719566063dSJacob Faibussowitsch M*/ 5729566063dSJacob Faibussowitsch #define CHKERRMPI(...) PetscCallMPI(__VA_ARGS__) 5739566063dSJacob Faibussowitsch 5749566063dSJacob Faibussowitsch /*MC 5759566063dSJacob Faibussowitsch PetscCallAbort - Checks error code returned from PETSc function, if non-zero it aborts immediately 5769566063dSJacob Faibussowitsch 5779566063dSJacob Faibussowitsch Synopsis: 5789566063dSJacob Faibussowitsch #include <petscerror.h> 5799566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm comm, PetscErrorCode ierr) 5809566063dSJacob Faibussowitsch 581c3339decSBarry Smith Collective 5829566063dSJacob Faibussowitsch 5839566063dSJacob Faibussowitsch Input Parameters: 5849566063dSJacob Faibussowitsch + comm - the MPI communicator on which to abort 5859566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 5869566063dSJacob Faibussowitsch 5879566063dSJacob Faibussowitsch Notes: 58887497f52SBarry Smith This macro has identical type and usage semantics to `PetscCall()` with the important caveat 5899566063dSJacob Faibussowitsch that this macro does not return. Instead, if ierr is nonzero it calls the PETSc error handler 59087497f52SBarry Smith and then immediately calls `MPI_Abort()`. It can therefore be used anywhere. 5919566063dSJacob Faibussowitsch 59287497f52SBarry Smith As per `MPI_Abort()` semantics the communicator passed must be valid, although there is currently 59387497f52SBarry Smith no attempt made at handling any potential errors from `MPI_Abort()`. Note that while 59487497f52SBarry Smith `MPI_Abort()` is required to terminate only those processes which reside on comm, it is often 59587497f52SBarry Smith the case that `MPI_Abort()` terminates *all* processes. 5969566063dSJacob Faibussowitsch 5979566063dSJacob Faibussowitsch Example Usage: 5989566063dSJacob Faibussowitsch .vb 5999566063dSJacob Faibussowitsch PetscErrorCode boom(void) { return PETSC_ERR_MEM; } 6009566063dSJacob Faibussowitsch 6019566063dSJacob Faibussowitsch void foo(void) 6029566063dSJacob Faibussowitsch { 6039566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type 6049566063dSJacob Faibussowitsch } 6059566063dSJacob Faibussowitsch 6069566063dSJacob Faibussowitsch double bar(void) 6079566063dSJacob Faibussowitsch { 6089566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type 6099566063dSJacob Faibussowitsch } 6109566063dSJacob Faibussowitsch 6119566063dSJacob Faibussowitsch PetscCallAbort(MPI_COMM_NULL,boom()); // ERROR, communicator should be valid 6129566063dSJacob Faibussowitsch 6139566063dSJacob Faibussowitsch struct baz 6149566063dSJacob Faibussowitsch { 6159566063dSJacob Faibussowitsch baz() 6169566063dSJacob Faibussowitsch { 6179566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,boom()); // OK 6189566063dSJacob Faibussowitsch } 6199566063dSJacob Faibussowitsch 6209566063dSJacob Faibussowitsch ~baz() 6219566063dSJacob Faibussowitsch { 6229566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,boom()); // OK (in fact the only way to handle PETSc errors) 6239566063dSJacob Faibussowitsch } 6249566063dSJacob Faibussowitsch }; 6259566063dSJacob Faibussowitsch .ve 6269566063dSJacob Faibussowitsch 6279566063dSJacob Faibussowitsch Level: intermediate 6289566063dSJacob Faibussowitsch 629db781477SPatrick Sanan .seealso: `SETERRABORT()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, 6305687f061SJacob Faibussowitsch `SETERRQ()`, `CHKMEMQ`, `PetscCallMPI()`, `PetscCallCXXAbort()` 6319566063dSJacob Faibussowitsch M*/ 6329566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 6339566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm, PetscErrorCode); 6349566063dSJacob Faibussowitsch void PetscCallContinue(PetscErrorCode); 6359566063dSJacob Faibussowitsch #else 6369371c9d4SSatish Balay #define PetscCallAbort(comm, ...) \ 6379371c9d4SSatish Balay do { \ 6389566063dSJacob Faibussowitsch PetscErrorCode ierr_abort_ = __VA_ARGS__; \ 6399566063dSJacob Faibussowitsch if (PetscUnlikely(ierr_abort_)) { \ 6409566063dSJacob Faibussowitsch PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_abort_, PETSC_ERROR_REPEAT, " "); \ 6419566063dSJacob Faibussowitsch MPI_Abort(comm, ierr_abort_); \ 6429566063dSJacob Faibussowitsch } \ 6439566063dSJacob Faibussowitsch } while (0) 6449371c9d4SSatish Balay #define PetscCallContinue(...) \ 6459371c9d4SSatish Balay do { \ 6469566063dSJacob Faibussowitsch PetscErrorCode ierr_continue_ = __VA_ARGS__; \ 6479566063dSJacob Faibussowitsch if (PetscUnlikely(ierr_continue_)) PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_continue_, PETSC_ERROR_REPEAT, " "); \ 6489566063dSJacob Faibussowitsch } while (0) 6499566063dSJacob Faibussowitsch #endif 6509566063dSJacob Faibussowitsch 6519566063dSJacob Faibussowitsch /*MC 6529566063dSJacob Faibussowitsch CHKERRABORT - Checks error code returned from PETSc function. If non-zero it aborts immediately. 6539566063dSJacob Faibussowitsch 6549566063dSJacob Faibussowitsch Synopsis: 6559566063dSJacob Faibussowitsch #include <petscerror.h> 6569566063dSJacob Faibussowitsch void CHKERRABORT(MPI_Comm comm, PetscErrorCode ierr) 6579566063dSJacob Faibussowitsch 6589566063dSJacob Faibussowitsch Not Collective 6599566063dSJacob Faibussowitsch 6609566063dSJacob Faibussowitsch Input Parameters: 6619566063dSJacob Faibussowitsch + comm - the MPI communicator 6629566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 6639566063dSJacob Faibussowitsch 6649566063dSJacob Faibussowitsch Notes: 66587497f52SBarry Smith Deprecated in favor of `PetscCallAbort()`. This routine behaves identically to it. 6669566063dSJacob Faibussowitsch 6679566063dSJacob Faibussowitsch Level: deprecated 6689566063dSJacob Faibussowitsch 669db781477SPatrick Sanan .seealso: `PetscCallAbort()` 6709566063dSJacob Faibussowitsch M*/ 6719566063dSJacob Faibussowitsch #define CHKERRABORT(comm, ...) PetscCallAbort(comm, __VA_ARGS__) 6729566063dSJacob Faibussowitsch #define CHKERRCONTINUE(...) PetscCallContinue(__VA_ARGS__) 6739566063dSJacob Faibussowitsch 6749566063dSJacob Faibussowitsch /*MC 67587497f52SBarry Smith CHKERRA - Fortran-only replacement for use of `CHKERRQ()` in the main program, which aborts immediately 676f388eb8bSPatrick Sanan 677f388eb8bSPatrick Sanan Synopsis: 678f388eb8bSPatrick Sanan #include <petscsys.h> 679f388eb8bSPatrick Sanan PetscErrorCode CHKERRA(PetscErrorCode ierr) 680f388eb8bSPatrick Sanan 681f388eb8bSPatrick Sanan Not Collective 682f388eb8bSPatrick Sanan 683f388eb8bSPatrick Sanan Input Parameters: 684f388eb8bSPatrick Sanan . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 685f388eb8bSPatrick Sanan 68687497f52SBarry Smith Level: deprecated 687f388eb8bSPatrick Sanan 68887497f52SBarry Smith Note: 68987497f52SBarry Smith This macro is rarely needed, normal usage is `PetscCallA()` in the main Fortran program. 690f388eb8bSPatrick Sanan 69187497f52SBarry Smith .seealso: `PetscCall()`, `PetscCallA()`, `PetscCallAbort()`, `CHKERRQ()`, `SETERRA()`, `SETERRQ()`, `SETERRABORT()` 692f388eb8bSPatrick Sanan M*/ 693f388eb8bSPatrick Sanan 6941e4f893aSSatish Balay PETSC_EXTERN PetscBool petscwaitonerrorflg; 6951e4f893aSSatish Balay PETSC_EXTERN PetscBool petscindebugger; 6967c66cc67SJunchao Zhang 6977c66cc67SJunchao Zhang /*MC 6987c66cc67SJunchao Zhang PETSCABORT - Call MPI_Abort with an informative error code 6997c66cc67SJunchao Zhang 7007c66cc67SJunchao Zhang Synopsis: 7017c66cc67SJunchao Zhang #include <petscsys.h> 7027c66cc67SJunchao Zhang PETSCABORT(MPI_Comm comm, PetscErrorCode ierr) 7037c66cc67SJunchao Zhang 7047c66cc67SJunchao Zhang Collective 7057c66cc67SJunchao Zhang 7067c66cc67SJunchao Zhang Input Parameters: 7077c66cc67SJunchao Zhang + comm - A communicator, so that the error can be collective 7087c66cc67SJunchao Zhang - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 7097c66cc67SJunchao Zhang 710bf4d2887SBarry Smith Level: advanced 7117c66cc67SJunchao Zhang 712bf4d2887SBarry Smith Notes: 713bf4d2887SBarry Smith If the option -start_in_debugger was used then this calls abort() to stop the program in the debugger. 714bf4d2887SBarry Smith 71587497f52SBarry Smith if `PetscCIEnabledPortableErrorOutput` is set it strives to exit cleanly without call `MPI_Abort()` 716660278c0SBarry Smith 7177c66cc67SJunchao Zhang M*/ 7189371c9d4SSatish Balay #define PETSCABORT(comm, ...) \ 7199371c9d4SSatish Balay do { \ 720baae8e41SSatish Balay if (petscwaitonerrorflg) PetscSleep(1000); \ 721bf4d2887SBarry Smith if (petscindebugger) abort(); \ 7223fcd9f07SJacob Faibussowitsch else { \ 7233fcd9f07SJacob Faibussowitsch PetscErrorCode ierr_petsc_abort_ = __VA_ARGS__; \ 724660278c0SBarry Smith PetscMPIInt size; \ 725660278c0SBarry Smith MPI_Comm_size(comm, &size); \ 726660278c0SBarry Smith if (PetscCIEnabledPortableErrorOutput && size == PetscGlobalSize && ierr_petsc_abort_ != PETSC_ERR_SIG) { \ 7279371c9d4SSatish Balay MPI_Finalize(); \ 7289371c9d4SSatish Balay exit(0); \ 729660278c0SBarry Smith } else if (PetscCIEnabledPortableErrorOutput && PetscGlobalSize == 1) { \ 730660278c0SBarry Smith exit(0); \ 731660278c0SBarry Smith } else { \ 732660278c0SBarry Smith MPI_Abort(comm, (PetscMPIInt)ierr_petsc_abort_); \ 733660278c0SBarry Smith } \ 7343fcd9f07SJacob Faibussowitsch } \ 7357c66cc67SJunchao Zhang } while (0) 736986eef2eSBarry Smith 7379566063dSJacob Faibussowitsch #ifdef PETSC_CLANGUAGE_CXX 738986eef2eSBarry Smith /*MC 7399566063dSJacob Faibussowitsch PetscCallThrow - Checks error code, if non-zero it calls the C++ error handler which throws 7409566063dSJacob Faibussowitsch an exception 741986eef2eSBarry Smith 742986eef2eSBarry Smith Synopsis: 7439566063dSJacob Faibussowitsch #include <petscerror.h> 7449566063dSJacob Faibussowitsch void PetscCallThrow(PetscErrorCode ierr) 745986eef2eSBarry Smith 746986eef2eSBarry Smith Not Collective 747986eef2eSBarry Smith 7489566063dSJacob Faibussowitsch Input Parameter: 749986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 750986eef2eSBarry Smith 751986eef2eSBarry Smith Notes: 7529566063dSJacob Faibussowitsch Requires PETSc to be configured with clanguage = c++. Throws a std::runtime_error() on error. 753986eef2eSBarry Smith 75487497f52SBarry Smith Once the error handler throws the exception you can use `PetscCallVoid()` which returns without 75587497f52SBarry Smith an error code (bad idea since the error is ignored) or `PetscCallAbort()` to have `MPI_Abort()` 7569566063dSJacob Faibussowitsch called immediately. 7579566063dSJacob Faibussowitsch 7589566063dSJacob Faibussowitsch Level: beginner 7599566063dSJacob Faibussowitsch 760db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscCall()`, `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, 761db781477SPatrick Sanan `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ` 762986eef2eSBarry Smith M*/ 7639371c9d4SSatish Balay #define PetscCallThrow(...) \ 7649371c9d4SSatish Balay do { \ 7659566063dSJacob Faibussowitsch PetscErrorCode ierr_cxx_ = __VA_ARGS__; \ 7669566063dSJacob Faibussowitsch if (PetscUnlikely(ierr_cxx_)) PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_cxx_, PETSC_ERROR_IN_CXX, PETSC_NULLPTR); \ 767ffc4695bSBarry Smith } while (0) 76885614651SBarry Smith 769cc26af49SMatthew Knepley /*MC 770cc26af49SMatthew Knepley CHKERRXX - Checks error code, if non-zero it calls the C++ error handler which throws an exception 771cc26af49SMatthew Knepley 772cc26af49SMatthew Knepley Synopsis: 7739566063dSJacob Faibussowitsch #include <petscerror.h> 7743af045c5SBarry Smith void CHKERRXX(PetscErrorCode ierr) 775cc26af49SMatthew Knepley 776eca87e8dSBarry Smith Not Collective 777cc26af49SMatthew Knepley 7789566063dSJacob Faibussowitsch Input Parameter: 7793af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 780cc26af49SMatthew Knepley 781cc26af49SMatthew Knepley Notes: 78287497f52SBarry Smith Deprecated in favor of `PetscCallThrow()`. This routine behaves identically to it. 783cc26af49SMatthew Knepley 7849566063dSJacob Faibussowitsch Level: deprecated 785cc26af49SMatthew Knepley 786db781477SPatrick Sanan .seealso: `PetscCallThrow()` 787cc26af49SMatthew Knepley M*/ 7889566063dSJacob Faibussowitsch #define CHKERRXX(...) PetscCallThrow(__VA_ARGS__) 789fd705b32SMatthew Knepley #endif 790fd705b32SMatthew Knepley 7913f520e80SJunchao Zhang /*MC 7929566063dSJacob Faibussowitsch PetscCallCXX - Checks C++ function calls and if they throw an exception, catch it and then 7939566063dSJacob Faibussowitsch return a PETSc error code 7943f520e80SJunchao Zhang 7953f520e80SJunchao Zhang Synopsis: 7969566063dSJacob Faibussowitsch #include <petscerror.h> 7975687f061SJacob Faibussowitsch void PetscCallCXX(...) noexcept; 7983f520e80SJunchao Zhang 7993f520e80SJunchao Zhang Not Collective 8003f520e80SJunchao Zhang 8019566063dSJacob Faibussowitsch Input Parameter: 8025687f061SJacob Faibussowitsch . __VA_ARGS__ - An arbitrary expression 8035687f061SJacob Faibussowitsch 8045687f061SJacob Faibussowitsch Level: beginner 8053f520e80SJunchao Zhang 8063f520e80SJunchao Zhang Notes: 8075687f061SJacob Faibussowitsch `PetscCallCXX(...)` is a macro replacement for 8089566063dSJacob Faibussowitsch .vb 8099566063dSJacob Faibussowitsch try { 8105687f061SJacob Faibussowitsch __VA_ARGS__; 8119566063dSJacob Faibussowitsch } catch (const std::exception& e) { 8129566063dSJacob Faibussowitsch return ConvertToPetscErrorCode(e); 8139566063dSJacob Faibussowitsch } 8149566063dSJacob Faibussowitsch .ve 8159566063dSJacob Faibussowitsch Due to the fact that it catches any (reasonable) exception, it is essentially noexcept. 8163f520e80SJunchao Zhang 8175687f061SJacob Faibussowitsch If you cannot return a `PetscErrorCode` use `PetscCallCXXAbort()` instead. 8185687f061SJacob Faibussowitsch 8199566063dSJacob Faibussowitsch Example Usage: 8209566063dSJacob Faibussowitsch .vb 8219566063dSJacob Faibussowitsch void foo(void) { throw std::runtime_error("error"); } 8223f520e80SJunchao Zhang 8239566063dSJacob Faibussowitsch void bar() 8249566063dSJacob Faibussowitsch { 825e8952933SJacob Faibussowitsch PetscCallCXX(foo()); // ERROR bar() does not return PetscErrorCode 8269566063dSJacob Faibussowitsch } 8279566063dSJacob Faibussowitsch 8289566063dSJacob Faibussowitsch PetscErrorCode baz() 8299566063dSJacob Faibussowitsch { 8309566063dSJacob Faibussowitsch PetscCallCXX(foo()); // OK 8319566063dSJacob Faibussowitsch 8329566063dSJacob Faibussowitsch PetscCallCXX( 8339566063dSJacob Faibussowitsch bar(); 834d5b43468SJose E. Roman foo(); // OK multiple statements allowed 8359566063dSJacob Faibussowitsch ); 8369566063dSJacob Faibussowitsch } 8379566063dSJacob Faibussowitsch 8389566063dSJacob Faibussowitsch struct bop 8399566063dSJacob Faibussowitsch { 8409566063dSJacob Faibussowitsch bop() 8419566063dSJacob Faibussowitsch { 842e8952933SJacob Faibussowitsch PetscCallCXX(foo()); // ERROR returns PetscErrorCode, cannot be used in constructors 8439566063dSJacob Faibussowitsch } 8449566063dSJacob Faibussowitsch }; 8459566063dSJacob Faibussowitsch 846e8952933SJacob Faibussowitsch // ERROR contains do-while, cannot be used as function-try block 8479566063dSJacob Faibussowitsch PetscErrorCode qux() PetscCallCXX( 8489566063dSJacob Faibussowitsch bar(); 8499566063dSJacob Faibussowitsch baz(); 8509566063dSJacob Faibussowitsch foo(); 8519566063dSJacob Faibussowitsch return 0; 8529566063dSJacob Faibussowitsch ) 8539566063dSJacob Faibussowitsch .ve 8549566063dSJacob Faibussowitsch 8555687f061SJacob Faibussowitsch .seealso: `PetscCallCXXAbort()`, `PetscCallThrow()`, `SETERRQ()`, `PetscCall()`, 8565687f061SJacob Faibussowitsch `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 8575687f061SJacob Faibussowitsch `PetscError()`, `CHKMEMQ` 8583f520e80SJunchao Zhang M*/ 8595687f061SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 8605687f061SJacob Faibussowitsch void PetscCallCXX(PetscErrorCode); 8615687f061SJacob Faibussowitsch #else 8629371c9d4SSatish Balay #define PetscCallCXX(...) \ 8639371c9d4SSatish Balay do { \ 86437154d10SBarry Smith PetscStackUpdateLine; \ 8653fcd9f07SJacob Faibussowitsch try { \ 8663fcd9f07SJacob Faibussowitsch __VA_ARGS__; \ 867d71ae5a4SJacob Faibussowitsch } catch (const std::exception &e) { \ 868d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "%s", e.what()); \ 869d71ae5a4SJacob Faibussowitsch } \ 8703fcd9f07SJacob Faibussowitsch } while (0) 8715687f061SJacob Faibussowitsch #endif 8725687f061SJacob Faibussowitsch 8735687f061SJacob Faibussowitsch /*MC 8745687f061SJacob Faibussowitsch PetscCallCXXAbort - Like `PetscCallCXX()` but calls `MPI_Abort()` instead of returning an 8755687f061SJacob Faibussowitsch error-code 8765687f061SJacob Faibussowitsch 8775687f061SJacob Faibussowitsch Synopsis: 8785687f061SJacob Faibussowitsch #include <petscerror.h> 8795687f061SJacob Faibussowitsch void PetscCallCXXAbort(MPI_Comm comm, ...) noexcept; 8805687f061SJacob Faibussowitsch 8815687f061SJacob Faibussowitsch Collective on `comm` 8825687f061SJacob Faibussowitsch 8835687f061SJacob Faibussowitsch Input Parameters: 8845687f061SJacob Faibussowitsch + comm - The MPI communicator to abort on 8855687f061SJacob Faibussowitsch - __VA_ARGS__ - An arbitrary expression 8865687f061SJacob Faibussowitsch 8875687f061SJacob Faibussowitsch Level: beginner 8885687f061SJacob Faibussowitsch 8895687f061SJacob Faibussowitsch Notes: 8905687f061SJacob Faibussowitsch This macro may be used to check C++ expressions for exceptions in cases where you cannot 8915687f061SJacob Faibussowitsch return an error code. This includes constructors, destructors, copy/move assignment functions 8925687f061SJacob Faibussowitsch or constructors among others. 8935687f061SJacob Faibussowitsch 8945687f061SJacob Faibussowitsch If an exception is caught, the macro calls `SETERRABORT()` on `comm`. The exception must 8955687f061SJacob Faibussowitsch derive from `std::exception` in order to be caught. 8965687f061SJacob Faibussowitsch 8975687f061SJacob Faibussowitsch If the routine _can_ return an error-code it is highly advised to use `PetscCallCXX()` 8985687f061SJacob Faibussowitsch instead. 8995687f061SJacob Faibussowitsch 9005687f061SJacob Faibussowitsch See `PetscCallCXX()` for additional discussion. 9015687f061SJacob Faibussowitsch 9025687f061SJacob Faibussowitsch Fortran Note: 9035687f061SJacob Faibussowitsch Not available from Fortran. 9045687f061SJacob Faibussowitsch 9055687f061SJacob Faibussowitsch Example Usage: 9065687f061SJacob Faibussowitsch .vb 9075687f061SJacob Faibussowitsch class Foo 9085687f061SJacob Faibussowitsch { 9095687f061SJacob Faibussowitsch std::vector<int> data_; 9105687f061SJacob Faibussowitsch 9115687f061SJacob Faibussowitsch public: 9125687f061SJacob Faibussowitsch // normally std::vector::reserve() may raise an exception, but since we handle it with 9135687f061SJacob Faibussowitsch // PetscCallCXXAbort() we may mark this routine as noexcept! 9145687f061SJacob Faibussowitsch Foo() noexcept 9155687f061SJacob Faibussowitsch { 9165687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, data_.reserve(10)); 9175687f061SJacob Faibussowitsch } 9185687f061SJacob Faibussowitsch }; 9195687f061SJacob Faibussowitsch 9205687f061SJacob Faibussowitsch std::vector<int> bar() 9215687f061SJacob Faibussowitsch { 9225687f061SJacob Faibussowitsch std::vector<int> v; 9235687f061SJacob Faibussowitsch 9245687f061SJacob Faibussowitsch PetscFunctionBegin; 9255687f061SJacob Faibussowitsch // OK! 9265687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1)); 9275687f061SJacob Faibussowitsch PetscFunctionReturn(v); 9285687f061SJacob Faibussowitsch } 9295687f061SJacob Faibussowitsch 9305687f061SJacob Faibussowitsch PetscErrorCode baz() 9315687f061SJacob Faibussowitsch { 9325687f061SJacob Faibussowitsch std::vector<int> v; 9335687f061SJacob Faibussowitsch 9345687f061SJacob Faibussowitsch PetscFunctionBegin; 9355687f061SJacob Faibussowitsch // WRONG! baz() returns a PetscErrorCode, prefer PetscCallCXX() instead 9365687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1)); 9375687f061SJacob Faibussowitsch PetscFunctionReturn(0); 9385687f061SJacob Faibussowitsch } 9395687f061SJacob Faibussowitsch .ve 9405687f061SJacob Faibussowitsch 9415687f061SJacob Faibussowitsch .seealso: `PetscCallCXX()`, `SETERRABORT()`, `PetscCallAbort()` 9425687f061SJacob Faibussowitsch M*/ 9435687f061SJacob Faibussowitsch #define PetscCallCXXAbort(comm, ...) \ 9445687f061SJacob Faibussowitsch do { \ 9455687f061SJacob Faibussowitsch PetscStackUpdateLine; \ 9465687f061SJacob Faibussowitsch try { \ 9475687f061SJacob Faibussowitsch __VA_ARGS__; \ 9485687f061SJacob Faibussowitsch } catch (const std::exception &e) { \ 9495687f061SJacob Faibussowitsch SETERRABORT(comm, PETSC_ERR_LIB, "%s", e.what()); \ 9505687f061SJacob Faibussowitsch } \ 9515687f061SJacob Faibussowitsch } while (0) 9525687f061SJacob Faibussowitsch 9535687f061SJacob Faibussowitsch #define PetscCallCXXAbort(comm, ...) \ 9545687f061SJacob Faibussowitsch do { \ 9555687f061SJacob Faibussowitsch PetscStackUpdateLine; \ 9565687f061SJacob Faibussowitsch try { \ 9575687f061SJacob Faibussowitsch __VA_ARGS__; \ 9585687f061SJacob Faibussowitsch } catch (const std::exception &e) { \ 9595687f061SJacob Faibussowitsch SETERRABORT(comm, PETSC_ERR_LIB, "%s", e.what()); \ 9605687f061SJacob Faibussowitsch } \ 9615687f061SJacob Faibussowitsch } while (0) 9623f520e80SJunchao Zhang 96330de9b25SBarry Smith /*MC 9649566063dSJacob Faibussowitsch CHKERRCXX - Checks C++ function calls and if they throw an exception, catch it and then 9659566063dSJacob Faibussowitsch return a PETSc error code 9669566063dSJacob Faibussowitsch 9679566063dSJacob Faibussowitsch Synopsis: 9689566063dSJacob Faibussowitsch #include <petscerror.h> 9699566063dSJacob Faibussowitsch void CHKERRCXX(func) noexcept; 9709566063dSJacob Faibussowitsch 9719566063dSJacob Faibussowitsch Not Collective 9729566063dSJacob Faibussowitsch 9739566063dSJacob Faibussowitsch Input Parameter: 9749566063dSJacob Faibussowitsch . func - C++ function calls 9759566063dSJacob Faibussowitsch 9769566063dSJacob Faibussowitsch Notes: 97787497f52SBarry Smith Deprecated in favor of `PetscCallCXX()`. This routine behaves identically to it. 9789566063dSJacob Faibussowitsch 9799566063dSJacob Faibussowitsch Level: deprecated 9809566063dSJacob Faibussowitsch 981db781477SPatrick Sanan .seealso: `PetscCallCXX()` 9829566063dSJacob Faibussowitsch M*/ 9839566063dSJacob Faibussowitsch #define CHKERRCXX(...) PetscCallCXX(__VA_ARGS__) 9849566063dSJacob Faibussowitsch 9859566063dSJacob Faibussowitsch /*MC 98630de9b25SBarry Smith CHKMEMQ - Checks the memory for corruption, calls error handler if any is detected 98730de9b25SBarry Smith 98830de9b25SBarry Smith Synopsis: 989aaa7dc30SBarry Smith #include <petscsys.h> 99091d3bdf4SKris Buschelman CHKMEMQ; 99130de9b25SBarry Smith 992eca87e8dSBarry Smith Not Collective 993eca87e8dSBarry Smith 99430de9b25SBarry Smith Level: beginner 99530de9b25SBarry Smith 99630de9b25SBarry Smith Notes: 997a17b96a8SKyle Gerard Felker We highly recommend using Valgrind https://petsc.org/release/faq/#valgrind or for NVIDIA CUDA systems 9985ed36255SBarry Smith https://docs.nvidia.com/cuda/cuda-memcheck/index.html for finding memory problems. The ``CHKMEMQ`` macro is useful on systems that 9995ed36255SBarry Smith do not have valgrind, but is not as good as valgrind or cuda-memcheck. 10001957e957SBarry Smith 100187497f52SBarry Smith Must run with the option -malloc_debug (-malloc_test in debug mode; or if `PetscMallocSetDebug()` called) to enable this option 100230de9b25SBarry Smith 100330de9b25SBarry Smith Once the error handler is called the calling function is then returned from with the given error code. 100430de9b25SBarry Smith 100530de9b25SBarry Smith By defaults prints location where memory that is corrupted was allocated. 100630de9b25SBarry Smith 100787497f52SBarry Smith Use `CHKMEMA` for functions that return void 1008f621e05eSBarry Smith 1009db781477SPatrick Sanan .seealso: `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `SETERRQ()`, `PetscMallocValidate()` 101030de9b25SBarry Smith M*/ 10116d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 1012064a246eSJacob Faibussowitsch #define CHKMEMQ 1013064a246eSJacob Faibussowitsch #define CHKMEMA 10146d210af2SJacob Faibussowitsch #else 10159371c9d4SSatish Balay #define CHKMEMQ \ 10169371c9d4SSatish Balay do { \ 101786d09637SLisandro Dalcin PetscErrorCode ierr_memq_ = PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__); \ 101886d09637SLisandro Dalcin if (PetscUnlikely(ierr_memq_)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_memq_, PETSC_ERROR_REPEAT, " "); \ 101986d09637SLisandro Dalcin } while (0) 10206d210af2SJacob Faibussowitsch #define CHKMEMA PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__) 1021064a246eSJacob Faibussowitsch #endif 10229566063dSJacob Faibussowitsch 1023668f157eSBarry Smith /*E 1024668f157eSBarry Smith PetscErrorType - passed to the PETSc error handling routines indicating if this is the first or a later call to the error handlers 1025668f157eSBarry Smith 1026668f157eSBarry Smith Level: advanced 1027668f157eSBarry Smith 102887497f52SBarry Smith `PETSC_ERROR_IN_CXX` indicates the error was detected in C++ and an exception should be generated 1029d736bfebSBarry Smith 103095452b02SPatrick Sanan Developer Notes: 103195452b02SPatrick Sanan This is currently used to decide when to print the detailed information about the run in PetscTraceBackErrorHandler() 1032668f157eSBarry Smith 103387497f52SBarry Smith .seealso: `PetscError()`, `SETERRQ()` 1034668f157eSBarry Smith E*/ 10359371c9d4SSatish Balay typedef enum { 10369371c9d4SSatish Balay PETSC_ERROR_INITIAL = 0, 10379371c9d4SSatish Balay PETSC_ERROR_REPEAT = 1, 10389371c9d4SSatish Balay PETSC_ERROR_IN_CXX = 2 10399371c9d4SSatish Balay } PetscErrorType; 10404b209cf6SBarry Smith 1041eb9e708aSLisandro Dalcin #if defined(__clang_analyzer__) 1042eb9e708aSLisandro Dalcin __attribute__((analyzer_noreturn)) 1043eb9e708aSLisandro Dalcin #endif 10449371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode 10459371c9d4SSatish Balay PetscError(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, ...) PETSC_ATTRIBUTE_COLD PETSC_ATTRIBUTE_FORMAT(7, 8); 1046eb9e708aSLisandro Dalcin 1047014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscErrorPrintfInitialize(void); 1048014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscErrorMessage(int, const char *[], char **); 1049d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscTraceBackErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1050d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscIgnoreErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1051d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscEmacsClientErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1052d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscMPIAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1053d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1054d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAttachDebuggerErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1055d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscReturnErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1056efca3c55SSatish Balay PETSC_EXTERN PetscErrorCode PetscPushErrorHandler(PetscErrorCode (*handler)(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *), void *); 1057014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopErrorHandler(void); 10588d359177SBarry Smith PETSC_EXTERN PetscErrorCode PetscSignalHandlerDefault(int, void *); 1059014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPushSignalHandler(PetscErrorCode (*)(int, void *), void *); 1060014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopSignalHandler(void); 106128559dc8SJed Brown PETSC_EXTERN PetscErrorCode PetscCheckPointerSetIntensity(PetscInt); 1062c2a741eeSJunchao Zhang PETSC_EXTERN void PetscSignalSegvCheckPointerOrMpi(void); 1063d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use PetscSignalSegvCheckPointerOrMpi() (since version 3.13)") static inline void PetscSignalSegvCheckPointer(void) 1064d71ae5a4SJacob Faibussowitsch { 10659371c9d4SSatish Balay PetscSignalSegvCheckPointerOrMpi(); 10669371c9d4SSatish Balay } 1067329f5518SBarry Smith 1068639ff905SBarry Smith /*MC 1069639ff905SBarry Smith PetscErrorPrintf - Prints error messages. 1070639ff905SBarry Smith 1071cf53795eSBarry Smith Not Collective; No Fortran Support 1072cf53795eSBarry Smith 1073639ff905SBarry Smith Synopsis: 1074aaa7dc30SBarry Smith #include <petscsys.h> 1075639ff905SBarry Smith PetscErrorCode (*PetscErrorPrintf)(const char format[],...); 1076639ff905SBarry Smith 1077f899ff85SJose E. Roman Input Parameter: 1078cd05f99aSJacob Faibussowitsch . format - the usual `printf()` format string 1079639ff905SBarry Smith 1080639ff905SBarry Smith Options Database Keys: 10811957e957SBarry Smith + -error_output_stdout - cause error messages to be printed to stdout instead of the (default) stderr 1082e1bc860dSBarry Smith - -error_output_none - to turn off all printing of error messages (does not change the way the error is handled.) 1083639ff905SBarry Smith 1084cf53795eSBarry Smith Level: developer 1085cf53795eSBarry Smith 108695452b02SPatrick Sanan Notes: 108795452b02SPatrick Sanan Use 1088639ff905SBarry Smith $ PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the 1089639ff905SBarry Smith $ error is handled.) and 10901957e957SBarry Smith $ PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on or you can use your own function 1091639ff905SBarry Smith 1092639ff905SBarry Smith Use 109387497f52SBarry Smith `PETSC_STDERR` = FILE* obtained from a file open etc. to have stderr printed to the file. 109487497f52SBarry Smith `PETSC_STDOUT` = FILE* obtained from a file open etc. to have stdout printed to the file. 1095639ff905SBarry Smith 1096639ff905SBarry Smith Use 109787497f52SBarry Smith `PetscPushErrorHandler()` to provide your own error handler that determines what kind of messages to print 1098639ff905SBarry Smith 1099db781477SPatrick Sanan .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscHelpPrintf()`, `PetscPrintf()`, `PetscPushErrorHandler()`, `PetscVFPrintf()`, `PetscHelpPrintf()` 1100639ff905SBarry Smith M*/ 11013ca90d2dSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode (*PetscErrorPrintf)(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2); 1102639ff905SBarry Smith 1103cf0818bdSBarry Smith /*E 1104cf0818bdSBarry Smith PetscFPTrap - types of floating point exceptions that may be trapped 1105cf0818bdSBarry Smith 110687497f52SBarry Smith Currently only `PETSC_FP_TRAP_OFF` and `PETSC_FP_TRAP_ON` are handled. All others are treated as `PETSC_FP_TRAP_ON`. 1107cf0818bdSBarry Smith 1108cf0818bdSBarry Smith Level: intermediate 1109cf0818bdSBarry Smith 1110cf0818bdSBarry Smith .seealso: `PetscSetFPTrap()`, `PetscPushFPTrap()` 1111cf0818bdSBarry Smith E*/ 11129371c9d4SSatish Balay typedef enum { 11139371c9d4SSatish Balay PETSC_FP_TRAP_OFF = 0, 11149371c9d4SSatish Balay PETSC_FP_TRAP_INDIV = 1, 11159371c9d4SSatish Balay PETSC_FP_TRAP_FLTOPERR = 2, 11169371c9d4SSatish Balay PETSC_FP_TRAP_FLTOVF = 4, 11179371c9d4SSatish Balay PETSC_FP_TRAP_FLTUND = 8, 11189371c9d4SSatish Balay PETSC_FP_TRAP_FLTDIV = 16, 11199371c9d4SSatish Balay PETSC_FP_TRAP_FLTINEX = 32 11209371c9d4SSatish Balay } PetscFPTrap; 1121bd2b07b1SBarry Smith #define PETSC_FP_TRAP_ON (PetscFPTrap)(PETSC_FP_TRAP_INDIV | PETSC_FP_TRAP_FLTOPERR | PETSC_FP_TRAP_FLTOVF | PETSC_FP_TRAP_FLTDIV | PETSC_FP_TRAP_FLTINEX) 1122014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSetFPTrap(PetscFPTrap); 1123014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPush(PetscFPTrap); 1124014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPop(void); 1125aba4c478SBarry Smith PETSC_EXTERN PetscErrorCode PetscDetermineInitialFPTrap(void); 112654a8ef01SBarry Smith 11273a40ed3dSBarry Smith /* 11283a40ed3dSBarry Smith Allows the code to build a stack frame as it runs 11293a40ed3dSBarry Smith */ 11303a40ed3dSBarry Smith 113199cd645aSJed Brown #define PETSCSTACKSIZE 64 11323a40ed3dSBarry Smith typedef struct { 11330e33f6ddSBarry Smith const char *function[PETSCSTACKSIZE]; 11340e33f6ddSBarry Smith const char *file[PETSCSTACKSIZE]; 1135184914b5SBarry Smith int line[PETSCSTACKSIZE]; 1136362febeeSStefano Zampini int petscroutine[PETSCSTACKSIZE]; /* 0 external called from petsc, 1 petsc functions, 2 petsc user functions */ 1137184914b5SBarry Smith int currentsize; 1138a2f94806SJed Brown int hotdepth; 11394be741a6SBarry Smith PetscBool check; /* option to check for correct Push/Pop semantics, true for default petscstack but not other stacks */ 11403a40ed3dSBarry Smith } PetscStack; 1141*dfb7d7afSStefano Zampini #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY) 114227104ee2SJacob Faibussowitsch PETSC_EXTERN PetscStack petscstack; 114327104ee2SJacob Faibussowitsch #endif 11443a40ed3dSBarry Smith 11455d12eec7SSatish Balay #if defined(PETSC_SERIALIZE_FUNCTIONS) 11465d12eec7SSatish Balay #include <petsc/private/petscfptimpl.h> 11475d12eec7SSatish Balay /* 11485d12eec7SSatish Balay Registers the current function into the global function pointer to function name table 11495d12eec7SSatish Balay 11505d12eec7SSatish Balay Have to fix this to handle errors but cannot return error since used in PETSC_VIEWER_DRAW_() etc 11515d12eec7SSatish Balay */ 11529371c9d4SSatish Balay #define PetscRegister__FUNCT__() \ 11539371c9d4SSatish Balay do { \ 11545d12eec7SSatish Balay static PetscBool __chked = PETSC_FALSE; \ 11555d12eec7SSatish Balay if (!__chked) { \ 11569371c9d4SSatish Balay void *ptr; \ 11579371c9d4SSatish Balay PetscDLSym(NULL, PETSC_FUNCTION_NAME, &ptr); \ 11585d12eec7SSatish Balay __chked = PETSC_TRUE; \ 11599371c9d4SSatish Balay } \ 11609371c9d4SSatish Balay } while (0) 11615d12eec7SSatish Balay #else 11625d12eec7SSatish Balay #define PetscRegister__FUNCT__() 11635d12eec7SSatish Balay #endif 11645d12eec7SSatish Balay 11656d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 11666d210af2SJacob Faibussowitsch #define PetscStackPushNoCheck(funct, petsc_routine, hot) 116737154d10SBarry Smith #define PetscStackUpdateLine 1168792fecdfSBarry Smith #define PetscStackPushExternal(funct) 11696d210af2SJacob Faibussowitsch #define PetscStackPopNoCheck 11706d210af2SJacob Faibussowitsch #define PetscStackClearTop 11716d210af2SJacob Faibussowitsch #define PetscFunctionBegin 11726d210af2SJacob Faibussowitsch #define PetscFunctionBeginUser 11736d210af2SJacob Faibussowitsch #define PetscFunctionBeginHot 11740a57284eSJacob Faibussowitsch #define PetscFunctionReturn(...) return __VA_ARGS__ 11756d210af2SJacob Faibussowitsch #define PetscFunctionReturnVoid() return 11766d210af2SJacob Faibussowitsch #define PetscStackPop 11776d210af2SJacob Faibussowitsch #define PetscStackPush(f) 1178*dfb7d7afSStefano Zampini #elif defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY) 1179660278c0SBarry Smith 11809371c9d4SSatish Balay #define PetscStackPush_Private(stack__, file__, func__, line__, petsc_routine__, hot__) \ 11819371c9d4SSatish Balay do { \ 11825a96b57dSJacob Faibussowitsch if (stack__.currentsize < PETSCSTACKSIZE) { \ 11835a96b57dSJacob Faibussowitsch stack__.function[stack__.currentsize] = func__; \ 1184ef1023bdSBarry Smith if (petsc_routine__) { \ 1185ef1023bdSBarry Smith stack__.file[stack__.currentsize] = file__; \ 11865a96b57dSJacob Faibussowitsch stack__.line[stack__.currentsize] = line__; \ 1187ef1023bdSBarry Smith } else { \ 1188648bc8c4SBarry Smith stack__.file[stack__.currentsize] = PETSC_NULLPTR; \ 1189ef1023bdSBarry Smith stack__.line[stack__.currentsize] = 0; \ 1190ef1023bdSBarry Smith } \ 11915a96b57dSJacob Faibussowitsch stack__.petscroutine[stack__.currentsize] = petsc_routine__; \ 11925a96b57dSJacob Faibussowitsch } \ 11935a96b57dSJacob Faibussowitsch ++stack__.currentsize; \ 11945a96b57dSJacob Faibussowitsch stack__.hotdepth += (hot__ || stack__.hotdepth); \ 11955a96b57dSJacob Faibussowitsch } while (0) 11965a96b57dSJacob Faibussowitsch 11974be741a6SBarry Smith /* uses PetscCheckAbort() because may be used in a function that does not return an error code */ 11989371c9d4SSatish Balay #define PetscStackPop_Private(stack__, func__) \ 11999371c9d4SSatish Balay do { \ 12004be741a6SBarry Smith PetscCheckAbort(!stack__.check || stack__.currentsize > 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid stack size %d, pop %s %s:%d.\n", stack__.currentsize, func__, __FILE__, __LINE__); \ 12015a96b57dSJacob Faibussowitsch if (--stack__.currentsize < PETSCSTACKSIZE) { \ 12029371c9d4SSatish Balay PetscCheckAbort(!stack__.check || stack__.petscroutine[stack__.currentsize] != 1 || stack__.function[stack__.currentsize] == (const char *)(func__), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid stack: push from %s %s:%d. Pop from %s %s:%d.\n", \ 12039371c9d4SSatish Balay stack__.function[stack__.currentsize], stack__.file[stack__.currentsize], stack__.line[stack__.currentsize], func__, __FILE__, __LINE__); \ 12045a96b57dSJacob Faibussowitsch stack__.function[stack__.currentsize] = PETSC_NULLPTR; \ 12055a96b57dSJacob Faibussowitsch stack__.file[stack__.currentsize] = PETSC_NULLPTR; \ 12065a96b57dSJacob Faibussowitsch stack__.line[stack__.currentsize] = 0; \ 12075a96b57dSJacob Faibussowitsch stack__.petscroutine[stack__.currentsize] = 0; \ 12085a96b57dSJacob Faibussowitsch } \ 12095a96b57dSJacob Faibussowitsch stack__.hotdepth = PetscMax(stack__.hotdepth - 1, 0); \ 12105a96b57dSJacob Faibussowitsch } while (0) 12115a96b57dSJacob Faibussowitsch 1212586f9135SBarry Smith /*MC 1213586f9135SBarry Smith PetscStackPushNoCheck - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is 1214586f9135SBarry Smith currently in the source code. 1215586f9135SBarry Smith 1216586f9135SBarry Smith Not Collective 1217586f9135SBarry Smith 1218586f9135SBarry Smith Synopsis: 1219586f9135SBarry Smith #include <petscsys.h> 1220586f9135SBarry Smith void PetscStackPushNoCheck(char *funct,int petsc_routine,PetscBool hot); 1221586f9135SBarry Smith 1222586f9135SBarry Smith Input Parameters: 1223586f9135SBarry Smith + funct - the function name 1224586f9135SBarry Smith . petsc_routine - 2 user function, 1 PETSc function, 0 some other function 1225586f9135SBarry Smith - hot - indicates that the function may be called often so expensive error checking should be turned off inside the function 1226586f9135SBarry Smith 1227586f9135SBarry Smith Level: developer 1228586f9135SBarry Smith 1229586f9135SBarry Smith Notes: 1230586f9135SBarry Smith In debug mode PETSc maintains a stack of the current function calls that can be used to help to quickly see where a problem has 123187497f52SBarry Smith occurred, for example, when a signal is received without running in the debugger. It is recommended to use the debugger if extensive information is needed to 1232586f9135SBarry Smith help debug the problem. 1233586f9135SBarry Smith 1234ef1023bdSBarry Smith This version does not check the memory corruption (an expensive operation), use `PetscStackPush()` to check the memory. 1235ef1023bdSBarry Smith 1236792fecdfSBarry Smith Use `PetscStackPushExternal()` for a function call that is about to be made to a non-PETSc or user function (such as BLAS etc). 1237ef1023bdSBarry Smith 1238586f9135SBarry Smith The default stack is a global variable called `petscstack`. 1239586f9135SBarry Smith 1240586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1241ef1023bdSBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPush()`, `PetscStackPop`, 1242792fecdfSBarry Smith `PetscStackPushExternal()` 1243586f9135SBarry Smith M*/ 12449371c9d4SSatish Balay #define PetscStackPushNoCheck(funct, petsc_routine, hot) \ 12459371c9d4SSatish Balay do { \ 1246e04113cfSBarry Smith PetscStackSAWsTakeAccess(); \ 12475a96b57dSJacob Faibussowitsch PetscStackPush_Private(petscstack, __FILE__, funct, __LINE__, petsc_routine, hot); \ 1248e04113cfSBarry Smith PetscStackSAWsGrantAccess(); \ 1249441dd030SJed Brown } while (0) 1250441dd030SJed Brown 1251586f9135SBarry Smith /*MC 125287497f52SBarry Smith PetscStackUpdateLine - in a function that has a `PetscFunctionBegin` or `PetscFunctionBeginUser` updates the stack line number to the 125337154d10SBarry Smith current line number. 125437154d10SBarry Smith 125537154d10SBarry Smith Not Collective 125637154d10SBarry Smith 125737154d10SBarry Smith Synopsis: 125837154d10SBarry Smith #include <petscsys.h> 125937154d10SBarry Smith void PetscStackUpdateLine 126037154d10SBarry Smith 126137154d10SBarry Smith Level: developer 126237154d10SBarry Smith 126337154d10SBarry Smith Notes: 126487497f52SBarry Smith Using `PetscCall()` and friends automatically handles this process 126587497f52SBarry Smith 126637154d10SBarry Smith In debug mode PETSc maintains a stack of the current function calls that can be used to help to quickly see where a problem has 126737154d10SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 126837154d10SBarry Smith help debug the problem. 126937154d10SBarry Smith 127037154d10SBarry Smith The default stack is a global variable called petscstack. 127137154d10SBarry Smith 127237154d10SBarry Smith This is used by `PetscCall()` and is otherwise not like to be needed 127337154d10SBarry Smith 127437154d10SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`, `PetscCall()` 127537154d10SBarry Smith M*/ 127637154d10SBarry Smith #define PetscStackUpdateLine \ 12779371c9d4SSatish Balay if (petscstack.currentsize > 0 && petscstack.function[petscstack.currentsize - 1] == PETSC_FUNCTION_NAME) { petscstack.line[petscstack.currentsize - 1] = __LINE__; } 127837154d10SBarry Smith 127937154d10SBarry Smith /*MC 1280792fecdfSBarry Smith PetscStackPushExternal - Pushes a new function name onto the PETSc default stack that tracks where the running program is 1281ef1023bdSBarry Smith currently in the source code. Does not include the filename or line number since this is called by the calling routine 1282ef1023bdSBarry Smith for non-PETSc or user functions. 1283ef1023bdSBarry Smith 1284ef1023bdSBarry Smith Not Collective 1285ef1023bdSBarry Smith 1286ef1023bdSBarry Smith Synopsis: 1287ef1023bdSBarry Smith #include <petscsys.h> 1288792fecdfSBarry Smith void PetscStackPushExternal(char *funct); 1289ef1023bdSBarry Smith 1290ef1023bdSBarry Smith Input Parameters: 1291ef1023bdSBarry Smith . funct - the function name 1292ef1023bdSBarry Smith 1293ef1023bdSBarry Smith Level: developer 1294ef1023bdSBarry Smith 1295ef1023bdSBarry Smith Notes: 129687497f52SBarry Smith Using `PetscCallExternal()` and friends automatically handles this process 129787497f52SBarry Smith 1298ef1023bdSBarry Smith In debug mode PETSc maintains a stack of the current function calls that can be used to help to quickly see where a problem has 1299ef1023bdSBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1300ef1023bdSBarry Smith help debug the problem. 1301ef1023bdSBarry Smith 1302ef1023bdSBarry Smith The default stack is a global variable called `petscstack`. 1303ef1023bdSBarry Smith 1304ef1023bdSBarry Smith This is to be used when calling an external package function such as a BLAS function. 1305ef1023bdSBarry Smith 1306ef1023bdSBarry Smith This also updates the stack line number for the current stack function. 1307ef1023bdSBarry Smith 1308ef1023bdSBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1309ef1023bdSBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1310ef1023bdSBarry Smith M*/ 13119371c9d4SSatish Balay #define PetscStackPushExternal(funct) \ 13129371c9d4SSatish Balay do { \ 13139371c9d4SSatish Balay PetscStackUpdateLine; \ 13149371c9d4SSatish Balay PetscStackPushNoCheck(funct, 0, PETSC_TRUE); \ 13159371c9d4SSatish Balay } while (0); 1316ef1023bdSBarry Smith 1317ef1023bdSBarry Smith /*MC 1318586f9135SBarry Smith PetscStackPopNoCheck - Pops a function name from the PETSc default stack that tracks where the running program is 1319586f9135SBarry Smith currently in the source code. 1320586f9135SBarry Smith 1321586f9135SBarry Smith Not Collective 1322586f9135SBarry Smith 1323586f9135SBarry Smith Synopsis: 1324586f9135SBarry Smith #include <petscsys.h> 1325586f9135SBarry Smith void PetscStackPopNoCheck(char *funct); 1326586f9135SBarry Smith 1327586f9135SBarry Smith Input Parameter: 1328586f9135SBarry Smith . funct - the function name 1329586f9135SBarry Smith 1330586f9135SBarry Smith Level: developer 1331586f9135SBarry Smith 1332586f9135SBarry Smith Notes: 133387497f52SBarry Smith Using `PetscCall()`, `PetscCallExternal()`, `PetscCallBack()` and friends negates the need to call this 133487497f52SBarry Smith 1335586f9135SBarry Smith In debug mode PETSc maintains a stack of the current function calls that can be used to help to quickly see where a problem has 1336586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1337586f9135SBarry Smith help debug the problem. 1338586f9135SBarry Smith 1339586f9135SBarry Smith The default stack is a global variable called petscstack. 1340586f9135SBarry Smith 1341586f9135SBarry Smith Developer Note: 1342586f9135SBarry Smith `PetscStackPopNoCheck()` takes a function argument while `PetscStackPop` does not, this difference is likely just historical. 1343586f9135SBarry Smith 1344586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1345586f9135SBarry Smith M*/ 13469371c9d4SSatish Balay #define PetscStackPopNoCheck(funct) \ 13479371c9d4SSatish Balay do { \ 1348362febeeSStefano Zampini PetscStackSAWsTakeAccess(); \ 13495a96b57dSJacob Faibussowitsch PetscStackPop_Private(petscstack, funct); \ 1350362febeeSStefano Zampini PetscStackSAWsGrantAccess(); \ 1351362febeeSStefano Zampini } while (0) 1352362febeeSStefano Zampini 13539371c9d4SSatish Balay #define PetscStackClearTop \ 13549371c9d4SSatish Balay do { \ 1355e04113cfSBarry Smith PetscStackSAWsTakeAccess(); \ 13569371c9d4SSatish Balay if (petscstack.currentsize > 0 && --petscstack.currentsize < PETSCSTACKSIZE) { \ 135727104ee2SJacob Faibussowitsch petscstack.function[petscstack.currentsize] = PETSC_NULLPTR; \ 135827104ee2SJacob Faibussowitsch petscstack.file[petscstack.currentsize] = PETSC_NULLPTR; \ 135927104ee2SJacob Faibussowitsch petscstack.line[petscstack.currentsize] = 0; \ 136027104ee2SJacob Faibussowitsch petscstack.petscroutine[petscstack.currentsize] = 0; \ 1361441dd030SJed Brown } \ 136227104ee2SJacob Faibussowitsch petscstack.hotdepth = PetscMax(petscstack.hotdepth - 1, 0); \ 1363e04113cfSBarry Smith PetscStackSAWsGrantAccess(); \ 1364441dd030SJed Brown } while (0) 1365441dd030SJed Brown 136630de9b25SBarry Smith /*MC 13671957e957SBarry Smith PetscFunctionBegin - First executable line of each PETSc function, used for error handling. Final 136887497f52SBarry Smith line of PETSc functions should be `PetscFunctionReturn`(0); 136930de9b25SBarry Smith 137030de9b25SBarry Smith Synopsis: 1371aaa7dc30SBarry Smith #include <petscsys.h> 137230de9b25SBarry Smith void PetscFunctionBegin; 137330de9b25SBarry Smith 1374eca87e8dSBarry Smith Not Collective 1375eca87e8dSBarry Smith 137630de9b25SBarry Smith Usage: 137730de9b25SBarry Smith .vb 137830de9b25SBarry Smith int something; 137930de9b25SBarry Smith 138030de9b25SBarry Smith PetscFunctionBegin; 138130de9b25SBarry Smith .ve 138230de9b25SBarry Smith 138330de9b25SBarry Smith Notes: 138487497f52SBarry Smith Use `PetscFunctionBeginUser` for application codes. 13851957e957SBarry Smith 138630de9b25SBarry Smith Not available in Fortran 138730de9b25SBarry Smith 138830de9b25SBarry Smith Level: developer 138930de9b25SBarry Smith 1390586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()` 139130de9b25SBarry Smith 139230de9b25SBarry Smith M*/ 13939371c9d4SSatish Balay #define PetscFunctionBegin \ 13949371c9d4SSatish Balay do { \ 1395362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_FALSE); \ 1396a2f94806SJed Brown PetscRegister__FUNCT__(); \ 1397a2f94806SJed Brown } while (0) 1398a2f94806SJed Brown 1399a2f94806SJed Brown /*MC 140087497f52SBarry Smith PetscFunctionBeginHot - Substitute for `PetscFunctionBegin` to be used in functions that are called in 1401a2f94806SJed Brown performance-critical circumstances. Use of this function allows for lighter profiling by default. 1402a2f94806SJed Brown 1403a2f94806SJed Brown Synopsis: 1404aaa7dc30SBarry Smith #include <petscsys.h> 1405a2f94806SJed Brown void PetscFunctionBeginHot; 1406a2f94806SJed Brown 1407a2f94806SJed Brown Not Collective 1408a2f94806SJed Brown 1409a2f94806SJed Brown Usage: 1410a2f94806SJed Brown .vb 1411a2f94806SJed Brown int something; 1412a2f94806SJed Brown 1413a2f94806SJed Brown PetscFunctionBeginHot; 1414a2f94806SJed Brown .ve 1415a2f94806SJed Brown 1416a2f94806SJed Brown Notes: 1417a2f94806SJed Brown Not available in Fortran 1418a2f94806SJed Brown 1419a2f94806SJed Brown Level: developer 1420a2f94806SJed Brown 1421586f9135SBarry Smith .seealso: `PetscFunctionBegin`, `PetscFunctionReturn()`, `PetscStackPushNoCheck()` 1422a2f94806SJed Brown 1423a2f94806SJed Brown M*/ 14249371c9d4SSatish Balay #define PetscFunctionBeginHot \ 14259371c9d4SSatish Balay do { \ 1426362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_TRUE); \ 14272d53ad75SBarry Smith PetscRegister__FUNCT__(); \ 142853c77d0aSJed Brown } while (0) 142953c77d0aSJed Brown 1430a8d2bbe5SBarry Smith /*MC 1431530d85c1SBarry Smith PetscFunctionBeginUser - First executable line of user provided routines 1432a8d2bbe5SBarry Smith 1433a8d2bbe5SBarry Smith Synopsis: 1434aaa7dc30SBarry Smith #include <petscsys.h> 1435a8d2bbe5SBarry Smith void PetscFunctionBeginUser; 1436a8d2bbe5SBarry Smith 1437a8d2bbe5SBarry Smith Not Collective 1438a8d2bbe5SBarry Smith 1439a8d2bbe5SBarry Smith Usage: 1440a8d2bbe5SBarry Smith .vb 1441a8d2bbe5SBarry Smith int something; 1442a8d2bbe5SBarry Smith 1443ac285190SBarry Smith PetscFunctionBeginUser; 1444a8d2bbe5SBarry Smith .ve 1445a8d2bbe5SBarry Smith 1446a8d2bbe5SBarry Smith Notes: 1447530d85c1SBarry Smith Functions that incorporate this must call `PetscFunctionReturn()` instead of return except for main(). 1448530d85c1SBarry Smith 1449530d85c1SBarry Smith May be used before `PetscInitialize()` 14501957e957SBarry Smith 1451a8d2bbe5SBarry Smith Not available in Fortran 1452a8d2bbe5SBarry Smith 1453530d85c1SBarry Smith This is identical to `PetscFunctionBegin` except it labels the routine as a user 1454ac285190SBarry Smith routine instead of as a PETSc library routine. 1455ac285190SBarry Smith 1456a2f94806SJed Brown Level: intermediate 1457a8d2bbe5SBarry Smith 1458586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, `PetscFunctionBeginHot`, `PetscStackPushNoCheck()` 1459a8d2bbe5SBarry Smith 1460a8d2bbe5SBarry Smith M*/ 14619371c9d4SSatish Balay #define PetscFunctionBeginUser \ 14629371c9d4SSatish Balay do { \ 1463362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 2, PETSC_FALSE); \ 1464a8d2bbe5SBarry Smith PetscRegister__FUNCT__(); \ 1465a8d2bbe5SBarry Smith } while (0) 1466a8d2bbe5SBarry Smith 1467586f9135SBarry Smith /*MC 1468586f9135SBarry Smith PetscStackPush - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is 1469586f9135SBarry Smith currently in the source code and verifies the memory is not corrupted. 1470586f9135SBarry Smith 1471586f9135SBarry Smith Not Collective 1472586f9135SBarry Smith 1473586f9135SBarry Smith Synopsis: 1474586f9135SBarry Smith #include <petscsys.h> 1475586f9135SBarry Smith void PetscStackPush(char *funct) 1476586f9135SBarry Smith 1477586f9135SBarry Smith Input Parameter: 1478586f9135SBarry Smith . funct - the function name 1479586f9135SBarry Smith 1480586f9135SBarry Smith Level: developer 1481586f9135SBarry Smith 1482586f9135SBarry Smith Notes: 1483586f9135SBarry Smith In debug mode PETSc maintains a stack of the current function calls that can be used to help to quickly see where a problem has 1484586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1485586f9135SBarry Smith help debug the problem. 1486586f9135SBarry Smith 1487586f9135SBarry Smith The default stack is a global variable called petscstack. 1488586f9135SBarry Smith 1489586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1490586f9135SBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1491586f9135SBarry Smith M*/ 14929371c9d4SSatish Balay #define PetscStackPush(n) \ 14939371c9d4SSatish Balay do { \ 1494362febeeSStefano Zampini PetscStackPushNoCheck(n, 0, PETSC_FALSE); \ 149515681b3cSBarry Smith CHKMEMQ; \ 149615681b3cSBarry Smith } while (0) 14973a40ed3dSBarry Smith 1498586f9135SBarry Smith /*MC 1499586f9135SBarry Smith PetscStackPop - Pops a function name from the PETSc default stack that tracks where the running program is 1500586f9135SBarry Smith currently in the source code and verifies the memory is not corrupted. 1501586f9135SBarry Smith 1502586f9135SBarry Smith Not Collective 1503586f9135SBarry Smith 1504586f9135SBarry Smith Synopsis: 1505586f9135SBarry Smith #include <petscsys.h> 1506586f9135SBarry Smith void PetscStackPop 1507586f9135SBarry Smith 1508586f9135SBarry Smith Level: developer 1509586f9135SBarry Smith 1510586f9135SBarry Smith Notes: 1511586f9135SBarry Smith In debug mode PETSc maintains a stack of the current function calls that can be used to help to quickly see where a problem has 1512586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1513586f9135SBarry Smith help debug the problem. 1514586f9135SBarry Smith 1515586f9135SBarry Smith The default stack is a global variable called petscstack. 1516586f9135SBarry Smith 1517586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPopNoCheck()`, `PetscStackPush()` 1518586f9135SBarry Smith M*/ 15199371c9d4SSatish Balay #define PetscStackPop \ 15209371c9d4SSatish Balay do { \ 1521441dd030SJed Brown CHKMEMQ; \ 1522362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 152315681b3cSBarry Smith } while (0) 1524d64ed03dSBarry Smith 152530de9b25SBarry Smith /*MC 15260a57284eSJacob Faibussowitsch PetscFunctionReturn - Last executable line of each PETSc function used for error 15270a57284eSJacob Faibussowitsch handling. Replaces `return()`. 152830de9b25SBarry Smith 152930de9b25SBarry Smith Synopsis: 15300a57284eSJacob Faibussowitsch #include <petscerror.h> 15310a57284eSJacob Faibussowitsch void PetscFunctionReturn(...) 153230de9b25SBarry Smith 1533eca87e8dSBarry Smith Not Collective 1534eca87e8dSBarry Smith 15355687f061SJacob Faibussowitsch Level: beginner 15365687f061SJacob Faibussowitsch 15370a57284eSJacob Faibussowitsch Notes: 15380a57284eSJacob Faibussowitsch This routine is a macro, so while it does not "return" anything itself, it does return from 15390a57284eSJacob Faibussowitsch the function in the literal sense. 15400a57284eSJacob Faibussowitsch 15410a57284eSJacob Faibussowitsch Usually the return value is the integer literal `0` (for example in any function returning 15420a57284eSJacob Faibussowitsch `PetscErrorCode`), however it is possible to return any arbitrary type. The arguments of 15430a57284eSJacob Faibussowitsch this macro are placed before the `return` statement as-is. 15440a57284eSJacob Faibussowitsch 15450a57284eSJacob Faibussowitsch Any routine which returns via `PetscFunctionReturn()` must begin with a corresponding 15460a57284eSJacob Faibussowitsch `PetscFunctionBegin`. 15470a57284eSJacob Faibussowitsch 15480a57284eSJacob Faibussowitsch For routines which return `void` use `PetscFunctionReturnVoid()` instead. 15490a57284eSJacob Faibussowitsch 15500a57284eSJacob Faibussowitsch Example Usage: 155130de9b25SBarry Smith .vb 15520a57284eSJacob Faibussowitsch PetscErrorCode foo(int *x) 15530a57284eSJacob Faibussowitsch { 15540a57284eSJacob Faibussowitsch PetscFunctionBegin; // don't forget the begin! 15550a57284eSJacob Faibussowitsch *x = 10; 155630de9b25SBarry Smith PetscFunctionReturn(0); 155730de9b25SBarry Smith } 155830de9b25SBarry Smith .ve 155930de9b25SBarry Smith 15600a57284eSJacob Faibussowitsch May return any arbitrary type\: 15610a57284eSJacob Faibussowitsch .vb 15620a57284eSJacob Faibussowitsch struct Foo 15630a57284eSJacob Faibussowitsch { 15640a57284eSJacob Faibussowitsch int x; 15650a57284eSJacob Faibussowitsch }; 15660a57284eSJacob Faibussowitsch 15670a57284eSJacob Faibussowitsch struct Foo make_foo(int value) 15680a57284eSJacob Faibussowitsch { 15690a57284eSJacob Faibussowitsch struct Foo f; 15700a57284eSJacob Faibussowitsch 15710a57284eSJacob Faibussowitsch PetscFunctionBegin; 15720a57284eSJacob Faibussowitsch f.x = value; 15730a57284eSJacob Faibussowitsch PetscFunctionReturn(f); 15740a57284eSJacob Faibussowitsch } 15750a57284eSJacob Faibussowitsch .ve 15760a57284eSJacob Faibussowitsch 15775687f061SJacob Faibussowitsch Fortran Note: 157830de9b25SBarry Smith Not available in Fortran 157930de9b25SBarry Smith 15800a57284eSJacob Faibussowitsch .seealso: `PetscFunctionBegin`, `PetscFunctionBeginUser`, `PetscFunctionReturnVoid()`, 15810a57284eSJacob Faibussowitsch `PetscStackPopNoCheck()` 158230de9b25SBarry Smith M*/ 15835687f061SJacob Faibussowitsch #define PetscFunctionReturn(...) \ 15849371c9d4SSatish Balay do { \ 1585362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 15865687f061SJacob Faibussowitsch return __VA_ARGS__; \ 158727104ee2SJacob Faibussowitsch } while (0) 1588d64ed03dSBarry Smith 15890a57284eSJacob Faibussowitsch /*MC 15900a57284eSJacob Faibussowitsch PetscFunctionReturnVoid - Like `PetscFunctionReturn()` but returns `void` 15910a57284eSJacob Faibussowitsch 15920a57284eSJacob Faibussowitsch Synopsis: 15930a57284eSJacob Faibussowitsch #include <petscerror.h> 15940a57284eSJacob Faibussowitsch void PetscFunctionReturnVoid() 15950a57284eSJacob Faibussowitsch 15960a57284eSJacob Faibussowitsch Not Collective 15970a57284eSJacob Faibussowitsch 15985687f061SJacob Faibussowitsch Level: beginner 15995687f061SJacob Faibussowitsch 16005687f061SJacob Faibussowitsch Note: 16010a57284eSJacob Faibussowitsch Behaves identically to `PetscFunctionReturn()` except that it returns `void`. That is, this 16025687f061SJacob Faibussowitsch macro culminates with `return`. 16030a57284eSJacob Faibussowitsch 16040a57284eSJacob Faibussowitsch Example Usage: 16050a57284eSJacob Faibussowitsch .vb 16060a57284eSJacob Faibussowitsch void foo() 16070a57284eSJacob Faibussowitsch { 16080a57284eSJacob Faibussowitsch PetscFunctionBegin; // must start with PetscFunctionBegin! 16090a57284eSJacob Faibussowitsch bar(); 16100a57284eSJacob Faibussowitsch baz(); 16110a57284eSJacob Faibussowitsch PetscFunctionReturnVoid(); 16120a57284eSJacob Faibussowitsch } 16130a57284eSJacob Faibussowitsch .ve 16140a57284eSJacob Faibussowitsch 16150a57284eSJacob Faibussowitsch .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, PetscFunctionBeginUser` 16160a57284eSJacob Faibussowitsch M*/ 16179371c9d4SSatish Balay #define PetscFunctionReturnVoid() \ 16189371c9d4SSatish Balay do { \ 1619362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 162027104ee2SJacob Faibussowitsch return; \ 162127104ee2SJacob Faibussowitsch } while (0) 162227104ee2SJacob Faibussowitsch #else /* PETSC_USE_DEBUG */ 162327104ee2SJacob Faibussowitsch #define PetscStackPushNoCheck(funct, petsc_routine, hot) 162437154d10SBarry Smith #define PetscStackUpdateLine 1625792fecdfSBarry Smith #define PetscStackPushExternal(funct) 162627104ee2SJacob Faibussowitsch #define PetscStackPopNoCheck 162727104ee2SJacob Faibussowitsch #define PetscStackClearTop 16283a40ed3dSBarry Smith #define PetscFunctionBegin 16290bdf7c52SPeter Brune #define PetscFunctionBeginUser 1630a2f94806SJed Brown #define PetscFunctionBeginHot 16310a57284eSJacob Faibussowitsch #define PetscFunctionReturn(...) return __VA_ARGS__ 16325665465eSBarry Smith #define PetscFunctionReturnVoid() return 1633812af9f3SBarry Smith #define PetscStackPop CHKMEMQ 1634812af9f3SBarry Smith #define PetscStackPush(f) CHKMEMQ 163527104ee2SJacob Faibussowitsch #endif /* PETSC_USE_DEBUG */ 16363a40ed3dSBarry Smith 16376d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 1638e77caa6dSBarry Smith #define PetscStackCallExternalVoid(name, routine) 1639792fecdfSBarry Smith #define PetscCallExternal(func, ...) 16406d210af2SJacob Faibussowitsch #else 1641586f9135SBarry Smith /*MC 1642e77caa6dSBarry Smith PetscStackCallExternalVoid - Calls an external library routine or user function after pushing the name of the routine on the stack. 1643eb6b5d47SBarry Smith 1644eb6b5d47SBarry Smith Input Parameters: 1645eb6b5d47SBarry Smith + name - string that gives the name of the function being called 1646586f9135SBarry Smith - routine - actual call to the routine, for example, functionname(a,b) 1647fd3f9acdSBarry Smith 1648586f9135SBarry Smith Level: developer 1649eb6b5d47SBarry Smith 1650586f9135SBarry Smith Note: 1651792fecdfSBarry Smith Often one should use `PetscCallExternal()` instead. This routine is intended for external library routines that DO NOT return error codes 1652eb6b5d47SBarry Smith 1653586f9135SBarry Smith In debug mode this also checks the memory for corruption at the end of the function call. 1654586f9135SBarry Smith 1655586f9135SBarry Smith Certain external packages, such as BLAS/LAPACK may have their own macros for managing the call, error checking, etc. 1656586f9135SBarry Smith 1657586f9135SBarry Smith Developer Note: 1658586f9135SBarry Smith This is so that when a user or external library routine results in a crash or corrupts memory, they get blamed instead of PETSc. 1659586f9135SBarry Smith 1660792fecdfSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscCallExternal()`, `PetscCallBLAS()` 1661586f9135SBarry Smith @*/ 16629371c9d4SSatish Balay #define PetscStackCallExternalVoid(name, routine) \ 16639371c9d4SSatish Balay do { \ 16649371c9d4SSatish Balay PetscStackPush(name); \ 16659371c9d4SSatish Balay routine; \ 16669371c9d4SSatish Balay PetscStackPop; \ 16679371c9d4SSatish Balay } while (0) 1668eb6b5d47SBarry Smith 1669586f9135SBarry Smith /*MC 1670792fecdfSBarry Smith PetscCallExternal - Calls an external library routine that returns an error code after pushing the name of the routine on the stack. 1671fd3f9acdSBarry Smith 1672fd3f9acdSBarry Smith Input Parameters: 1673fd3f9acdSBarry Smith + func- name of the routine 1674586f9135SBarry Smith - args - arguments to the routine 1675586f9135SBarry Smith 1676586f9135SBarry Smith Level: developer 1677fd3f9acdSBarry Smith 167895452b02SPatrick Sanan Notes: 1679e77caa6dSBarry Smith This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not. 1680dbf62e16SBarry Smith 1681586f9135SBarry Smith In debug mode this also checks the memory for corruption at the end of the function call. 1682fd3f9acdSBarry Smith 168387497f52SBarry Smith Assumes the error return code of the function is an integer and that a value of 0 indicates success 168487497f52SBarry Smith 1685586f9135SBarry Smith Developer Note: 1686d5b43468SJose E. Roman This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc. 1687586f9135SBarry Smith 1688e77caa6dSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()` 1689586f9135SBarry Smith M*/ 16909371c9d4SSatish Balay #define PetscCallExternal(func, ...) \ 16919371c9d4SSatish Balay do { \ 1692a74df02fSJacob Faibussowitsch PetscStackPush(PetscStringize(func)); \ 1693a74df02fSJacob Faibussowitsch PetscErrorCode __ierr = func(__VA_ARGS__); \ 16941d4906efSStefano Zampini PetscStackPop; \ 169532771afcSJacob Faibussowitsch PetscCheck(!__ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in %s(): error code %d", PetscStringize(func), __ierr); \ 1696fd3f9acdSBarry Smith } while (0) 16976d210af2SJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */ 169806d1fe2cSBarry Smith 169906d1fe2cSBarry Smith #endif 1700