154a8ef01SBarry Smith /* 2f621e05eSBarry Smith Contains all error handling interfaces for PETSc. 354a8ef01SBarry Smith */ 4a4963045SJacob Faibussowitsch #pragma once 56c7e564aSBarry Smith 6e1bf4ed2SJacob Faibussowitsch #include <petscmacros.h> 7e1bf4ed2SJacob Faibussowitsch #include <petscsystypes.h> 8e1bf4ed2SJacob Faibussowitsch 9af75f258SJacob Faibussowitsch #if defined(__cplusplus) 10af75f258SJacob Faibussowitsch #include <exception> // std::exception 11af75f258SJacob Faibussowitsch #endif 12af75f258SJacob Faibussowitsch 13ac09b921SBarry Smith /* SUBMANSEC = Sys */ 14ac09b921SBarry Smith 15edd03b47SJacob Faibussowitsch #define SETERRQ1(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 16edd03b47SJacob Faibussowitsch #define SETERRQ2(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 17edd03b47SJacob Faibussowitsch #define SETERRQ3(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 18edd03b47SJacob Faibussowitsch #define SETERRQ4(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 19edd03b47SJacob Faibussowitsch #define SETERRQ5(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 20edd03b47SJacob Faibussowitsch #define SETERRQ6(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 21edd03b47SJacob Faibussowitsch #define SETERRQ7(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 22edd03b47SJacob Faibussowitsch #define SETERRQ8(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 23edd03b47SJacob Faibussowitsch #define SETERRQ9(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 2498921bdaSJacob Faibussowitsch 2530de9b25SBarry Smith /*MC 261957e957SBarry Smith SETERRQ - Macro to be called when an error has been detected, 2730de9b25SBarry Smith 2830de9b25SBarry Smith Synopsis: 29aaa7dc30SBarry Smith #include <petscsys.h> 3098921bdaSJacob Faibussowitsch PetscErrorCode SETERRQ(MPI_Comm comm, PetscErrorCode ierr, char *message, ...) 3130de9b25SBarry Smith 32d083f849SBarry Smith Collective 3330de9b25SBarry Smith 3430de9b25SBarry Smith Input Parameters: 3595bd0b28SBarry Smith + comm - An MPI communicator, use `PETSC_COMM_SELF` unless you know all ranks of another communicator will detect the error 363af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 3730de9b25SBarry Smith - message - error message 3830de9b25SBarry Smith 3930de9b25SBarry Smith Level: beginner 4030de9b25SBarry Smith 4130de9b25SBarry Smith Notes: 4287497f52SBarry Smith This is rarely needed, one should use `PetscCheck()` and `PetscCall()` and friends to automatically handle error conditions. 4330de9b25SBarry Smith Once the error handler is called the calling function is then returned from with the given error code. 4430de9b25SBarry Smith 4587497f52SBarry Smith Experienced users can set the error handler with `PetscPushErrorHandler()`. 4630de9b25SBarry Smith 4716a05f60SBarry Smith Fortran Note: 48989c49a2SBarry Smith `SETERRQ()` may be called from Fortran subroutines but `SETERRA()` must be called from the 493b1008b8SBarry Smith Fortran main program. 503b1008b8SBarry Smith 51db781477SPatrick Sanan .seealso: `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 52af27ebaaSBarry Smith `PetscError()`, `PetscCall()`, `CHKMEMQ`, `CHKERRA()`, `PetscCallMPI()`, `PetscErrorCode` 5330de9b25SBarry Smith M*/ 54e809461dSJacob Faibussowitsch #define SETERRQ(comm, ierr, ...) \ 55e809461dSJacob Faibussowitsch do { \ 56e809461dSJacob Faibussowitsch PetscErrorCode ierr_seterrq_petsc_ = PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \ 57e809461dSJacob Faibussowitsch return ierr_seterrq_petsc_ ? ierr_seterrq_petsc_ : PETSC_ERR_RETURN; \ 58e809461dSJacob Faibussowitsch } while (0) 59986eef2eSBarry Smith 60648c30bcSBarry Smith #define SETERRQNULL(comm, ierr, ...) \ 61648c30bcSBarry Smith do { \ 62648c30bcSBarry Smith (void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \ 63648c30bcSBarry Smith return NULL; \ 64648c30bcSBarry Smith } while (0) 65648c30bcSBarry Smith 66ffc4695bSBarry Smith /* 67ffc4695bSBarry Smith Returned from PETSc functions that are called from MPI, such as related to attributes 68ffc4695bSBarry Smith Do not confuse PETSC_MPI_ERROR_CODE and PETSC_ERR_MPI, the first is registered with MPI and returned to MPI as 69888a1fb3SBarry 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. 70ffc4695bSBarry Smith */ 71ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CLASS; 72ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CODE; 73ffc4695bSBarry Smith 74986eef2eSBarry Smith /*MC 75986eef2eSBarry Smith SETERRMPI - Macro to be called when an error has been detected within an MPI callback function 76986eef2eSBarry Smith 77989c49a2SBarry Smith No Fortran Support 78989c49a2SBarry Smith 79986eef2eSBarry Smith Synopsis: 80986eef2eSBarry Smith #include <petscsys.h> 8198921bdaSJacob Faibussowitsch PetscErrorCode SETERRMPI(MPI_Comm comm, PetscErrorCode ierr, char *message, ...) 82986eef2eSBarry Smith 83d083f849SBarry Smith Collective 84986eef2eSBarry Smith 85986eef2eSBarry Smith Input Parameters: 8695bd0b28SBarry Smith + comm - An MPI communicator, use `PETSC_COMM_SELF` unless you know all ranks of another communicator will detect the error 87986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 88986eef2eSBarry Smith - message - error message 89986eef2eSBarry Smith 90986eef2eSBarry Smith Level: developer 91986eef2eSBarry Smith 9295bd0b28SBarry Smith Note: 9387497f52SBarry 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` 9487497f52SBarry Smith which is registered with `MPI_Add_error_code()` when PETSc is initialized. 95986eef2eSBarry Smith 96af27ebaaSBarry Smith .seealso: `SETERRQ()`, `PetscCall()`, `PetscCallMPI()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `PetscErrorCode` 97986eef2eSBarry Smith M*/ 983ba16761SJacob Faibussowitsch #define SETERRMPI(comm, ierr, ...) return ((void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__), PETSC_MPI_ERROR_CODE) 99ee8199e6SMatthew G. Knepley 100ee8199e6SMatthew G. Knepley /*MC 101f388eb8bSPatrick Sanan SETERRA - Fortran-only macro that can be called when an error has been detected from the main program 102f388eb8bSPatrick Sanan 103f388eb8bSPatrick Sanan Synopsis: 104f388eb8bSPatrick Sanan #include <petscsys.h> 105f388eb8bSPatrick Sanan PetscErrorCode SETERRA(MPI_Comm comm, PetscErrorCode ierr, char *message) 106f388eb8bSPatrick Sanan 107f388eb8bSPatrick Sanan Collective 108f388eb8bSPatrick Sanan 109f388eb8bSPatrick Sanan Input Parameters: 11095bd0b28SBarry Smith + comm - An MPI communicator, so that the error can be collective 111f388eb8bSPatrick Sanan . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 112e55b9c91SBarry Smith - message - error message in the `printf()` format 113f388eb8bSPatrick Sanan 114f388eb8bSPatrick Sanan Level: beginner 115f388eb8bSPatrick Sanan 116f388eb8bSPatrick Sanan Notes: 11787497f52SBarry Smith This should only be used with Fortran. With C/C++, use `SETERRQ()`. 118f388eb8bSPatrick Sanan 11987497f52SBarry Smith `SETERRQ()` may be called from Fortran subroutines but `SETERRA()` must be called from the 120f388eb8bSPatrick Sanan Fortran main program. 121f388eb8bSPatrick Sanan 122af27ebaaSBarry Smith .seealso: `SETERRQ()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()`, `PetscErrorCode` 123f388eb8bSPatrick Sanan M*/ 124f388eb8bSPatrick Sanan 125f388eb8bSPatrick Sanan /*MC 126fa190f98SMatthew G. Knepley SETERRABORT - Macro that can be called when an error has been detected, 127fa190f98SMatthew G. Knepley 128fa190f98SMatthew G. Knepley Synopsis: 129fa190f98SMatthew G. Knepley #include <petscsys.h> 13098921bdaSJacob Faibussowitsch PetscErrorCode SETERRABORT(MPI_Comm comm, PetscErrorCode ierr, char *message, ...) 131fa190f98SMatthew G. Knepley 132d083f849SBarry Smith Collective 133fa190f98SMatthew G. Knepley 134fa190f98SMatthew G. Knepley Input Parameters: 13595bd0b28SBarry Smith + comm - An MPI communicator, so that the error can be collective 1363af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 137e55b9c91SBarry Smith - message - error message in the `printf()` format 138fa190f98SMatthew G. Knepley 139fa190f98SMatthew G. Knepley Level: beginner 140fa190f98SMatthew G. Knepley 141fa190f98SMatthew G. Knepley Notes: 14287497f52SBarry Smith This function just calls `MPI_Abort()`. 14387497f52SBarry Smith 14487497f52SBarry Smith This should only be called in routines that cannot return an error code, such as in C++ constructors. 145fa190f98SMatthew G. Knepley 146989c49a2SBarry Smith Fortran Note: 147989c49a2SBarry Smith Use `SETERRA()` in Fortran main program and `SETERRQ()` in Fortran subroutines 148989c49a2SBarry Smith 149989c49a2SBarry Smith Developer Note: 150989c49a2SBarry Smith In Fortran `SETERRA()` could be called `SETERRABORT()` since they serve the same purpose 151989c49a2SBarry Smith 152af27ebaaSBarry Smith .seealso: `SETERRQ()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `PetscCall()`, `CHKMEMQ`, `PetscErrorCode` 153fa190f98SMatthew G. Knepley M*/ 1549371c9d4SSatish Balay #define SETERRABORT(comm, ierr, ...) \ 1559371c9d4SSatish Balay do { \ 1563ba16761SJacob Faibussowitsch (void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \ 157dd460d27SBarry Smith (void)MPI_Abort(comm, ierr); \ 15898921bdaSJacob Faibussowitsch } while (0) 1599a00fa46SSatish Balay 16030de9b25SBarry Smith /*MC 161139c8099SSatish Balay PetscCheck - Checks that a particular condition is true; if not true, then returns the provided error code 1622c71b3e2SJacob Faibussowitsch 1632c71b3e2SJacob Faibussowitsch Synopsis: 1642c71b3e2SJacob Faibussowitsch #include <petscerror.h> 1652c71b3e2SJacob Faibussowitsch void PetscCheck(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 1662c71b3e2SJacob Faibussowitsch 1677cdbe19fSJose E. Roman Collective; No Fortran Support 1682c71b3e2SJacob Faibussowitsch 1692c71b3e2SJacob Faibussowitsch Input Parameters: 1702c71b3e2SJacob Faibussowitsch + cond - The boolean condition 1712c71b3e2SJacob Faibussowitsch . comm - The communicator on which the check can be collective on 1722c71b3e2SJacob Faibussowitsch . ierr - A nonzero error code, see include/petscerror.h for the complete list 173e55b9c91SBarry Smith - message - Error message in the `printf()` format 1742c71b3e2SJacob Faibussowitsch 175667f096bSBarry Smith Level: beginner 176667f096bSBarry Smith 1772c71b3e2SJacob Faibussowitsch Notes: 1782c71b3e2SJacob Faibussowitsch Enabled in both optimized and debug builds. 1792c71b3e2SJacob Faibussowitsch 180a2454668SBarry Smith As a general rule, `PetscCheck()` is used to check "usage error" (for example, passing an incorrect value as a function argument), 181a2454668SBarry Smith `PetscAssert()` is used to "check for bugs in PETSc" (for example, is a value in a PETSc data structure nonsensical). 182a2454668SBarry Smith However, for functions that are called in a "hot spot", for example, thousands of times in a loop, `PetscAssert()` should be used instead 183a2454668SBarry Smith of `PetscCheck()` since the former is compiled out in PETSc's optimization code. 184a2454668SBarry Smith 18587497f52SBarry Smith Calls `SETERRQ()` if the assertion fails, so can only be called from functions returning a 18687497f52SBarry Smith `PetscErrorCode` (or equivalent type after conversion). 1872c71b3e2SJacob Faibussowitsch 188af27ebaaSBarry Smith .seealso: `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheckAbort()`, `PetscErrorCode` 1899566063dSJacob Faibussowitsch M*/ 1909371c9d4SSatish Balay #define PetscCheck(cond, comm, ierr, ...) \ 1919371c9d4SSatish Balay do { \ 1929371c9d4SSatish Balay if (PetscUnlikely(!(cond))) SETERRQ(comm, ierr, __VA_ARGS__); \ 1939371c9d4SSatish Balay } while (0) 1942c71b3e2SJacob Faibussowitsch 1952c71b3e2SJacob Faibussowitsch /*MC 1964be741a6SBarry Smith PetscCheckAbort - Check that a particular condition is true, otherwise prints error and aborts 1974be741a6SBarry Smith 1984be741a6SBarry Smith Synopsis: 1994be741a6SBarry Smith #include <petscerror.h> 2004be741a6SBarry Smith void PetscCheckAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 2014be741a6SBarry Smith 2027cdbe19fSJose E. Roman Collective; No Fortran Support 2034be741a6SBarry Smith 2044be741a6SBarry Smith Input Parameters: 2054be741a6SBarry Smith + cond - The boolean condition 2064be741a6SBarry Smith . comm - The communicator on which the check can be collective on 2074be741a6SBarry Smith . ierr - A nonzero error code, see include/petscerror.h for the complete list 208e55b9c91SBarry Smith - message - Error message in the `printf()` format 2094be741a6SBarry Smith 210667f096bSBarry Smith Level: developer 211667f096bSBarry Smith 2124be741a6SBarry Smith Notes: 2134be741a6SBarry Smith Enabled in both optimized and debug builds. 2144be741a6SBarry Smith 21587497f52SBarry Smith Calls `SETERRABORT()` if the assertion fails, can be called from a function that does not return an 21687497f52SBarry Smith error code, such as a C++ constructor. usually `PetscCheck()` should be used. 2174be741a6SBarry Smith 218af27ebaaSBarry Smith .seealso: `PetscAssertAbort()`, `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheck()`, `SETERRABORT()`, `PetscErrorCode` 2194be741a6SBarry Smith M*/ 2209371c9d4SSatish Balay #define PetscCheckAbort(cond, comm, ierr, ...) \ 2210e6b6b59SJacob Faibussowitsch do { \ 2220e6b6b59SJacob Faibussowitsch if (PetscUnlikely(!(cond))) SETERRABORT(comm, ierr, __VA_ARGS__); \ 223c7481402SJacob Faibussowitsch } while (0) 2244be741a6SBarry Smith 2254be741a6SBarry Smith /*MC 2269ace16cdSJacob Faibussowitsch PetscAssert - Assert that a particular condition is true 2279ace16cdSJacob Faibussowitsch 2289ace16cdSJacob Faibussowitsch Synopsis: 2299ace16cdSJacob Faibussowitsch #include <petscerror.h> 2309ace16cdSJacob Faibussowitsch void PetscAssert(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 2319ace16cdSJacob Faibussowitsch 2327cdbe19fSJose E. Roman Collective; No Fortran Support 2339ace16cdSJacob Faibussowitsch 2349ace16cdSJacob Faibussowitsch Input Parameters: 2359ace16cdSJacob Faibussowitsch + cond - The boolean condition 2369ace16cdSJacob Faibussowitsch . comm - The communicator on which the check can be collective on 2379ace16cdSJacob Faibussowitsch . ierr - A nonzero error code, see include/petscerror.h for the complete list 238af27ebaaSBarry Smith - message - Error message in `printf()` format 2399ace16cdSJacob Faibussowitsch 240667f096bSBarry Smith Level: beginner 241667f096bSBarry Smith 2429ace16cdSJacob Faibussowitsch Notes: 243c7481402SJacob Faibussowitsch Equivalent to `PetscCheck()` if debugging is enabled, and `PetscAssume(cond)` otherwise. 2442c71b3e2SJacob Faibussowitsch 24587497f52SBarry Smith See `PetscCheck()` for usage and behaviour. 24687497f52SBarry Smith 24787497f52SBarry Smith This is needed instead of simply using `assert()` because this correctly handles the collective nature of errors under MPI 2489ace16cdSJacob Faibussowitsch 249af27ebaaSBarry Smith .seealso: `PetscCheck()`, `SETERRQ()`, `PetscError()`, `PetscAssertAbort()`, `PetscErrorCode` 2509566063dSJacob Faibussowitsch M*/ 251c7481402SJacob Faibussowitsch #if PetscDefined(USE_DEBUG) 252c7481402SJacob Faibussowitsch #define PetscAssert(cond, comm, ierr, ...) PetscCheck(cond, comm, ierr, __VA_ARGS__) 253c7481402SJacob Faibussowitsch #else 254c7481402SJacob Faibussowitsch #define PetscAssert(cond, ...) PetscAssume(cond) 255c7481402SJacob Faibussowitsch #endif 2569ace16cdSJacob Faibussowitsch 2579ace16cdSJacob Faibussowitsch /*MC 2580e6b6b59SJacob Faibussowitsch PetscAssertAbort - Assert that a particular condition is true, otherwise prints error and aborts 2590e6b6b59SJacob Faibussowitsch 2600e6b6b59SJacob Faibussowitsch Synopsis: 2610e6b6b59SJacob Faibussowitsch #include <petscerror.h> 2620e6b6b59SJacob Faibussowitsch void PetscAssertAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 2630e6b6b59SJacob Faibussowitsch 2647cdbe19fSJose E. Roman Collective; No Fortran Support 2650e6b6b59SJacob Faibussowitsch 2660e6b6b59SJacob Faibussowitsch Input Parameters: 2670e6b6b59SJacob Faibussowitsch + cond - The boolean condition 2680e6b6b59SJacob Faibussowitsch . comm - The communicator on which the check can be collective on 2690e6b6b59SJacob Faibussowitsch . ierr - A nonzero error code, see include/petscerror.h for the complete list 270e55b9c91SBarry Smith - message - Error message in the `printf()` format 2710e6b6b59SJacob Faibussowitsch 272667f096bSBarry Smith Level: beginner 273667f096bSBarry Smith 27495bd0b28SBarry Smith Note: 2750e6b6b59SJacob Faibussowitsch Enabled only in debug builds. See `PetscCheckAbort()` for usage. 2760e6b6b59SJacob Faibussowitsch 2770e6b6b59SJacob Faibussowitsch .seealso: `PetscCheckAbort()`, `PetscAssert()`, `PetscCheck()`, `SETERRABORT()`, `PetscError()` 2780e6b6b59SJacob Faibussowitsch M*/ 2793ba16761SJacob Faibussowitsch #if PetscDefined(USE_DEBUG) 2803ba16761SJacob Faibussowitsch #define PetscAssertAbort(cond, comm, ierr, ...) PetscCheckAbort(cond, comm, ierr, __VA_ARGS__) 2813ba16761SJacob Faibussowitsch #else 2823ba16761SJacob Faibussowitsch #define PetscAssertAbort(cond, comm, ierr, ...) PetscAssume(cond) 2833ba16761SJacob Faibussowitsch #endif 2840e6b6b59SJacob Faibussowitsch 2850e6b6b59SJacob Faibussowitsch /*MC 2863ba16761SJacob Faibussowitsch PetscCall - Calls a PETSc function and then checks the resulting error code, if it is 2873ba16761SJacob Faibussowitsch non-zero it calls the error handler and returns from the current function with the error 2883ba16761SJacob Faibussowitsch code. 2899566063dSJacob Faibussowitsch 2909566063dSJacob Faibussowitsch Synopsis: 2919566063dSJacob Faibussowitsch #include <petscerror.h> 29249c86fc7SBarry Smith void PetscCall(PetscFunction(args)) 2939566063dSJacob Faibussowitsch 2949566063dSJacob Faibussowitsch Not Collective 2959566063dSJacob Faibussowitsch 2969566063dSJacob Faibussowitsch Input Parameter: 29749c86fc7SBarry Smith . PetscFunction - any PETSc function that returns an error code 2989566063dSJacob Faibussowitsch 299667f096bSBarry Smith Level: beginner 300667f096bSBarry Smith 3019566063dSJacob Faibussowitsch Notes: 3029566063dSJacob Faibussowitsch Once the error handler is called the calling function is then returned from with the given 30387497f52SBarry Smith error code. Experienced users can set the error handler with `PetscPushErrorHandler()`. 3049566063dSJacob Faibussowitsch 30587497f52SBarry Smith `PetscCall()` cannot be used in functions returning a datatype not convertible to 3068af9f190SBarry Smith `PetscErrorCode`. For example, `PetscCall()` may not be used in functions returning `void`, use 3073ba16761SJacob Faibussowitsch `PetscCallAbort()` or `PetscCallVoid()` in this case. 3089566063dSJacob Faibussowitsch 3099566063dSJacob Faibussowitsch Example Usage: 3109566063dSJacob Faibussowitsch .vb 3119566063dSJacob Faibussowitsch PetscCall(PetscInitiailize(...)); // OK to call even when PETSc is not yet initialized! 3129566063dSJacob Faibussowitsch 3139566063dSJacob Faibussowitsch struct my_struct 3149566063dSJacob Faibussowitsch { 3159566063dSJacob Faibussowitsch void *data; 3169566063dSJacob Faibussowitsch } my_complex_type; 3179566063dSJacob Faibussowitsch 3189566063dSJacob Faibussowitsch struct my_struct bar(void) 3199566063dSJacob Faibussowitsch { 3206aad120cSJose E. Roman PetscCall(foo(15)); // ERROR PetscErrorCode not convertible to struct my_struct! 3219566063dSJacob Faibussowitsch } 3229566063dSJacob Faibussowitsch 3236aad120cSJose E. Roman PetscCall(bar()) // ERROR input not convertible to PetscErrorCode 3249566063dSJacob Faibussowitsch .ve 3259566063dSJacob Faibussowitsch 32687497f52SBarry Smith It is also possible to call this directly on a `PetscErrorCode` variable 32749c86fc7SBarry Smith .vb 32849c86fc7SBarry Smith PetscCall(ierr); // check if ierr is nonzero 32949c86fc7SBarry Smith .ve 33049c86fc7SBarry Smith 331792fecdfSBarry Smith Should not be used to call callback functions provided by users, `PetscCallBack()` should be used in that situation. 332ef1023bdSBarry Smith 3336a8be23eSBarry Smith `PetscUseTypeMethod()` or `PetscTryTypeMethod()` should be used when calling functions pointers contained in a PETSc object's `ops` array 3346a8be23eSBarry Smith 33549c86fc7SBarry Smith Fortran Notes: 33698d50953SPierre Jolivet The Fortran function in which this is used must declare a `PetscErrorCode` variable necessarily named `ierr`, and `ierr` must be 33787497f52SBarry Smith the final argument to the PETSc function being called. 33849c86fc7SBarry Smith 33998d50953SPierre Jolivet In the main program and in Fortran subroutines that do not have `ierr` as the final return parameter, one 34087497f52SBarry Smith should use `PetscCallA()` 34149c86fc7SBarry Smith 34249c86fc7SBarry Smith Example Fortran Usage: 34349c86fc7SBarry Smith .vb 34449c86fc7SBarry Smith PetscErrorCode ierr 34549c86fc7SBarry Smith Vec v 34649c86fc7SBarry Smith 34749c86fc7SBarry Smith ... 34849c86fc7SBarry Smith PetscCall(VecShift(v, 1.0, ierr)) 34949c86fc7SBarry Smith PetscCallA(VecShift(v, 1.0, ierr)) 35049c86fc7SBarry Smith .ve 35149c86fc7SBarry Smith 352667f096bSBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`, 353667f096bSBarry Smith `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, 354648c30bcSBarry Smith `CHKERRMPI()`, `PetscCallBack()`, `PetscCallAbort()`, `PetscCallVoid()`, `PetscCallNull()` 355648c30bcSBarry Smith M*/ 356648c30bcSBarry Smith 357648c30bcSBarry Smith /*MC 358648c30bcSBarry Smith PetscCallNull - Calls a PETSc function and then checks the resulting error code, if it is 359648c30bcSBarry Smith non-zero it calls the error handler and returns a `NULL` 360648c30bcSBarry Smith 361648c30bcSBarry Smith Synopsis: 362648c30bcSBarry Smith #include <petscerror.h> 363648c30bcSBarry Smith void PetscCallNull(PetscFunction(args)) 364648c30bcSBarry Smith 365648c30bcSBarry Smith Not Collective; No Fortran Support 366648c30bcSBarry Smith 367648c30bcSBarry Smith Input Parameter: 368648c30bcSBarry Smith . PetscFunction - any PETSc function that returns something that can be returned as a `NULL` 369648c30bcSBarry Smith 370648c30bcSBarry Smith Level: developer 371648c30bcSBarry Smith 372648c30bcSBarry Smith .seealso: `PetscCall()`, `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`, 373648c30bcSBarry Smith `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, 374648c30bcSBarry Smith `CHKERRMPI()`, `PetscCallBack()`, `PetscCallAbort()`, `PetscCallVoid()`, `PetscCall()` 3759566063dSJacob Faibussowitsch M*/ 376ef1023bdSBarry Smith 377ef1023bdSBarry Smith /*MC 37898d50953SPierre Jolivet PetscCallA - Fortran-only macro that should be used in the main program and subroutines that do not have `ierr` as the final return parameter, to call PETSc functions instead of using 37998d50953SPierre Jolivet `PetscCall()` which should be used in other Fortran subroutines 380989c49a2SBarry Smith 381989c49a2SBarry Smith Synopsis: 382989c49a2SBarry Smith #include <petscsys.h> 383989c49a2SBarry Smith PetscErrorCode PetscCallA(PetscFunction(arguments, ierr)) 384989c49a2SBarry Smith 385989c49a2SBarry Smith Collective 386989c49a2SBarry Smith 387989c49a2SBarry Smith Input Parameter: 388989c49a2SBarry Smith . PetscFunction(arguments,ierr) - the call to the function 389989c49a2SBarry Smith 390989c49a2SBarry Smith Level: beginner 391989c49a2SBarry Smith 392989c49a2SBarry Smith Notes: 393989c49a2SBarry Smith This should only be used with Fortran. With C/C++, use `PetscCall()` always. 394989c49a2SBarry Smith 39598d50953SPierre Jolivet The Fortran function in which this is used must declare a `PetscErrorCode` variable necessarily named `ierr` 396989c49a2SBarry Smith Use `SETERRA()` to set an error in a Fortran main program and `SETERRQ()` in Fortran subroutines 397989c49a2SBarry Smith 398989c49a2SBarry Smith .seealso: `SETERRQ()`, `SETERRA()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()` 399989c49a2SBarry Smith M*/ 400989c49a2SBarry Smith 401989c49a2SBarry Smith /*MC 402792fecdfSBarry 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 403ef1023bdSBarry Smith handler and returns from the current function with the error code. 404ef1023bdSBarry Smith 405ef1023bdSBarry Smith Synopsis: 406ef1023bdSBarry Smith #include <petscerror.h> 407792fecdfSBarry Smith void PetscCallBack(const char *functionname, PetscFunction(args)) 408ef1023bdSBarry Smith 4097cdbe19fSJose E. Roman Not Collective; No Fortran Support 410ef1023bdSBarry Smith 411ef1023bdSBarry Smith Input Parameters: 412ef1023bdSBarry Smith + functionname - the name of the function being called, this can be a string with spaces that describes the meaning of the callback 41387497f52SBarry Smith - PetscFunction - user provided callback function that returns an error code 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 422667f096bSBarry Smith Notes: 4237e6f8dd6SBarry Smith `PetscUseTypeMethod()` and ` PetscTryTypeMethod()` are the preferred API for this functionality. But when the callback functions are associated with a 4247e6f8dd6SBarry Smith `DMSNES` or `DMTS` this API must be used. 4257e6f8dd6SBarry Smith 426667f096bSBarry Smith Once the error handler is called the calling function is then returned from with the given 427667f096bSBarry Smith error code. Experienced users can set the error handler with `PetscPushErrorHandler()`. 428667f096bSBarry Smith 429667f096bSBarry Smith `PetscCallBack()` should only be called in PETSc when a call is being made to a user provided call-back routine. 430667f096bSBarry Smith 4317e6f8dd6SBarry Smith Developer Note: 4327e6f8dd6SBarry Smith It would be good to provide a new API for when the callbacks are associated with `DMSNES` or `DMTS` so this routine could be used less 4337e6f8dd6SBarry Smith 43487497f52SBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscCall()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()` 4357e6f8dd6SBarry Smith `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, `CHKERRMPI()`, `PetscCall()`, `PetscUseTypeMethod()`, `PetscTryTypeMethod()` 436ef1023bdSBarry Smith M*/ 437ef1023bdSBarry Smith 4383ba16761SJacob Faibussowitsch /*MC 4398af9f190SBarry Smith PetscCallVoid - Like `PetscCall()` but for use in functions that return `void` 4403ba16761SJacob Faibussowitsch 4413ba16761SJacob Faibussowitsch Synopsis: 4423ba16761SJacob Faibussowitsch #include <petscerror.h> 4438af9f190SBarry Smith void PetscCallVoid(PetscFunction(args)) 4443ba16761SJacob Faibussowitsch 4457cdbe19fSJose E. Roman Not Collective; No Fortran Support 4463ba16761SJacob Faibussowitsch 4473ba16761SJacob Faibussowitsch Input Parameter: 4483ba16761SJacob Faibussowitsch . PetscFunction - any PETSc function that returns an error code 4493ba16761SJacob Faibussowitsch 4503ba16761SJacob Faibussowitsch Example Usage: 4513ba16761SJacob Faibussowitsch .vb 4523ba16761SJacob Faibussowitsch void foo() 4533ba16761SJacob Faibussowitsch { 4543ba16761SJacob Faibussowitsch KSP ksp; 4553ba16761SJacob Faibussowitsch 4563ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 4573ba16761SJacob Faibussowitsch // OK, properly handles PETSc error codes 4583ba16761SJacob Faibussowitsch PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp)); 4598af9f190SBarry Smith PetscFunctionReturnVoid(); 4603ba16761SJacob Faibussowitsch } 4613ba16761SJacob Faibussowitsch 4623ba16761SJacob Faibussowitsch PetscErrorCode bar() 4633ba16761SJacob Faibussowitsch { 4643ba16761SJacob Faibussowitsch KSP ksp; 4653ba16761SJacob Faibussowitsch 4663ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 4673ba16761SJacob Faibussowitsch // ERROR, Non-void function 'bar' should return a value 4683ba16761SJacob Faibussowitsch PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp)); 4693ba16761SJacob Faibussowitsch // OK, returning PetscErrorCode 4703ba16761SJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp)); 4713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4723ba16761SJacob Faibussowitsch } 473667f096bSBarry Smith .ve 4743ba16761SJacob Faibussowitsch 4753ba16761SJacob Faibussowitsch Level: beginner 4763ba16761SJacob Faibussowitsch 477667f096bSBarry Smith Notes: 478667f096bSBarry Smith Has identical usage to `PetscCall()`, except that it returns `void` on error instead of a 479667f096bSBarry Smith `PetscErrorCode`. See `PetscCall()` for more detailed discussion. 480667f096bSBarry Smith 481667f096bSBarry Smith Note that users should prefer `PetscCallAbort()` to this routine. While this routine does 482667f096bSBarry Smith "handle" errors by returning from the enclosing function, it effectively gobbles the 483667f096bSBarry Smith error. Since the enclosing function itself returns `void`, its callers have no way of knowing 484667f096bSBarry Smith that the routine returned early due to an error. `PetscCallAbort()` at least ensures that the 485667f096bSBarry Smith program crashes gracefully. 486667f096bSBarry Smith 487648c30bcSBarry Smith .seealso: `PetscCall()`, `PetscErrorCode`, `PetscCallAbort()`, `PetscCallNull()` 4883ba16761SJacob Faibussowitsch M*/ 4899566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 4909566063dSJacob Faibussowitsch void PetscCall(PetscErrorCode); 491792fecdfSBarry Smith void PetscCallBack(const char *, PetscErrorCode); 4929566063dSJacob Faibussowitsch void PetscCallVoid(PetscErrorCode); 493648c30bcSBarry Smith void PetscCallNull(PetscErrorCode); 4949566063dSJacob Faibussowitsch #else 4959371c9d4SSatish Balay #define PetscCall(...) \ 4969371c9d4SSatish Balay do { \ 4973ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_q_; \ 49837154d10SBarry Smith PetscStackUpdateLine; \ 4993ba16761SJacob Faibussowitsch ierr_petsc_call_q_ = __VA_ARGS__; \ 5003ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_q_ != PETSC_SUCCESS)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_q_, PETSC_ERROR_REPEAT, " "); \ 5019566063dSJacob Faibussowitsch } while (0) 502648c30bcSBarry Smith #define PetscCallNull(...) \ 503648c30bcSBarry Smith do { \ 504648c30bcSBarry Smith PetscErrorCode ierr_petsc_call_q_; \ 505648c30bcSBarry Smith PetscStackUpdateLine; \ 506648c30bcSBarry Smith ierr_petsc_call_q_ = __VA_ARGS__; \ 507648c30bcSBarry Smith if (PetscUnlikely(ierr_petsc_call_q_ != PETSC_SUCCESS)) { \ 508648c30bcSBarry Smith (void)PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); \ 509648c30bcSBarry Smith PetscFunctionReturn(NULL); \ 510648c30bcSBarry Smith } \ 511648c30bcSBarry Smith } while (0) 5129371c9d4SSatish Balay #define PetscCallBack(function, ...) \ 5139371c9d4SSatish Balay do { \ 5143ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_q_; \ 515e33ced7fSLisandro Dalcin PetscStackUpdateLine; \ 516e33ced7fSLisandro Dalcin PetscStackPushExternal(function); \ 5173ba16761SJacob Faibussowitsch ierr_petsc_call_q_ = __VA_ARGS__; \ 518ef1023bdSBarry Smith PetscStackPop; \ 5193ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_q_ != PETSC_SUCCESS)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_q_, PETSC_ERROR_REPEAT, " "); \ 520ef1023bdSBarry Smith } while (0) 5219371c9d4SSatish Balay #define PetscCallVoid(...) \ 5229371c9d4SSatish Balay do { \ 5233ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_void_; \ 524e33ced7fSLisandro Dalcin PetscStackUpdateLine; \ 5253ba16761SJacob Faibussowitsch ierr_petsc_call_void_ = __VA_ARGS__; \ 5263ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_void_ != PETSC_SUCCESS)) { \ 527dd460d27SBarry Smith (void)PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_void_, PETSC_ERROR_REPEAT, " "); \ 5289371c9d4SSatish Balay return; \ 5299371c9d4SSatish Balay } \ 5309566063dSJacob Faibussowitsch } while (0) 5319566063dSJacob Faibussowitsch #endif 5329566063dSJacob Faibussowitsch 5339566063dSJacob Faibussowitsch /*MC 5349566063dSJacob Faibussowitsch CHKERRQ - Checks error code returned from PETSc function 53530de9b25SBarry Smith 53630de9b25SBarry Smith Synopsis: 537aaa7dc30SBarry Smith #include <petscsys.h> 5389566063dSJacob Faibussowitsch void CHKERRQ(PetscErrorCode ierr) 5399566063dSJacob Faibussowitsch 5409566063dSJacob Faibussowitsch Not Collective 5419566063dSJacob Faibussowitsch 5422fe279fdSBarry Smith Input Parameter: 5439566063dSJacob Faibussowitsch . ierr - nonzero error code 5449566063dSJacob Faibussowitsch 5459566063dSJacob Faibussowitsch Level: deprecated 5469566063dSJacob Faibussowitsch 547667f096bSBarry Smith Note: 548667f096bSBarry Smith Deprecated in favor of `PetscCall()`. This routine behaves identically to it. 549667f096bSBarry Smith 550db781477SPatrick Sanan .seealso: `PetscCall()` 5519566063dSJacob Faibussowitsch M*/ 5529566063dSJacob Faibussowitsch #define CHKERRQ(...) PetscCall(__VA_ARGS__) 5539566063dSJacob Faibussowitsch #define CHKERRV(...) PetscCallVoid(__VA_ARGS__) 5549566063dSJacob Faibussowitsch 555a1fd7ae3SBarry Smith PETSC_EXTERN void PetscMPIErrorString(PetscMPIInt, size_t, char *); 556db9cea48SBarry Smith 5579566063dSJacob Faibussowitsch /*MC 5589566063dSJacob Faibussowitsch PetscCallMPI - Checks error code returned from MPI calls, if non-zero it calls the error 559648c30bcSBarry Smith handler and then returns a `PetscErrorCode` 5609566063dSJacob Faibussowitsch 5619566063dSJacob Faibussowitsch Synopsis: 5629566063dSJacob Faibussowitsch #include <petscerror.h> 56349c86fc7SBarry Smith void PetscCallMPI(MPI_Function(args)) 56430de9b25SBarry Smith 565eca87e8dSBarry Smith Not Collective 56630de9b25SBarry Smith 5672fe279fdSBarry Smith Input Parameter: 56849c86fc7SBarry Smith . MPI_Function - an MPI function that returns an MPI error code 56930de9b25SBarry Smith 570667f096bSBarry Smith Level: beginner 571667f096bSBarry Smith 5729566063dSJacob Faibussowitsch Notes: 57387497f52SBarry Smith Always returns the error code `PETSC_ERR_MPI`; the MPI error code and string are embedded in 5749566063dSJacob Faibussowitsch the string error message. Do not use this to call any other routines (for example PETSc 5753ba16761SJacob Faibussowitsch routines), it should only be used for direct MPI calls. The user may configure PETSc with the 5763ba16761SJacob Faibussowitsch `--with-strict-petscerrorcode` option to check this at compile-time, otherwise they must 5779566063dSJacob Faibussowitsch check this themselves. 5789566063dSJacob Faibussowitsch 579aaa8cc7dSPierre Jolivet This routine can only be used in functions returning `PetscErrorCode` themselves. If the 5803ba16761SJacob Faibussowitsch calling function returns a different type, use `PetscCallMPIAbort()` instead. 5813ba16761SJacob Faibussowitsch 5829566063dSJacob Faibussowitsch Example Usage: 5839566063dSJacob Faibussowitsch .vb 5849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(...)); // OK, calling MPI function 5859566063dSJacob Faibussowitsch 5869566063dSJacob Faibussowitsch PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead! 5879566063dSJacob Faibussowitsch .ve 5889566063dSJacob Faibussowitsch 58949c86fc7SBarry Smith Fortran Notes: 59087497f52SBarry Smith The Fortran function from which this is used must declare a variable `PetscErrorCode` ierr and ierr must be 59149c86fc7SBarry Smith the final argument to the MPI function being called. 59249c86fc7SBarry Smith 59349c86fc7SBarry Smith In the main program and in Fortran subroutines that do not have ierr as the final return parameter one 59487497f52SBarry Smith should use `PetscCallMPIA()` 59549c86fc7SBarry Smith 59649c86fc7SBarry Smith Fortran Usage: 59749c86fc7SBarry Smith .vb 59849c86fc7SBarry Smith PetscErrorCode ierr or integer ierr 59949c86fc7SBarry Smith ... 60049c86fc7SBarry Smith PetscCallMPI(MPI_Comm_size(...,ierr)) 60149c86fc7SBarry Smith PetscCallMPIA(MPI_Comm_size(...,ierr)) ! Will abort after calling error handler 60249c86fc7SBarry Smith 60349c86fc7SBarry Smith PetscCallMPI(MPI_Comm_size(...,eflag)) ! ERROR, final argument must be ierr 60449c86fc7SBarry Smith .ve 60549c86fc7SBarry Smith 606db781477SPatrick Sanan .seealso: `SETERRMPI()`, `PetscCall()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`, 6073ba16761SJacob Faibussowitsch `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 608648c30bcSBarry Smith `PetscError()`, `CHKMEMQ`, `PetscCallMPINull()` 609648c30bcSBarry Smith M*/ 610648c30bcSBarry Smith 611648c30bcSBarry Smith /*MC 612648c30bcSBarry Smith PetscCallMPINull - Checks error code returned from MPI calls, if non-zero it calls the error 613648c30bcSBarry Smith handler and then returns a `NULL` 614648c30bcSBarry Smith 615648c30bcSBarry Smith Synopsis: 616648c30bcSBarry Smith #include <petscerror.h> 617648c30bcSBarry Smith void PetscCallMPINull(MPI_Function(args)) 618648c30bcSBarry Smith 619648c30bcSBarry Smith Not Collective; No Fortran Support 620648c30bcSBarry Smith 621648c30bcSBarry Smith Input Parameter: 622648c30bcSBarry Smith . MPI_Function - an MPI function that returns an MPI error code 623648c30bcSBarry Smith 624648c30bcSBarry Smith Level: beginner 625648c30bcSBarry Smith 626648c30bcSBarry Smith Notes: 627648c30bcSBarry Smith Always passes the error code `PETSC_ERR_MPI` to the error handler `PetscError()`; the MPI error code and string are embedded in 628648c30bcSBarry Smith the string error message. Do not use this to call any other routines (for example PETSc 629648c30bcSBarry Smith routines), it should only be used for direct MPI calls. 630648c30bcSBarry Smith 631648c30bcSBarry Smith This routine can only be used in functions returning anything that can be returned as a `NULL` themselves. If the 632648c30bcSBarry Smith calling function returns a different type, use `PetscCallMPIAbort()` instead. 633648c30bcSBarry Smith 634648c30bcSBarry Smith Example Usage: 635648c30bcSBarry Smith .vb 636648c30bcSBarry Smith PetscCallMPINull(MPI_Comm_size(...)); // OK, calling MPI function 637648c30bcSBarry Smith 638648c30bcSBarry Smith PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead! 639648c30bcSBarry Smith .ve 640648c30bcSBarry Smith 641648c30bcSBarry Smith .seealso: `SETERRMPI()`, `PetscCall()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`, 642648c30bcSBarry Smith `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 643648c30bcSBarry Smith `PetscError()`, `CHKMEMQ`, `PetscCallMPI()` 6443ba16761SJacob Faibussowitsch M*/ 6453ba16761SJacob Faibussowitsch 6463ba16761SJacob Faibussowitsch /*MC 6473ba16761SJacob Faibussowitsch PetscCallMPIAbort - Like `PetscCallMPI()` but calls `MPI_Abort()` on error 6483ba16761SJacob Faibussowitsch 6493ba16761SJacob Faibussowitsch Synopsis: 6503ba16761SJacob Faibussowitsch #include <petscerror.h> 6513ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm comm, MPI_Function(args)) 6523ba16761SJacob Faibussowitsch 6533ba16761SJacob Faibussowitsch Not Collective 6543ba16761SJacob Faibussowitsch 6553ba16761SJacob Faibussowitsch Input Parameters: 6563ba16761SJacob Faibussowitsch + comm - the MPI communicator to abort on 6573ba16761SJacob Faibussowitsch - MPI_Function - an MPI function that returns an MPI error code 6583ba16761SJacob Faibussowitsch 659667f096bSBarry Smith Level: beginner 660667f096bSBarry Smith 6613ba16761SJacob Faibussowitsch Notes: 6623ba16761SJacob Faibussowitsch Usage is identical to `PetscCallMPI()`. See `PetscCallMPI()` for detailed discussion. 6633ba16761SJacob Faibussowitsch 6643ba16761SJacob Faibussowitsch This routine may be used in functions returning `void` or other non-`PetscErrorCode` types. 6653ba16761SJacob Faibussowitsch 666989c49a2SBarry Smith Fortran Note: 667989c49a2SBarry Smith In Fortran this is called `PetscCallMPIA()` and is intended to be used in the main program while `PetscCallMPI()` is 668989c49a2SBarry Smith used in Fortran subroutines. 669989c49a2SBarry Smith 670989c49a2SBarry Smith Developer Note: 671989c49a2SBarry Smith This should have the same name in Fortran. 672989c49a2SBarry Smith 6733ba16761SJacob Faibussowitsch .seealso: `PetscCallMPI()`, `PetscCallAbort()`, `SETERRABORT()` 67430de9b25SBarry Smith M*/ 6753fcd9f07SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 67663a3b9bcSJacob Faibussowitsch void PetscCallMPI(PetscMPIInt); 6773ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm, PetscMPIInt); 678648c30bcSBarry Smith void PetscCallMPINull(PetscMPIInt); 679064a246eSJacob Faibussowitsch #else 6803ba16761SJacob Faibussowitsch #define PetscCallMPI_Private(__PETSC_STACK_POP_FUNC__, __SETERR_FUNC__, __COMM__, ...) \ 6819371c9d4SSatish Balay do { \ 6823ba16761SJacob Faibussowitsch PetscMPIInt ierr_petsc_call_mpi_; \ 683ef1023bdSBarry Smith PetscStackUpdateLine; \ 684792fecdfSBarry Smith PetscStackPushExternal("MPI function"); \ 685d71ae5a4SJacob Faibussowitsch { \ 6863ba16761SJacob Faibussowitsch ierr_petsc_call_mpi_ = __VA_ARGS__; \ 687d71ae5a4SJacob Faibussowitsch } \ 6883ba16761SJacob Faibussowitsch __PETSC_STACK_POP_FUNC__; \ 6893ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_mpi_ != MPI_SUCCESS)) { \ 6903ba16761SJacob Faibussowitsch char petsc_mpi_7_errorstring[2 * MPI_MAX_ERROR_STRING]; \ 691a1fd7ae3SBarry Smith PetscMPIErrorString(ierr_petsc_call_mpi_, 2 * MPI_MAX_ERROR_STRING, (char *)petsc_mpi_7_errorstring); \ 692a1fd7ae3SBarry Smith __SETERR_FUNC__(__COMM__, PETSC_ERR_MPI, "MPI error %d %s", ierr_petsc_call_mpi_, petsc_mpi_7_errorstring); \ 6933fcd9f07SJacob Faibussowitsch } \ 6943fcd9f07SJacob Faibussowitsch } while (0) 6953ba16761SJacob Faibussowitsch 6963ba16761SJacob Faibussowitsch #define PetscCallMPI(...) PetscCallMPI_Private(PetscStackPop, SETERRQ, PETSC_COMM_SELF, __VA_ARGS__) 6973ba16761SJacob Faibussowitsch #define PetscCallMPIAbort(comm, ...) PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRABORT, comm, __VA_ARGS__) 698648c30bcSBarry Smith #define PetscCallMPINull(...) PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRQNULL, PETSC_COMM_SELF, __VA_ARGS__) 6999566063dSJacob Faibussowitsch #endif 700064a246eSJacob Faibussowitsch 7017037b0edSPatrick Sanan /*MC 7029566063dSJacob Faibussowitsch CHKERRMPI - Checks error code returned from MPI calls, if non-zero it calls the error 7039566063dSJacob Faibussowitsch handler and then returns 7049566063dSJacob Faibussowitsch 7059566063dSJacob Faibussowitsch Synopsis: 7069566063dSJacob Faibussowitsch #include <petscerror.h> 7079566063dSJacob Faibussowitsch void CHKERRMPI(PetscErrorCode ierr) 7089566063dSJacob Faibussowitsch 7099566063dSJacob Faibussowitsch Not Collective 7109566063dSJacob Faibussowitsch 7119566063dSJacob Faibussowitsch Input Parameter: 7129566063dSJacob Faibussowitsch . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 7139566063dSJacob Faibussowitsch 7149566063dSJacob Faibussowitsch Level: deprecated 7159566063dSJacob Faibussowitsch 716667f096bSBarry Smith Note: 717667f096bSBarry Smith Deprecated in favor of `PetscCallMPI()`. This routine behaves identically to it. 718667f096bSBarry Smith 719db781477SPatrick Sanan .seealso: `PetscCallMPI()` 7209566063dSJacob Faibussowitsch M*/ 7219566063dSJacob Faibussowitsch #define CHKERRMPI(...) PetscCallMPI(__VA_ARGS__) 7229566063dSJacob Faibussowitsch 7239566063dSJacob Faibussowitsch /*MC 724989c49a2SBarry Smith PetscCallAbort - Checks error code returned from PETSc function, if non-zero it aborts immediately by calling `MPI_Abort()` 7259566063dSJacob Faibussowitsch 7269566063dSJacob Faibussowitsch Synopsis: 7279566063dSJacob Faibussowitsch #include <petscerror.h> 7289566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm comm, PetscErrorCode ierr) 7299566063dSJacob Faibussowitsch 730c3339decSBarry Smith Collective 7319566063dSJacob Faibussowitsch 7329566063dSJacob Faibussowitsch Input Parameters: 7339566063dSJacob Faibussowitsch + comm - the MPI communicator on which to abort 7349566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 7359566063dSJacob Faibussowitsch 736667f096bSBarry Smith Level: intermediate 737667f096bSBarry Smith 7389566063dSJacob Faibussowitsch Notes: 73987497f52SBarry Smith This macro has identical type and usage semantics to `PetscCall()` with the important caveat 7409566063dSJacob Faibussowitsch that this macro does not return. Instead, if ierr is nonzero it calls the PETSc error handler 74187497f52SBarry Smith and then immediately calls `MPI_Abort()`. It can therefore be used anywhere. 7429566063dSJacob Faibussowitsch 74387497f52SBarry Smith As per `MPI_Abort()` semantics the communicator passed must be valid, although there is currently 74487497f52SBarry Smith no attempt made at handling any potential errors from `MPI_Abort()`. Note that while 74587497f52SBarry Smith `MPI_Abort()` is required to terminate only those processes which reside on comm, it is often 74687497f52SBarry Smith the case that `MPI_Abort()` terminates *all* processes. 7479566063dSJacob Faibussowitsch 7489566063dSJacob Faibussowitsch Example Usage: 7499566063dSJacob Faibussowitsch .vb 7509566063dSJacob Faibussowitsch PetscErrorCode boom(void) { return PETSC_ERR_MEM; } 7519566063dSJacob Faibussowitsch 7529566063dSJacob Faibussowitsch void foo(void) 7539566063dSJacob Faibussowitsch { 7549566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type 7559566063dSJacob Faibussowitsch } 7569566063dSJacob Faibussowitsch 7579566063dSJacob Faibussowitsch double bar(void) 7589566063dSJacob Faibussowitsch { 7599566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type 7609566063dSJacob Faibussowitsch } 7619566063dSJacob Faibussowitsch 7629566063dSJacob Faibussowitsch PetscCallAbort(MPI_COMM_NULL,boom()); // ERROR, communicator should be valid 7639566063dSJacob Faibussowitsch 7649566063dSJacob Faibussowitsch struct baz 7659566063dSJacob Faibussowitsch { 7669566063dSJacob Faibussowitsch baz() 7679566063dSJacob Faibussowitsch { 7689566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,boom()); // OK 7699566063dSJacob Faibussowitsch } 7709566063dSJacob Faibussowitsch 7719566063dSJacob Faibussowitsch ~baz() 7729566063dSJacob Faibussowitsch { 7739566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,boom()); // OK (in fact the only way to handle PETSc errors) 7749566063dSJacob Faibussowitsch } 7759566063dSJacob Faibussowitsch }; 7769566063dSJacob Faibussowitsch .ve 7779566063dSJacob Faibussowitsch 778989c49a2SBarry Smith Fortran Note: 779989c49a2SBarry Smith Use `PetscCallA()`. 780989c49a2SBarry Smith 781989c49a2SBarry Smith Developer Note: 782989c49a2SBarry Smith This should have the same name in Fortran as in C. 783989c49a2SBarry Smith 784db781477SPatrick Sanan .seealso: `SETERRABORT()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, 7855687f061SJacob Faibussowitsch `SETERRQ()`, `CHKMEMQ`, `PetscCallMPI()`, `PetscCallCXXAbort()` 7869566063dSJacob Faibussowitsch M*/ 7879566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 7889566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm, PetscErrorCode); 7899566063dSJacob Faibussowitsch void PetscCallContinue(PetscErrorCode); 7909566063dSJacob Faibussowitsch #else 7919371c9d4SSatish Balay #define PetscCallAbort(comm, ...) \ 7929371c9d4SSatish Balay do { \ 7937da6b3ccSLisandro Dalcin PetscErrorCode ierr_petsc_call_abort_; \ 7943ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 7957da6b3ccSLisandro Dalcin ierr_petsc_call_abort_ = __VA_ARGS__; \ 7963ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_abort_ != PETSC_SUCCESS)) { \ 7973ba16761SJacob Faibussowitsch ierr_petsc_call_abort_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_abort_, PETSC_ERROR_REPEAT, " "); \ 7983ba16761SJacob Faibussowitsch (void)MPI_Abort(comm, (PetscMPIInt)ierr_petsc_call_abort_); \ 7999566063dSJacob Faibussowitsch } \ 8009566063dSJacob Faibussowitsch } while (0) 8019371c9d4SSatish Balay #define PetscCallContinue(...) \ 8029371c9d4SSatish Balay do { \ 8037da6b3ccSLisandro Dalcin PetscErrorCode ierr_petsc_call_continue_; \ 8043ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 8057da6b3ccSLisandro Dalcin ierr_petsc_call_continue_ = __VA_ARGS__; \ 8063ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_continue_ != PETSC_SUCCESS)) { \ 8073ba16761SJacob Faibussowitsch ierr_petsc_call_continue_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_continue_, PETSC_ERROR_REPEAT, " "); \ 8083ba16761SJacob Faibussowitsch (void)ierr_petsc_call_continue_; \ 8093ba16761SJacob Faibussowitsch } \ 8109566063dSJacob Faibussowitsch } while (0) 8119566063dSJacob Faibussowitsch #endif 8129566063dSJacob Faibussowitsch 8139566063dSJacob Faibussowitsch /*MC 8149566063dSJacob Faibussowitsch CHKERRABORT - Checks error code returned from PETSc function. If non-zero it aborts immediately. 8159566063dSJacob Faibussowitsch 8169566063dSJacob Faibussowitsch Synopsis: 8179566063dSJacob Faibussowitsch #include <petscerror.h> 8189566063dSJacob Faibussowitsch void CHKERRABORT(MPI_Comm comm, PetscErrorCode ierr) 8199566063dSJacob Faibussowitsch 8209566063dSJacob Faibussowitsch Not Collective 8219566063dSJacob Faibussowitsch 8229566063dSJacob Faibussowitsch Input Parameters: 8239566063dSJacob Faibussowitsch + comm - the MPI communicator 8249566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 8259566063dSJacob Faibussowitsch 8269566063dSJacob Faibussowitsch Level: deprecated 8279566063dSJacob Faibussowitsch 828667f096bSBarry Smith Note: 829667f096bSBarry Smith Deprecated in favor of `PetscCallAbort()`. This routine behaves identically to it. 830667f096bSBarry Smith 831af27ebaaSBarry Smith .seealso: `PetscCallAbort()`, `PetscErrorCode` 8329566063dSJacob Faibussowitsch M*/ 8339566063dSJacob Faibussowitsch #define CHKERRABORT(comm, ...) PetscCallAbort(comm, __VA_ARGS__) 8349566063dSJacob Faibussowitsch #define CHKERRCONTINUE(...) PetscCallContinue(__VA_ARGS__) 8359566063dSJacob Faibussowitsch 8369566063dSJacob Faibussowitsch /*MC 83787497f52SBarry Smith CHKERRA - Fortran-only replacement for use of `CHKERRQ()` in the main program, which aborts immediately 838f388eb8bSPatrick Sanan 839f388eb8bSPatrick Sanan Synopsis: 840f388eb8bSPatrick Sanan #include <petscsys.h> 841f388eb8bSPatrick Sanan PetscErrorCode CHKERRA(PetscErrorCode ierr) 842f388eb8bSPatrick Sanan 843f388eb8bSPatrick Sanan Not Collective 844f388eb8bSPatrick Sanan 8452fe279fdSBarry Smith Input Parameter: 846f388eb8bSPatrick Sanan . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 847f388eb8bSPatrick Sanan 84887497f52SBarry Smith Level: deprecated 849f388eb8bSPatrick Sanan 85087497f52SBarry Smith Note: 85187497f52SBarry Smith This macro is rarely needed, normal usage is `PetscCallA()` in the main Fortran program. 852f388eb8bSPatrick Sanan 853989c49a2SBarry Smith Developer Note: 854989c49a2SBarry Smith Why isn't this named `CHKERRABORT()` in Fortran? 855989c49a2SBarry Smith 85687497f52SBarry Smith .seealso: `PetscCall()`, `PetscCallA()`, `PetscCallAbort()`, `CHKERRQ()`, `SETERRA()`, `SETERRQ()`, `SETERRABORT()` 857f388eb8bSPatrick Sanan M*/ 858f388eb8bSPatrick Sanan 8591e4f893aSSatish Balay PETSC_EXTERN PetscBool petscwaitonerrorflg; 8601e4f893aSSatish Balay PETSC_EXTERN PetscBool petscindebugger; 86135f00c14SToby Isaac PETSC_EXTERN PetscBool petscabortmpifinalize; 8627c66cc67SJunchao Zhang 8631fc6d13cSAlexander #if defined(PETSC_CLANG_STATIC_ANALYZER) 8641fc6d13cSAlexander void PETSCABORTWITHERR_Private(MPI_Comm, PetscErrorCode); 8651fc6d13cSAlexander #else 8661fc6d13cSAlexander #define PETSCABORTWITHIERR_Private(comm, ierr) \ 8671fc6d13cSAlexander do { \ 8681fc6d13cSAlexander PetscMPIInt size_; \ 869dd460d27SBarry Smith (void)MPI_Comm_size(comm, &size_); \ 8701fc6d13cSAlexander if (PetscCIEnabledPortableErrorOutput && (size_ == PetscGlobalSize || petscabortmpifinalize) && ierr != PETSC_ERR_SIG) { \ 871dd460d27SBarry Smith (void)MPI_Finalize(); \ 8721fc6d13cSAlexander exit(0); \ 8731fc6d13cSAlexander } else if (PetscCIEnabledPortableErrorOutput && PetscGlobalSize == 1) { \ 8741fc6d13cSAlexander exit(0); \ 8751fc6d13cSAlexander } else { \ 876dd460d27SBarry Smith (void)MPI_Abort(comm, (PetscMPIInt)ierr); \ 8771fc6d13cSAlexander } \ 8781fc6d13cSAlexander } while (0) 8791fc6d13cSAlexander #endif 8801fc6d13cSAlexander 8817c66cc67SJunchao Zhang /*MC 882667f096bSBarry Smith PETSCABORT - Call `MPI_Abort()` with an informative error code 8837c66cc67SJunchao Zhang 8847c66cc67SJunchao Zhang Synopsis: 8857c66cc67SJunchao Zhang #include <petscsys.h> 8867c66cc67SJunchao Zhang PETSCABORT(MPI_Comm comm, PetscErrorCode ierr) 8877c66cc67SJunchao Zhang 8887cdbe19fSJose E. Roman Collective; No Fortran Support 8897c66cc67SJunchao Zhang 8907c66cc67SJunchao Zhang Input Parameters: 89195bd0b28SBarry Smith + comm - An MPI communicator, so that the error can be collective 8927c66cc67SJunchao Zhang - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 8937c66cc67SJunchao Zhang 894bf4d2887SBarry Smith Level: advanced 8957c66cc67SJunchao Zhang 896bf4d2887SBarry Smith Notes: 897989c49a2SBarry Smith If the option `-start_in_debugger` was used then this calls `abort()` to stop the program in the debugger. 898bf4d2887SBarry Smith 899989c49a2SBarry Smith if `PetscCIEnabledPortableErrorOutput` is set, which means the code is running in the PETSc test harness (make test), 900989c49a2SBarry Smith and `comm` is `MPI_COMM_WORLD` it strives to exit cleanly without calling `MPI_Abort()` and instead calling `MPI_Finalize()`. 901660278c0SBarry Smith 902989c49a2SBarry Smith This is currently only used when an error propagates up to the C `main()` program and is detected by a `PetscCall()`, `PetscCallMPI()`, 903989c49a2SBarry Smith or is set in `main()` with `SETERRQ()`. Abort calls such as `SETERRABORT()`, 904989c49a2SBarry Smith `PetscCheckAbort()`, `PetscCallMPIAbort()`, and `PetscCallAbort()` always call `MPI_Abort()` and do not have any special 905989c49a2SBarry Smith handling for the test harness. 906989c49a2SBarry Smith 907989c49a2SBarry Smith Developer Note: 908989c49a2SBarry Smith Should the other abort calls also pass through this call instead of calling `MPI_Abort()` directly? 909989c49a2SBarry Smith 910989c49a2SBarry Smith .seealso: `PetscError()`, `PetscCall()`, `SETERRABORT()`, `PetscCheckAbort()`, `PetscCallMPIAbort()`, `PetscCall()`, `PetscCallMPI()`, 911af27ebaaSBarry Smith `PetscCallAbort()`, `MPI_Abort()`, `PetscErrorCode` 9127c66cc67SJunchao Zhang M*/ 91329f88feaSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 91429f88feaSJacob Faibussowitsch void PETSCABORT(MPI_Comm, PetscErrorCode); 91529f88feaSJacob Faibussowitsch #else 9169371c9d4SSatish Balay #define PETSCABORT(comm, ...) \ 9179371c9d4SSatish Balay do { \ 9183ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_abort_; \ 9193ba16761SJacob Faibussowitsch if (petscwaitonerrorflg) { ierr_petsc_abort_ = PetscSleep(1000); } \ 9203ba16761SJacob Faibussowitsch if (petscindebugger) { \ 9213ba16761SJacob Faibussowitsch abort(); \ 9223ba16761SJacob Faibussowitsch } else { \ 9233ba16761SJacob Faibussowitsch ierr_petsc_abort_ = __VA_ARGS__; \ 9241fc6d13cSAlexander PETSCABORTWITHIERR_Private(comm, ierr_petsc_abort_); \ 9253fcd9f07SJacob Faibussowitsch } \ 9267c66cc67SJunchao Zhang } while (0) 92729f88feaSJacob Faibussowitsch #endif 928986eef2eSBarry Smith 9299566063dSJacob Faibussowitsch #ifdef PETSC_CLANGUAGE_CXX 930986eef2eSBarry Smith /*MC 9319566063dSJacob Faibussowitsch PetscCallThrow - Checks error code, if non-zero it calls the C++ error handler which throws 9329566063dSJacob Faibussowitsch an exception 933986eef2eSBarry Smith 934986eef2eSBarry Smith Synopsis: 9359566063dSJacob Faibussowitsch #include <petscerror.h> 9369566063dSJacob Faibussowitsch void PetscCallThrow(PetscErrorCode ierr) 937986eef2eSBarry Smith 938986eef2eSBarry Smith Not Collective 939986eef2eSBarry Smith 9409566063dSJacob Faibussowitsch Input Parameter: 941986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 942986eef2eSBarry Smith 943667f096bSBarry Smith Level: beginner 944667f096bSBarry Smith 945986eef2eSBarry Smith Notes: 946af27ebaaSBarry Smith Requires PETSc to be configured with clanguage of c++. Throws a std::runtime_error() on error. 947986eef2eSBarry Smith 94887497f52SBarry Smith Once the error handler throws the exception you can use `PetscCallVoid()` which returns without 94987497f52SBarry Smith an error code (bad idea since the error is ignored) or `PetscCallAbort()` to have `MPI_Abort()` 9509566063dSJacob Faibussowitsch called immediately. 9519566063dSJacob Faibussowitsch 952db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscCall()`, `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, 953db781477SPatrick Sanan `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ` 954986eef2eSBarry Smith M*/ 9559371c9d4SSatish Balay #define PetscCallThrow(...) \ 9569371c9d4SSatish Balay do { \ 9573ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 9583ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_throw_ = __VA_ARGS__; \ 9593ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_throw_ != PETSC_SUCCESS)) PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_throw_, PETSC_ERROR_IN_CXX, PETSC_NULLPTR); \ 960ffc4695bSBarry Smith } while (0) 96185614651SBarry Smith 962cc26af49SMatthew Knepley /*MC 963cc26af49SMatthew Knepley CHKERRXX - Checks error code, if non-zero it calls the C++ error handler which throws an exception 964cc26af49SMatthew Knepley 965cc26af49SMatthew Knepley Synopsis: 9669566063dSJacob Faibussowitsch #include <petscerror.h> 9673af045c5SBarry Smith void CHKERRXX(PetscErrorCode ierr) 968cc26af49SMatthew Knepley 969eca87e8dSBarry Smith Not Collective 970cc26af49SMatthew Knepley 9719566063dSJacob Faibussowitsch Input Parameter: 9723af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 973cc26af49SMatthew Knepley 9749566063dSJacob Faibussowitsch Level: deprecated 975cc26af49SMatthew Knepley 976667f096bSBarry Smith Note: 977667f096bSBarry Smith Deprecated in favor of `PetscCallThrow()`. This routine behaves identically to it. 978667f096bSBarry Smith 979db781477SPatrick Sanan .seealso: `PetscCallThrow()` 980cc26af49SMatthew Knepley M*/ 9819566063dSJacob Faibussowitsch #define CHKERRXX(...) PetscCallThrow(__VA_ARGS__) 982fd705b32SMatthew Knepley #endif 983fd705b32SMatthew Knepley 9843ba16761SJacob Faibussowitsch #define PetscCallCXX_Private(__SETERR_FUNC__, __COMM__, ...) \ 9853ba16761SJacob Faibussowitsch do { \ 9863ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 9873ba16761SJacob Faibussowitsch try { \ 9883ba16761SJacob Faibussowitsch __VA_ARGS__; \ 9893ba16761SJacob Faibussowitsch } catch (const std::exception &e) { \ 9903ba16761SJacob Faibussowitsch __SETERR_FUNC__(__COMM__, PETSC_ERR_LIB, "%s", e.what()); \ 9913ba16761SJacob Faibussowitsch } \ 9923ba16761SJacob Faibussowitsch } while (0) 9933ba16761SJacob Faibussowitsch 9943f520e80SJunchao Zhang /*MC 9959566063dSJacob Faibussowitsch PetscCallCXX - Checks C++ function calls and if they throw an exception, catch it and then 9969566063dSJacob Faibussowitsch return a PETSc error code 9973f520e80SJunchao Zhang 9983f520e80SJunchao Zhang Synopsis: 9999566063dSJacob Faibussowitsch #include <petscerror.h> 10005687f061SJacob Faibussowitsch void PetscCallCXX(...) noexcept; 10013f520e80SJunchao Zhang 10023f520e80SJunchao Zhang Not Collective 10033f520e80SJunchao Zhang 10049566063dSJacob Faibussowitsch Input Parameter: 10055687f061SJacob Faibussowitsch . __VA_ARGS__ - An arbitrary expression 10065687f061SJacob Faibussowitsch 10075687f061SJacob Faibussowitsch Level: beginner 10083f520e80SJunchao Zhang 10093f520e80SJunchao Zhang Notes: 10105687f061SJacob Faibussowitsch `PetscCallCXX(...)` is a macro replacement for 10119566063dSJacob Faibussowitsch .vb 10129566063dSJacob Faibussowitsch try { 10135687f061SJacob Faibussowitsch __VA_ARGS__; 10149566063dSJacob Faibussowitsch } catch (const std::exception& e) { 10159566063dSJacob Faibussowitsch return ConvertToPetscErrorCode(e); 10169566063dSJacob Faibussowitsch } 10179566063dSJacob Faibussowitsch .ve 10189566063dSJacob Faibussowitsch Due to the fact that it catches any (reasonable) exception, it is essentially noexcept. 10193f520e80SJunchao Zhang 10205687f061SJacob Faibussowitsch If you cannot return a `PetscErrorCode` use `PetscCallCXXAbort()` instead. 10215687f061SJacob Faibussowitsch 10229566063dSJacob Faibussowitsch Example Usage: 10239566063dSJacob Faibussowitsch .vb 10249566063dSJacob Faibussowitsch void foo(void) { throw std::runtime_error("error"); } 10253f520e80SJunchao Zhang 10269566063dSJacob Faibussowitsch void bar() 10279566063dSJacob Faibussowitsch { 1028e8952933SJacob Faibussowitsch PetscCallCXX(foo()); // ERROR bar() does not return PetscErrorCode 10299566063dSJacob Faibussowitsch } 10309566063dSJacob Faibussowitsch 10319566063dSJacob Faibussowitsch PetscErrorCode baz() 10329566063dSJacob Faibussowitsch { 10339566063dSJacob Faibussowitsch PetscCallCXX(foo()); // OK 10349566063dSJacob Faibussowitsch 10359566063dSJacob Faibussowitsch PetscCallCXX( 10369566063dSJacob Faibussowitsch bar(); 1037d5b43468SJose E. Roman foo(); // OK multiple statements allowed 10389566063dSJacob Faibussowitsch ); 10399566063dSJacob Faibussowitsch } 10409566063dSJacob Faibussowitsch 10419566063dSJacob Faibussowitsch struct bop 10429566063dSJacob Faibussowitsch { 10439566063dSJacob Faibussowitsch bop() 10449566063dSJacob Faibussowitsch { 1045e8952933SJacob Faibussowitsch PetscCallCXX(foo()); // ERROR returns PetscErrorCode, cannot be used in constructors 10469566063dSJacob Faibussowitsch } 10479566063dSJacob Faibussowitsch }; 10489566063dSJacob Faibussowitsch 1049e8952933SJacob Faibussowitsch // ERROR contains do-while, cannot be used as function-try block 10509566063dSJacob Faibussowitsch PetscErrorCode qux() PetscCallCXX( 10519566063dSJacob Faibussowitsch bar(); 10529566063dSJacob Faibussowitsch baz(); 10539566063dSJacob Faibussowitsch foo(); 10549566063dSJacob Faibussowitsch return 0; 10559566063dSJacob Faibussowitsch ) 10569566063dSJacob Faibussowitsch .ve 10579566063dSJacob Faibussowitsch 10585687f061SJacob Faibussowitsch .seealso: `PetscCallCXXAbort()`, `PetscCallThrow()`, `SETERRQ()`, `PetscCall()`, 10595687f061SJacob Faibussowitsch `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 10605687f061SJacob Faibussowitsch `PetscError()`, `CHKMEMQ` 10613f520e80SJunchao Zhang M*/ 10623ba16761SJacob Faibussowitsch #define PetscCallCXX(...) PetscCallCXX_Private(SETERRQ, PETSC_COMM_SELF, __VA_ARGS__) 10635687f061SJacob Faibussowitsch 10645687f061SJacob Faibussowitsch /*MC 10655687f061SJacob Faibussowitsch PetscCallCXXAbort - Like `PetscCallCXX()` but calls `MPI_Abort()` instead of returning an 10665687f061SJacob Faibussowitsch error-code 10675687f061SJacob Faibussowitsch 10685687f061SJacob Faibussowitsch Synopsis: 10695687f061SJacob Faibussowitsch #include <petscerror.h> 10705687f061SJacob Faibussowitsch void PetscCallCXXAbort(MPI_Comm comm, ...) noexcept; 10715687f061SJacob Faibussowitsch 107220f4b53cSBarry Smith Collective; No Fortran Support 10735687f061SJacob Faibussowitsch 10745687f061SJacob Faibussowitsch Input Parameters: 10755687f061SJacob Faibussowitsch + comm - The MPI communicator to abort on 10765687f061SJacob Faibussowitsch - __VA_ARGS__ - An arbitrary expression 10775687f061SJacob Faibussowitsch 10785687f061SJacob Faibussowitsch Level: beginner 10795687f061SJacob Faibussowitsch 10805687f061SJacob Faibussowitsch Notes: 10815687f061SJacob Faibussowitsch This macro may be used to check C++ expressions for exceptions in cases where you cannot 10825687f061SJacob Faibussowitsch return an error code. This includes constructors, destructors, copy/move assignment functions 10835687f061SJacob Faibussowitsch or constructors among others. 10845687f061SJacob Faibussowitsch 10855687f061SJacob Faibussowitsch If an exception is caught, the macro calls `SETERRABORT()` on `comm`. The exception must 10865687f061SJacob Faibussowitsch derive from `std::exception` in order to be caught. 10875687f061SJacob Faibussowitsch 10885687f061SJacob Faibussowitsch If the routine _can_ return an error-code it is highly advised to use `PetscCallCXX()` 10895687f061SJacob Faibussowitsch instead. 10905687f061SJacob Faibussowitsch 10915687f061SJacob Faibussowitsch See `PetscCallCXX()` for additional discussion. 10925687f061SJacob Faibussowitsch 10935687f061SJacob Faibussowitsch Example Usage: 10945687f061SJacob Faibussowitsch .vb 10955687f061SJacob Faibussowitsch class Foo 10965687f061SJacob Faibussowitsch { 10975687f061SJacob Faibussowitsch std::vector<int> data_; 10985687f061SJacob Faibussowitsch 10995687f061SJacob Faibussowitsch public: 11005687f061SJacob Faibussowitsch // normally std::vector::reserve() may raise an exception, but since we handle it with 11015687f061SJacob Faibussowitsch // PetscCallCXXAbort() we may mark this routine as noexcept! 11025687f061SJacob Faibussowitsch Foo() noexcept 11035687f061SJacob Faibussowitsch { 11045687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, data_.reserve(10)); 11055687f061SJacob Faibussowitsch } 11065687f061SJacob Faibussowitsch }; 11075687f061SJacob Faibussowitsch 11085687f061SJacob Faibussowitsch std::vector<int> bar() 11095687f061SJacob Faibussowitsch { 11105687f061SJacob Faibussowitsch std::vector<int> v; 11115687f061SJacob Faibussowitsch 11125687f061SJacob Faibussowitsch PetscFunctionBegin; 11135687f061SJacob Faibussowitsch // OK! 11145687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1)); 11155687f061SJacob Faibussowitsch PetscFunctionReturn(v); 11165687f061SJacob Faibussowitsch } 11175687f061SJacob Faibussowitsch 11185687f061SJacob Faibussowitsch PetscErrorCode baz() 11195687f061SJacob Faibussowitsch { 11205687f061SJacob Faibussowitsch std::vector<int> v; 11215687f061SJacob Faibussowitsch 11225687f061SJacob Faibussowitsch PetscFunctionBegin; 11235687f061SJacob Faibussowitsch // WRONG! baz() returns a PetscErrorCode, prefer PetscCallCXX() instead 11245687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1)); 11253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11265687f061SJacob Faibussowitsch } 11275687f061SJacob Faibussowitsch .ve 11285687f061SJacob Faibussowitsch 11295687f061SJacob Faibussowitsch .seealso: `PetscCallCXX()`, `SETERRABORT()`, `PetscCallAbort()` 11305687f061SJacob Faibussowitsch M*/ 11313ba16761SJacob Faibussowitsch #define PetscCallCXXAbort(comm, ...) PetscCallCXX_Private(SETERRABORT, comm, __VA_ARGS__) 11323f520e80SJunchao Zhang 113330de9b25SBarry Smith /*MC 11349566063dSJacob Faibussowitsch CHKERRCXX - Checks C++ function calls and if they throw an exception, catch it and then 11359566063dSJacob Faibussowitsch return a PETSc error code 11369566063dSJacob Faibussowitsch 11379566063dSJacob Faibussowitsch Synopsis: 11389566063dSJacob Faibussowitsch #include <petscerror.h> 11399566063dSJacob Faibussowitsch void CHKERRCXX(func) noexcept; 11409566063dSJacob Faibussowitsch 11419566063dSJacob Faibussowitsch Not Collective 11429566063dSJacob Faibussowitsch 11439566063dSJacob Faibussowitsch Input Parameter: 11449566063dSJacob Faibussowitsch . func - C++ function calls 11459566063dSJacob Faibussowitsch 11469566063dSJacob Faibussowitsch Level: deprecated 11479566063dSJacob Faibussowitsch 1148667f096bSBarry Smith Note: 1149667f096bSBarry Smith Deprecated in favor of `PetscCallCXX()`. This routine behaves identically to it. 1150667f096bSBarry Smith 1151db781477SPatrick Sanan .seealso: `PetscCallCXX()` 11529566063dSJacob Faibussowitsch M*/ 11539566063dSJacob Faibussowitsch #define CHKERRCXX(...) PetscCallCXX(__VA_ARGS__) 11549566063dSJacob Faibussowitsch 11559566063dSJacob Faibussowitsch /*MC 115630de9b25SBarry Smith CHKMEMQ - Checks the memory for corruption, calls error handler if any is detected 115730de9b25SBarry Smith 115830de9b25SBarry Smith Synopsis: 1159aaa7dc30SBarry Smith #include <petscsys.h> 116091d3bdf4SKris Buschelman CHKMEMQ; 116130de9b25SBarry Smith 1162eca87e8dSBarry Smith Not Collective 1163eca87e8dSBarry Smith 116430de9b25SBarry Smith Level: beginner 116530de9b25SBarry Smith 116630de9b25SBarry Smith Notes: 1167af27ebaaSBarry Smith We recommend using Valgrind <https://petsc.org/release/faq/#valgrind> or for NVIDIA CUDA systems 1168af27ebaaSBarry Smith <https://docs.nvidia.com/cuda/cuda-memcheck/index.html> for finding memory problems. The ``CHKMEMQ`` macro is useful on systems that 11695ed36255SBarry Smith do not have valgrind, but is not as good as valgrind or cuda-memcheck. 11701957e957SBarry Smith 1171667f096bSBarry Smith Must run with the option `-malloc_debug` (`-malloc_test` in debug mode; or if `PetscMallocSetDebug()` called) to enable this option 117230de9b25SBarry Smith 117330de9b25SBarry Smith Once the error handler is called the calling function is then returned from with the given error code. 117430de9b25SBarry Smith 117530de9b25SBarry Smith By defaults prints location where memory that is corrupted was allocated. 117630de9b25SBarry Smith 11778af9f190SBarry Smith Use `CHKMEMA` for functions that return `void` 1178f621e05eSBarry Smith 1179db781477SPatrick Sanan .seealso: `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `SETERRQ()`, `PetscMallocValidate()` 118030de9b25SBarry Smith M*/ 11816d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 1182064a246eSJacob Faibussowitsch #define CHKMEMQ 1183064a246eSJacob Faibussowitsch #define CHKMEMA 11846d210af2SJacob Faibussowitsch #else 11859371c9d4SSatish Balay #define CHKMEMQ \ 11869371c9d4SSatish Balay do { \ 11873ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_memq_ = PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__); \ 11883ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_memq_ != PETSC_SUCCESS)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_memq_, PETSC_ERROR_REPEAT, " "); \ 118986d09637SLisandro Dalcin } while (0) 11906d210af2SJacob Faibussowitsch #define CHKMEMA PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__) 1191064a246eSJacob Faibussowitsch #endif 11929566063dSJacob Faibussowitsch 1193668f157eSBarry Smith /*E 1194668f157eSBarry Smith PetscErrorType - passed to the PETSc error handling routines indicating if this is the first or a later call to the error handlers 1195668f157eSBarry Smith 1196668f157eSBarry Smith Level: advanced 1197668f157eSBarry Smith 1198667f096bSBarry Smith Note: 119987497f52SBarry Smith `PETSC_ERROR_IN_CXX` indicates the error was detected in C++ and an exception should be generated 1200d736bfebSBarry Smith 1201667f096bSBarry Smith Developer Note: 1202667f096bSBarry Smith This is currently used to decide when to print the detailed information about the run in `PetscTraceBackErrorHandler()` 1203668f157eSBarry Smith 120487497f52SBarry Smith .seealso: `PetscError()`, `SETERRQ()` 1205668f157eSBarry Smith E*/ 12069371c9d4SSatish Balay typedef enum { 12079371c9d4SSatish Balay PETSC_ERROR_INITIAL = 0, 12089371c9d4SSatish Balay PETSC_ERROR_REPEAT = 1, 12099371c9d4SSatish Balay PETSC_ERROR_IN_CXX = 2 12109371c9d4SSatish Balay } PetscErrorType; 12114b209cf6SBarry Smith 1212eb9e708aSLisandro Dalcin #if defined(__clang_analyzer__) 1213eb9e708aSLisandro Dalcin __attribute__((analyzer_noreturn)) 1214eb9e708aSLisandro Dalcin #endif 12159371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode 12169371c9d4SSatish Balay PetscError(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, ...) PETSC_ATTRIBUTE_COLD PETSC_ATTRIBUTE_FORMAT(7, 8); 1217eb9e708aSLisandro Dalcin 1218014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscErrorPrintfInitialize(void); 12193ba16761SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscErrorMessage(PetscErrorCode, const char *[], char **); 1220d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscTraceBackErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1221d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscIgnoreErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1222d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscEmacsClientErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1223d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscMPIAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1224d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1225d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAttachDebuggerErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1226d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscReturnErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1227efca3c55SSatish Balay PETSC_EXTERN PetscErrorCode PetscPushErrorHandler(PetscErrorCode (*handler)(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *), void *); 1228014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopErrorHandler(void); 12298d359177SBarry Smith PETSC_EXTERN PetscErrorCode PetscSignalHandlerDefault(int, void *); 1230014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPushSignalHandler(PetscErrorCode (*)(int, void *), void *); 1231014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopSignalHandler(void); 123228559dc8SJed Brown PETSC_EXTERN PetscErrorCode PetscCheckPointerSetIntensity(PetscInt); 1233c2a741eeSJunchao Zhang PETSC_EXTERN void PetscSignalSegvCheckPointerOrMpi(void); 1234edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 13, 0, "PetscSignalSegvCheckPointerOrMpi()", ) static inline void PetscSignalSegvCheckPointer(void) 1235d71ae5a4SJacob Faibussowitsch { 12369371c9d4SSatish Balay PetscSignalSegvCheckPointerOrMpi(); 12379371c9d4SSatish Balay } 1238329f5518SBarry Smith 1239639ff905SBarry Smith /*MC 1240639ff905SBarry Smith PetscErrorPrintf - Prints error messages. 1241639ff905SBarry Smith 1242639ff905SBarry Smith Synopsis: 1243aaa7dc30SBarry Smith #include <petscsys.h> 1244639ff905SBarry Smith PetscErrorCode (*PetscErrorPrintf)(const char format[], ...); 1245639ff905SBarry Smith 12467cdbe19fSJose E. Roman Not Collective; No Fortran Support 12477cdbe19fSJose E. Roman 1248f899ff85SJose E. Roman Input Parameter: 1249cd05f99aSJacob Faibussowitsch . format - the usual `printf()` format string 1250639ff905SBarry Smith 1251639ff905SBarry Smith Options Database Keys: 12521957e957SBarry Smith + -error_output_stdout - cause error messages to be printed to stdout instead of the (default) stderr 1253e1bc860dSBarry Smith - -error_output_none - to turn off all printing of error messages (does not change the way the error is handled.) 1254639ff905SBarry Smith 1255cf53795eSBarry Smith Level: developer 1256cf53795eSBarry Smith 125795452b02SPatrick Sanan Notes: 125895452b02SPatrick Sanan Use 1259667f096bSBarry Smith .vb 126095bd0b28SBarry Smith PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the error is handled) and 1261667f096bSBarry Smith PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on or you can use your own function 1262667f096bSBarry Smith .ve 1263639ff905SBarry Smith Use 1264667f096bSBarry Smith .vb 126587497f52SBarry Smith `PETSC_STDERR` = FILE* obtained from a file open etc. to have stderr printed to the file. 126687497f52SBarry Smith `PETSC_STDOUT` = FILE* obtained from a file open etc. to have stdout printed to the file. 1267667f096bSBarry Smith .ve 1268639ff905SBarry Smith Use 1269af27ebaaSBarry Smith .vb 127087497f52SBarry Smith `PetscPushErrorHandler()` to provide your own error handler that determines what kind of messages to print 1271af27ebaaSBarry Smith .ve 1272639ff905SBarry Smith 1273db781477SPatrick Sanan .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscHelpPrintf()`, `PetscPrintf()`, `PetscPushErrorHandler()`, `PetscVFPrintf()`, `PetscHelpPrintf()` 1274639ff905SBarry Smith M*/ 12753ca90d2dSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode (*PetscErrorPrintf)(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2); 1276639ff905SBarry Smith 1277cf0818bdSBarry Smith /*E 1278cf0818bdSBarry Smith PetscFPTrap - types of floating point exceptions that may be trapped 1279cf0818bdSBarry Smith 128087497f52SBarry Smith Currently only `PETSC_FP_TRAP_OFF` and `PETSC_FP_TRAP_ON` are handled. All others are treated as `PETSC_FP_TRAP_ON`. 1281cf0818bdSBarry Smith 1282cf0818bdSBarry Smith Level: intermediate 1283cf0818bdSBarry Smith 12847de69702SBarry Smith .seealso: `PetscSetFPTrap()`, `PetscFPTrapPush()` 1285cf0818bdSBarry Smith E*/ 12869371c9d4SSatish Balay typedef enum { 12879371c9d4SSatish Balay PETSC_FP_TRAP_OFF = 0, 12889371c9d4SSatish Balay PETSC_FP_TRAP_INDIV = 1, 12899371c9d4SSatish Balay PETSC_FP_TRAP_FLTOPERR = 2, 12909371c9d4SSatish Balay PETSC_FP_TRAP_FLTOVF = 4, 12919371c9d4SSatish Balay PETSC_FP_TRAP_FLTUND = 8, 12929371c9d4SSatish Balay PETSC_FP_TRAP_FLTDIV = 16, 1293dd460d27SBarry Smith PETSC_FP_TRAP_FLTINEX = 32, 1294dd460d27SBarry Smith PETSC_FP_TRAP_ON = (PETSC_FP_TRAP_INDIV | PETSC_FP_TRAP_FLTOPERR | PETSC_FP_TRAP_FLTOVF | PETSC_FP_TRAP_FLTDIV | PETSC_FP_TRAP_FLTINEX) 12959371c9d4SSatish Balay } PetscFPTrap; 1296dd460d27SBarry Smith 1297014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSetFPTrap(PetscFPTrap); 1298014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPush(PetscFPTrap); 1299014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPop(void); 1300aba4c478SBarry Smith PETSC_EXTERN PetscErrorCode PetscDetermineInitialFPTrap(void); 130154a8ef01SBarry Smith 13023a40ed3dSBarry Smith /* 13033a40ed3dSBarry Smith Allows the code to build a stack frame as it runs 13043a40ed3dSBarry Smith */ 13053a40ed3dSBarry Smith 130699cd645aSJed Brown #define PETSCSTACKSIZE 64 13073a40ed3dSBarry Smith typedef struct { 13080e33f6ddSBarry Smith const char *function[PETSCSTACKSIZE]; 13090e33f6ddSBarry Smith const char *file[PETSCSTACKSIZE]; 1310184914b5SBarry Smith int line[PETSCSTACKSIZE]; 1311362febeeSStefano Zampini int petscroutine[PETSCSTACKSIZE]; /* 0 external called from petsc, 1 petsc functions, 2 petsc user functions */ 1312184914b5SBarry Smith int currentsize; 1313a2f94806SJed Brown int hotdepth; 13144be741a6SBarry Smith PetscBool check; /* option to check for correct Push/Pop semantics, true for default petscstack but not other stacks */ 13153a40ed3dSBarry Smith } PetscStack; 1316dfb7d7afSStefano Zampini #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY) 131727104ee2SJacob Faibussowitsch PETSC_EXTERN PetscStack petscstack; 131827104ee2SJacob Faibussowitsch #endif 13193a40ed3dSBarry Smith 13205d12eec7SSatish Balay #if defined(PETSC_SERIALIZE_FUNCTIONS) 13215d12eec7SSatish Balay #include <petsc/private/petscfptimpl.h> 13225d12eec7SSatish Balay /* 13235d12eec7SSatish Balay Registers the current function into the global function pointer to function name table 13245d12eec7SSatish Balay 13255d12eec7SSatish Balay Have to fix this to handle errors but cannot return error since used in PETSC_VIEWER_DRAW_() etc 13265d12eec7SSatish Balay */ 13279371c9d4SSatish Balay #define PetscRegister__FUNCT__() \ 13289371c9d4SSatish Balay do { \ 13295d12eec7SSatish Balay static PetscBool __chked = PETSC_FALSE; \ 13305d12eec7SSatish Balay if (!__chked) { \ 13319371c9d4SSatish Balay void *ptr; \ 13323ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscDLSym(NULL, PETSC_FUNCTION_NAME, &ptr)); \ 13335d12eec7SSatish Balay __chked = PETSC_TRUE; \ 13349371c9d4SSatish Balay } \ 13359371c9d4SSatish Balay } while (0) 13365d12eec7SSatish Balay #else 13375d12eec7SSatish Balay #define PetscRegister__FUNCT__() 13385d12eec7SSatish Balay #endif 13395d12eec7SSatish Balay 1340eae3dc7dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) || defined(__clang_analyzer__) 13416d210af2SJacob Faibussowitsch #define PetscStackPushNoCheck(funct, petsc_routine, hot) 134237154d10SBarry Smith #define PetscStackUpdateLine 1343792fecdfSBarry Smith #define PetscStackPushExternal(funct) 1344dd460d27SBarry Smith #define PetscStackPopNoCheck(funct) 13456d210af2SJacob Faibussowitsch #define PetscStackClearTop 13466d210af2SJacob Faibussowitsch #define PetscFunctionBegin 13476d210af2SJacob Faibussowitsch #define PetscFunctionBeginUser 13486d210af2SJacob Faibussowitsch #define PetscFunctionBeginHot 13490a57284eSJacob Faibussowitsch #define PetscFunctionReturn(...) return __VA_ARGS__ 13506d210af2SJacob Faibussowitsch #define PetscFunctionReturnVoid() return 13516d210af2SJacob Faibussowitsch #define PetscStackPop 13526d210af2SJacob Faibussowitsch #define PetscStackPush(f) 1353*63bfac88SBarry Smith #define PetscStackPush_Private(stack__, file__, func__, line__, petsc_routine__, hot__) \ 1354*63bfac88SBarry Smith (void)file__; \ 1355*63bfac88SBarry Smith (void)func__; \ 1356*63bfac88SBarry Smith (void)line__ 1357*63bfac88SBarry Smith #define PetscStackPop_Private(stack__, func__) (void)func__ 1358dfb7d7afSStefano Zampini #elif defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY) 1359660278c0SBarry Smith 13609371c9d4SSatish Balay #define PetscStackPush_Private(stack__, file__, func__, line__, petsc_routine__, hot__) \ 13619371c9d4SSatish Balay do { \ 13625a96b57dSJacob Faibussowitsch if (stack__.currentsize < PETSCSTACKSIZE) { \ 13635a96b57dSJacob Faibussowitsch stack__.function[stack__.currentsize] = func__; \ 1364ef1023bdSBarry Smith if (petsc_routine__) { \ 1365ef1023bdSBarry Smith stack__.file[stack__.currentsize] = file__; \ 13665a96b57dSJacob Faibussowitsch stack__.line[stack__.currentsize] = line__; \ 1367ef1023bdSBarry Smith } else { \ 1368648bc8c4SBarry Smith stack__.file[stack__.currentsize] = PETSC_NULLPTR; \ 1369ef1023bdSBarry Smith stack__.line[stack__.currentsize] = 0; \ 1370ef1023bdSBarry Smith } \ 13715a96b57dSJacob Faibussowitsch stack__.petscroutine[stack__.currentsize] = petsc_routine__; \ 13725a96b57dSJacob Faibussowitsch } \ 13735a96b57dSJacob Faibussowitsch ++stack__.currentsize; \ 13745a96b57dSJacob Faibussowitsch stack__.hotdepth += (hot__ || stack__.hotdepth); \ 13755a96b57dSJacob Faibussowitsch } while (0) 13765a96b57dSJacob Faibussowitsch 13774be741a6SBarry Smith /* uses PetscCheckAbort() because may be used in a function that does not return an error code */ 13789371c9d4SSatish Balay #define PetscStackPop_Private(stack__, func__) \ 13799371c9d4SSatish Balay do { \ 13804be741a6SBarry 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__); \ 13815a96b57dSJacob Faibussowitsch if (--stack__.currentsize < PETSCSTACKSIZE) { \ 13829371c9d4SSatish 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", \ 13839371c9d4SSatish Balay stack__.function[stack__.currentsize], stack__.file[stack__.currentsize], stack__.line[stack__.currentsize], func__, __FILE__, __LINE__); \ 13845a96b57dSJacob Faibussowitsch stack__.function[stack__.currentsize] = PETSC_NULLPTR; \ 13855a96b57dSJacob Faibussowitsch stack__.file[stack__.currentsize] = PETSC_NULLPTR; \ 13865a96b57dSJacob Faibussowitsch stack__.line[stack__.currentsize] = 0; \ 13875a96b57dSJacob Faibussowitsch stack__.petscroutine[stack__.currentsize] = 0; \ 13885a96b57dSJacob Faibussowitsch } \ 13895a96b57dSJacob Faibussowitsch stack__.hotdepth = PetscMax(stack__.hotdepth - 1, 0); \ 13905a96b57dSJacob Faibussowitsch } while (0) 13915a96b57dSJacob Faibussowitsch 1392586f9135SBarry Smith /*MC 1393586f9135SBarry Smith PetscStackPushNoCheck - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is 1394586f9135SBarry Smith currently in the source code. 1395586f9135SBarry Smith 1396586f9135SBarry Smith Synopsis: 1397586f9135SBarry Smith #include <petscsys.h> 1398586f9135SBarry Smith void PetscStackPushNoCheck(char *funct,int petsc_routine,PetscBool hot); 1399586f9135SBarry Smith 14007cdbe19fSJose E. Roman Not Collective 14017cdbe19fSJose E. Roman 1402586f9135SBarry Smith Input Parameters: 1403586f9135SBarry Smith + funct - the function name 1404586f9135SBarry Smith . petsc_routine - 2 user function, 1 PETSc function, 0 some other function 1405586f9135SBarry Smith - hot - indicates that the function may be called often so expensive error checking should be turned off inside the function 1406586f9135SBarry Smith 1407586f9135SBarry Smith Level: developer 1408586f9135SBarry Smith 1409586f9135SBarry Smith Notes: 1410586f9135SBarry 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 141187497f52SBarry 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 1412586f9135SBarry Smith help debug the problem. 1413586f9135SBarry Smith 1414ef1023bdSBarry Smith This version does not check the memory corruption (an expensive operation), use `PetscStackPush()` to check the memory. 1415ef1023bdSBarry Smith 1416792fecdfSBarry Smith Use `PetscStackPushExternal()` for a function call that is about to be made to a non-PETSc or user function (such as BLAS etc). 1417ef1023bdSBarry Smith 1418586f9135SBarry Smith The default stack is a global variable called `petscstack`. 1419586f9135SBarry Smith 1420586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1421ef1023bdSBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPush()`, `PetscStackPop`, 1422792fecdfSBarry Smith `PetscStackPushExternal()` 1423586f9135SBarry Smith M*/ 14249371c9d4SSatish Balay #define PetscStackPushNoCheck(funct, petsc_routine, hot) \ 14259371c9d4SSatish Balay do { \ 1426e04113cfSBarry Smith PetscStackSAWsTakeAccess(); \ 14275a96b57dSJacob Faibussowitsch PetscStackPush_Private(petscstack, __FILE__, funct, __LINE__, petsc_routine, hot); \ 1428e04113cfSBarry Smith PetscStackSAWsGrantAccess(); \ 1429441dd030SJed Brown } while (0) 1430441dd030SJed Brown 1431586f9135SBarry Smith /*MC 143287497f52SBarry Smith PetscStackUpdateLine - in a function that has a `PetscFunctionBegin` or `PetscFunctionBeginUser` updates the stack line number to the 143337154d10SBarry Smith current line number. 143437154d10SBarry Smith 143537154d10SBarry Smith Synopsis: 143637154d10SBarry Smith #include <petscsys.h> 143737154d10SBarry Smith void PetscStackUpdateLine 143837154d10SBarry Smith 14397cdbe19fSJose E. Roman Not Collective 14407cdbe19fSJose E. Roman 144137154d10SBarry Smith Level: developer 144237154d10SBarry Smith 144337154d10SBarry Smith Notes: 144487497f52SBarry Smith Using `PetscCall()` and friends automatically handles this process 144587497f52SBarry Smith 144637154d10SBarry 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 144737154d10SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 144837154d10SBarry Smith help debug the problem. 144937154d10SBarry Smith 145095bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 145137154d10SBarry Smith 145237154d10SBarry Smith This is used by `PetscCall()` and is otherwise not like to be needed 145337154d10SBarry Smith 145437154d10SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`, `PetscCall()` 145537154d10SBarry Smith M*/ 145637154d10SBarry Smith #define PetscStackUpdateLine \ 14573ba16761SJacob Faibussowitsch do { \ 1458f1e99121SPierre Jolivet if (petscstack.currentsize > 0 && petscstack.currentsize < PETSCSTACKSIZE && petscstack.function[petscstack.currentsize - 1] == PETSC_FUNCTION_NAME) { petscstack.line[petscstack.currentsize - 1] = __LINE__; } \ 14593ba16761SJacob Faibussowitsch } while (0) 146037154d10SBarry Smith 146137154d10SBarry Smith /*MC 1462792fecdfSBarry Smith PetscStackPushExternal - Pushes a new function name onto the PETSc default stack that tracks where the running program is 1463ef1023bdSBarry Smith currently in the source code. Does not include the filename or line number since this is called by the calling routine 1464ef1023bdSBarry Smith for non-PETSc or user functions. 1465ef1023bdSBarry Smith 1466ef1023bdSBarry Smith Synopsis: 1467ef1023bdSBarry Smith #include <petscsys.h> 1468792fecdfSBarry Smith void PetscStackPushExternal(char *funct); 1469ef1023bdSBarry Smith 14707cdbe19fSJose E. Roman Not Collective 14717cdbe19fSJose E. Roman 14722fe279fdSBarry Smith Input Parameter: 1473ef1023bdSBarry Smith . funct - the function name 1474ef1023bdSBarry Smith 1475ef1023bdSBarry Smith Level: developer 1476ef1023bdSBarry Smith 1477ef1023bdSBarry Smith Notes: 147887497f52SBarry Smith Using `PetscCallExternal()` and friends automatically handles this process 147987497f52SBarry Smith 1480ef1023bdSBarry 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 1481ef1023bdSBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1482ef1023bdSBarry Smith help debug the problem. 1483ef1023bdSBarry Smith 1484ef1023bdSBarry Smith The default stack is a global variable called `petscstack`. 1485ef1023bdSBarry Smith 1486ef1023bdSBarry Smith This is to be used when calling an external package function such as a BLAS function. 1487ef1023bdSBarry Smith 1488ef1023bdSBarry Smith This also updates the stack line number for the current stack function. 1489ef1023bdSBarry Smith 1490ef1023bdSBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1491ef1023bdSBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1492ef1023bdSBarry Smith M*/ 14939371c9d4SSatish Balay #define PetscStackPushExternal(funct) \ 14949371c9d4SSatish Balay do { \ 14959371c9d4SSatish Balay PetscStackUpdateLine; \ 14969371c9d4SSatish Balay PetscStackPushNoCheck(funct, 0, PETSC_TRUE); \ 1497a8f51744SPierre Jolivet } while (0) 1498ef1023bdSBarry Smith 1499ef1023bdSBarry Smith /*MC 1500586f9135SBarry Smith PetscStackPopNoCheck - Pops a function name from the PETSc default stack that tracks where the running program is 1501586f9135SBarry Smith currently in the source code. 1502586f9135SBarry Smith 1503586f9135SBarry Smith Synopsis: 1504586f9135SBarry Smith #include <petscsys.h> 1505586f9135SBarry Smith void PetscStackPopNoCheck(char *funct); 1506586f9135SBarry Smith 15077cdbe19fSJose E. Roman Not Collective 15087cdbe19fSJose E. Roman 1509586f9135SBarry Smith Input Parameter: 1510586f9135SBarry Smith . funct - the function name 1511586f9135SBarry Smith 1512586f9135SBarry Smith Level: developer 1513586f9135SBarry Smith 1514586f9135SBarry Smith Notes: 151587497f52SBarry Smith Using `PetscCall()`, `PetscCallExternal()`, `PetscCallBack()` and friends negates the need to call this 151687497f52SBarry Smith 1517586f9135SBarry 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 1518586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1519586f9135SBarry Smith help debug the problem. 1520586f9135SBarry Smith 152195bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 1522586f9135SBarry Smith 1523586f9135SBarry Smith Developer Note: 1524586f9135SBarry Smith `PetscStackPopNoCheck()` takes a function argument while `PetscStackPop` does not, this difference is likely just historical. 1525586f9135SBarry Smith 1526586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1527586f9135SBarry Smith M*/ 15289371c9d4SSatish Balay #define PetscStackPopNoCheck(funct) \ 15299371c9d4SSatish Balay do { \ 1530362febeeSStefano Zampini PetscStackSAWsTakeAccess(); \ 15315a96b57dSJacob Faibussowitsch PetscStackPop_Private(petscstack, funct); \ 1532362febeeSStefano Zampini PetscStackSAWsGrantAccess(); \ 1533362febeeSStefano Zampini } while (0) 1534362febeeSStefano Zampini 15359371c9d4SSatish Balay #define PetscStackClearTop \ 15369371c9d4SSatish Balay do { \ 1537e04113cfSBarry Smith PetscStackSAWsTakeAccess(); \ 15389371c9d4SSatish Balay if (petscstack.currentsize > 0 && --petscstack.currentsize < PETSCSTACKSIZE) { \ 153927104ee2SJacob Faibussowitsch petscstack.function[petscstack.currentsize] = PETSC_NULLPTR; \ 154027104ee2SJacob Faibussowitsch petscstack.file[petscstack.currentsize] = PETSC_NULLPTR; \ 154127104ee2SJacob Faibussowitsch petscstack.line[petscstack.currentsize] = 0; \ 154227104ee2SJacob Faibussowitsch petscstack.petscroutine[petscstack.currentsize] = 0; \ 1543441dd030SJed Brown } \ 154427104ee2SJacob Faibussowitsch petscstack.hotdepth = PetscMax(petscstack.hotdepth - 1, 0); \ 1545e04113cfSBarry Smith PetscStackSAWsGrantAccess(); \ 1546441dd030SJed Brown } while (0) 1547441dd030SJed Brown 154830de9b25SBarry Smith /*MC 15491957e957SBarry Smith PetscFunctionBegin - First executable line of each PETSc function, used for error handling. Final 155087497f52SBarry Smith line of PETSc functions should be `PetscFunctionReturn`(0); 155130de9b25SBarry Smith 155230de9b25SBarry Smith Synopsis: 1553aaa7dc30SBarry Smith #include <petscsys.h> 155430de9b25SBarry Smith void PetscFunctionBegin; 155530de9b25SBarry Smith 155620f4b53cSBarry Smith Not Collective; No Fortran Support 1557eca87e8dSBarry Smith 155830de9b25SBarry Smith Usage: 155930de9b25SBarry Smith .vb 156030de9b25SBarry Smith int something; 156130de9b25SBarry Smith 156230de9b25SBarry Smith PetscFunctionBegin; 156330de9b25SBarry Smith .ve 156430de9b25SBarry Smith 156530de9b25SBarry Smith Level: developer 156630de9b25SBarry Smith 156720f4b53cSBarry Smith Note: 156820f4b53cSBarry Smith Use `PetscFunctionBeginUser` for application codes. 156920f4b53cSBarry Smith 1570586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()` 157130de9b25SBarry Smith 157230de9b25SBarry Smith M*/ 15739371c9d4SSatish Balay #define PetscFunctionBegin \ 15749371c9d4SSatish Balay do { \ 1575362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_FALSE); \ 1576a2f94806SJed Brown PetscRegister__FUNCT__(); \ 1577a2f94806SJed Brown } while (0) 1578a2f94806SJed Brown 1579a2f94806SJed Brown /*MC 158087497f52SBarry Smith PetscFunctionBeginHot - Substitute for `PetscFunctionBegin` to be used in functions that are called in 1581a2f94806SJed Brown performance-critical circumstances. Use of this function allows for lighter profiling by default. 1582a2f94806SJed Brown 1583a2f94806SJed Brown Synopsis: 1584aaa7dc30SBarry Smith #include <petscsys.h> 1585a2f94806SJed Brown void PetscFunctionBeginHot; 1586a2f94806SJed Brown 158720f4b53cSBarry Smith Not Collective; No Fortran Support 1588a2f94806SJed Brown 1589a2f94806SJed Brown Usage: 1590a2f94806SJed Brown .vb 1591a2f94806SJed Brown int something; 1592a2f94806SJed Brown 1593a2f94806SJed Brown PetscFunctionBeginHot; 1594a2f94806SJed Brown .ve 1595a2f94806SJed Brown 1596a2f94806SJed Brown Level: developer 1597a2f94806SJed Brown 1598586f9135SBarry Smith .seealso: `PetscFunctionBegin`, `PetscFunctionReturn()`, `PetscStackPushNoCheck()` 1599a2f94806SJed Brown 1600a2f94806SJed Brown M*/ 16019371c9d4SSatish Balay #define PetscFunctionBeginHot \ 16029371c9d4SSatish Balay do { \ 1603362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_TRUE); \ 16042d53ad75SBarry Smith PetscRegister__FUNCT__(); \ 160553c77d0aSJed Brown } while (0) 160653c77d0aSJed Brown 1607a8d2bbe5SBarry Smith /*MC 1608530d85c1SBarry Smith PetscFunctionBeginUser - First executable line of user provided routines 1609a8d2bbe5SBarry Smith 1610a8d2bbe5SBarry Smith Synopsis: 1611aaa7dc30SBarry Smith #include <petscsys.h> 1612a8d2bbe5SBarry Smith void PetscFunctionBeginUser; 1613a8d2bbe5SBarry Smith 161420f4b53cSBarry Smith Not Collective; No Fortran Support 1615a8d2bbe5SBarry Smith 1616a8d2bbe5SBarry Smith Usage: 1617a8d2bbe5SBarry Smith .vb 1618a8d2bbe5SBarry Smith int something; 1619a8d2bbe5SBarry Smith 1620ac285190SBarry Smith PetscFunctionBeginUser; 1621a8d2bbe5SBarry Smith .ve 1622a8d2bbe5SBarry Smith 162320f4b53cSBarry Smith Level: intermediate 162420f4b53cSBarry Smith 1625a8d2bbe5SBarry Smith Notes: 1626530d85c1SBarry Smith Functions that incorporate this must call `PetscFunctionReturn()` instead of return except for main(). 1627530d85c1SBarry Smith 1628530d85c1SBarry Smith May be used before `PetscInitialize()` 16291957e957SBarry Smith 1630530d85c1SBarry Smith This is identical to `PetscFunctionBegin` except it labels the routine as a user 1631ac285190SBarry Smith routine instead of as a PETSc library routine. 1632ac285190SBarry Smith 1633586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, `PetscFunctionBeginHot`, `PetscStackPushNoCheck()` 1634a8d2bbe5SBarry Smith M*/ 16359371c9d4SSatish Balay #define PetscFunctionBeginUser \ 16369371c9d4SSatish Balay do { \ 1637362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 2, PETSC_FALSE); \ 1638a8d2bbe5SBarry Smith PetscRegister__FUNCT__(); \ 1639a8d2bbe5SBarry Smith } while (0) 1640a8d2bbe5SBarry Smith 1641586f9135SBarry Smith /*MC 1642586f9135SBarry Smith PetscStackPush - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is 1643586f9135SBarry Smith currently in the source code and verifies the memory is not corrupted. 1644586f9135SBarry Smith 1645586f9135SBarry Smith Synopsis: 1646586f9135SBarry Smith #include <petscsys.h> 1647586f9135SBarry Smith void PetscStackPush(char *funct) 1648586f9135SBarry Smith 16497cdbe19fSJose E. Roman Not Collective 16507cdbe19fSJose E. Roman 1651586f9135SBarry Smith Input Parameter: 1652586f9135SBarry Smith . funct - the function name 1653586f9135SBarry Smith 1654586f9135SBarry Smith Level: developer 1655586f9135SBarry Smith 1656586f9135SBarry Smith Notes: 1657586f9135SBarry 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 1658586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1659586f9135SBarry Smith help debug the problem. 1660586f9135SBarry Smith 166195bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 1662586f9135SBarry Smith 1663586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1664586f9135SBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1665586f9135SBarry Smith M*/ 16669371c9d4SSatish Balay #define PetscStackPush(n) \ 16679371c9d4SSatish Balay do { \ 1668362febeeSStefano Zampini PetscStackPushNoCheck(n, 0, PETSC_FALSE); \ 166915681b3cSBarry Smith CHKMEMQ; \ 167015681b3cSBarry Smith } while (0) 16713a40ed3dSBarry Smith 1672586f9135SBarry Smith /*MC 1673586f9135SBarry Smith PetscStackPop - Pops a function name from the PETSc default stack that tracks where the running program is 1674586f9135SBarry Smith currently in the source code and verifies the memory is not corrupted. 1675586f9135SBarry Smith 1676586f9135SBarry Smith Synopsis: 1677586f9135SBarry Smith #include <petscsys.h> 1678586f9135SBarry Smith void PetscStackPop 1679586f9135SBarry Smith 16807cdbe19fSJose E. Roman Not Collective 16817cdbe19fSJose E. Roman 1682586f9135SBarry Smith Level: developer 1683586f9135SBarry Smith 1684586f9135SBarry Smith Notes: 1685586f9135SBarry 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 1686586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1687586f9135SBarry Smith help debug the problem. 1688586f9135SBarry Smith 168995bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 1690586f9135SBarry Smith 1691586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPopNoCheck()`, `PetscStackPush()` 1692586f9135SBarry Smith M*/ 16939371c9d4SSatish Balay #define PetscStackPop \ 16949371c9d4SSatish Balay do { \ 1695441dd030SJed Brown CHKMEMQ; \ 1696362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 169715681b3cSBarry Smith } while (0) 1698d64ed03dSBarry Smith 169930de9b25SBarry Smith /*MC 17000a57284eSJacob Faibussowitsch PetscFunctionReturn - Last executable line of each PETSc function used for error 17010a57284eSJacob Faibussowitsch handling. Replaces `return()`. 170230de9b25SBarry Smith 170330de9b25SBarry Smith Synopsis: 17040a57284eSJacob Faibussowitsch #include <petscerror.h> 17050a57284eSJacob Faibussowitsch void PetscFunctionReturn(...) 170630de9b25SBarry Smith 1707667f096bSBarry Smith Not Collective; No Fortran Support 1708eca87e8dSBarry Smith 17095687f061SJacob Faibussowitsch Level: beginner 17105687f061SJacob Faibussowitsch 17110a57284eSJacob Faibussowitsch Notes: 17120a57284eSJacob Faibussowitsch This routine is a macro, so while it does not "return" anything itself, it does return from 17130a57284eSJacob Faibussowitsch the function in the literal sense. 17140a57284eSJacob Faibussowitsch 17150a57284eSJacob Faibussowitsch Usually the return value is the integer literal `0` (for example in any function returning 17160a57284eSJacob Faibussowitsch `PetscErrorCode`), however it is possible to return any arbitrary type. The arguments of 17170a57284eSJacob Faibussowitsch this macro are placed before the `return` statement as-is. 17180a57284eSJacob Faibussowitsch 17190a57284eSJacob Faibussowitsch Any routine which returns via `PetscFunctionReturn()` must begin with a corresponding 17200a57284eSJacob Faibussowitsch `PetscFunctionBegin`. 17210a57284eSJacob Faibussowitsch 17220a57284eSJacob Faibussowitsch For routines which return `void` use `PetscFunctionReturnVoid()` instead. 17230a57284eSJacob Faibussowitsch 17240a57284eSJacob Faibussowitsch Example Usage: 172530de9b25SBarry Smith .vb 17260a57284eSJacob Faibussowitsch PetscErrorCode foo(int *x) 17270a57284eSJacob Faibussowitsch { 17280a57284eSJacob Faibussowitsch PetscFunctionBegin; // don't forget the begin! 17290a57284eSJacob Faibussowitsch *x = 10; 17303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173130de9b25SBarry Smith } 173230de9b25SBarry Smith .ve 173330de9b25SBarry Smith 17340a57284eSJacob Faibussowitsch May return any arbitrary type\: 17350a57284eSJacob Faibussowitsch .vb 17360a57284eSJacob Faibussowitsch struct Foo 17370a57284eSJacob Faibussowitsch { 17380a57284eSJacob Faibussowitsch int x; 17390a57284eSJacob Faibussowitsch }; 17400a57284eSJacob Faibussowitsch 17410a57284eSJacob Faibussowitsch struct Foo make_foo(int value) 17420a57284eSJacob Faibussowitsch { 17430a57284eSJacob Faibussowitsch struct Foo f; 17440a57284eSJacob Faibussowitsch 17450a57284eSJacob Faibussowitsch PetscFunctionBegin; 17460a57284eSJacob Faibussowitsch f.x = value; 17470a57284eSJacob Faibussowitsch PetscFunctionReturn(f); 17480a57284eSJacob Faibussowitsch } 17490a57284eSJacob Faibussowitsch .ve 17500a57284eSJacob Faibussowitsch 17510a57284eSJacob Faibussowitsch .seealso: `PetscFunctionBegin`, `PetscFunctionBeginUser`, `PetscFunctionReturnVoid()`, 17520a57284eSJacob Faibussowitsch `PetscStackPopNoCheck()` 175330de9b25SBarry Smith M*/ 17545687f061SJacob Faibussowitsch #define PetscFunctionReturn(...) \ 17559371c9d4SSatish Balay do { \ 1756362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 17575687f061SJacob Faibussowitsch return __VA_ARGS__; \ 175827104ee2SJacob Faibussowitsch } while (0) 1759d64ed03dSBarry Smith 17600a57284eSJacob Faibussowitsch /*MC 17610a57284eSJacob Faibussowitsch PetscFunctionReturnVoid - Like `PetscFunctionReturn()` but returns `void` 17620a57284eSJacob Faibussowitsch 17630a57284eSJacob Faibussowitsch Synopsis: 17640a57284eSJacob Faibussowitsch #include <petscerror.h> 17650a57284eSJacob Faibussowitsch void PetscFunctionReturnVoid() 17660a57284eSJacob Faibussowitsch 17670a57284eSJacob Faibussowitsch Not Collective 17680a57284eSJacob Faibussowitsch 17695687f061SJacob Faibussowitsch Level: beginner 17705687f061SJacob Faibussowitsch 17715687f061SJacob Faibussowitsch Note: 17720a57284eSJacob Faibussowitsch Behaves identically to `PetscFunctionReturn()` except that it returns `void`. That is, this 17735687f061SJacob Faibussowitsch macro culminates with `return`. 17740a57284eSJacob Faibussowitsch 17750a57284eSJacob Faibussowitsch Example Usage: 17760a57284eSJacob Faibussowitsch .vb 17770a57284eSJacob Faibussowitsch void foo() 17780a57284eSJacob Faibussowitsch { 17790a57284eSJacob Faibussowitsch PetscFunctionBegin; // must start with PetscFunctionBegin! 17800a57284eSJacob Faibussowitsch bar(); 17810a57284eSJacob Faibussowitsch baz(); 17820a57284eSJacob Faibussowitsch PetscFunctionReturnVoid(); 17830a57284eSJacob Faibussowitsch } 17840a57284eSJacob Faibussowitsch .ve 17850a57284eSJacob Faibussowitsch 17860a57284eSJacob Faibussowitsch .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, PetscFunctionBeginUser` 17870a57284eSJacob Faibussowitsch M*/ 17889371c9d4SSatish Balay #define PetscFunctionReturnVoid() \ 17899371c9d4SSatish Balay do { \ 1790362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 179127104ee2SJacob Faibussowitsch return; \ 179227104ee2SJacob Faibussowitsch } while (0) 179327104ee2SJacob Faibussowitsch #else /* PETSC_USE_DEBUG */ 179427104ee2SJacob Faibussowitsch #define PetscStackPushNoCheck(funct, petsc_routine, hot) 179537154d10SBarry Smith #define PetscStackUpdateLine 1796792fecdfSBarry Smith #define PetscStackPushExternal(funct) 17973ba16761SJacob Faibussowitsch #define PetscStackPopNoCheck(...) 179827104ee2SJacob Faibussowitsch #define PetscStackClearTop 17993a40ed3dSBarry Smith #define PetscFunctionBegin 18000bdf7c52SPeter Brune #define PetscFunctionBeginUser 1801a2f94806SJed Brown #define PetscFunctionBeginHot 18020a57284eSJacob Faibussowitsch #define PetscFunctionReturn(...) return __VA_ARGS__ 18035665465eSBarry Smith #define PetscFunctionReturnVoid() return 1804812af9f3SBarry Smith #define PetscStackPop CHKMEMQ 1805812af9f3SBarry Smith #define PetscStackPush(f) CHKMEMQ 180627104ee2SJacob Faibussowitsch #endif /* PETSC_USE_DEBUG */ 18073a40ed3dSBarry Smith 18086d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 18093ba16761SJacob Faibussowitsch #define PetscStackCallExternalVoid(...) 18103ba16761SJacob Faibussowitsch template <typename F, typename... Args> 18113ba16761SJacob Faibussowitsch void PetscCallExternal(F, Args...); 18121fc6d13cSAlexander template <typename F, typename... Args> 18131fc6d13cSAlexander void PetscCallExternalAbort(F, Args...); 18146d210af2SJacob Faibussowitsch #else 1815586f9135SBarry Smith /*MC 1816e77caa6dSBarry Smith PetscStackCallExternalVoid - Calls an external library routine or user function after pushing the name of the routine on the stack. 1817eb6b5d47SBarry Smith 1818eb6b5d47SBarry Smith Input Parameters: 1819eb6b5d47SBarry Smith + name - string that gives the name of the function being called 1820586f9135SBarry Smith - routine - actual call to the routine, for example, functionname(a,b) 1821fd3f9acdSBarry Smith 1822586f9135SBarry Smith Level: developer 1823eb6b5d47SBarry Smith 182495bd0b28SBarry Smith Notes: 1825792fecdfSBarry Smith Often one should use `PetscCallExternal()` instead. This routine is intended for external library routines that DO NOT return error codes 1826eb6b5d47SBarry Smith 1827586f9135SBarry Smith In debug mode this also checks the memory for corruption at the end of the function call. 1828586f9135SBarry Smith 182995bd0b28SBarry Smith Certain external packages, such as BLAS/LAPACK may have their own macros, `PetscCallBLAS()` for managing the call, error checking, etc. 1830586f9135SBarry Smith 1831586f9135SBarry Smith Developer Note: 1832586f9135SBarry 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. 1833586f9135SBarry Smith 1834792fecdfSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscCallExternal()`, `PetscCallBLAS()` 1835586f9135SBarry Smith @*/ 18363ba16761SJacob Faibussowitsch #define PetscStackCallExternalVoid(name, ...) \ 18379371c9d4SSatish Balay do { \ 183828770333SStefano Zampini PetscStackPushExternal(name); \ 18393ba16761SJacob Faibussowitsch __VA_ARGS__; \ 18409371c9d4SSatish Balay PetscStackPop; \ 18419371c9d4SSatish Balay } while (0) 1842eb6b5d47SBarry Smith 1843586f9135SBarry Smith /*MC 1844792fecdfSBarry Smith PetscCallExternal - Calls an external library routine that returns an error code after pushing the name of the routine on the stack. 1845fd3f9acdSBarry Smith 1846fd3f9acdSBarry Smith Input Parameters: 1847fd3f9acdSBarry Smith + func - name of the routine 1848586f9135SBarry Smith - args - arguments to the routine 1849586f9135SBarry Smith 1850586f9135SBarry Smith Level: developer 1851fd3f9acdSBarry Smith 185295452b02SPatrick Sanan Notes: 1853e77caa6dSBarry Smith This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not. 1854dbf62e16SBarry Smith 1855586f9135SBarry Smith In debug mode this also checks the memory for corruption at the end of the function call. 1856fd3f9acdSBarry Smith 185787497f52SBarry Smith Assumes the error return code of the function is an integer and that a value of 0 indicates success 185887497f52SBarry Smith 1859586f9135SBarry Smith Developer Note: 1860d5b43468SJose E. Roman This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc. 1861586f9135SBarry Smith 18621fc6d13cSAlexander .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()`, `PetscCallExternalAbort()` 1863586f9135SBarry Smith M*/ 18649371c9d4SSatish Balay #define PetscCallExternal(func, ...) \ 18659371c9d4SSatish Balay do { \ 1866a74df02fSJacob Faibussowitsch PetscStackPush(PetscStringize(func)); \ 18673ba16761SJacob Faibussowitsch int ierr_petsc_call_external_ = func(__VA_ARGS__); \ 18681d4906efSStefano Zampini PetscStackPop; \ 18693ba16761SJacob Faibussowitsch PetscCheck(ierr_petsc_call_external_ == 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in %s(): error code %d", PetscStringize(func), ierr_petsc_call_external_); \ 1870fd3f9acdSBarry Smith } while (0) 18711fc6d13cSAlexander 18721fc6d13cSAlexander /*MC 18731fc6d13cSAlexander PetscCallExternalAbort - Calls an external library routine that returns an error code after pushing the name of the routine on the stack. If the external library function return code indicates an error, this prints the error and aborts 18741fc6d13cSAlexander 18751fc6d13cSAlexander Input Parameters: 18761fc6d13cSAlexander + func - name of the routine 18771fc6d13cSAlexander - args - arguments to the routine 18781fc6d13cSAlexander 18791fc6d13cSAlexander Level: developer 18801fc6d13cSAlexander 18811fc6d13cSAlexander Notes: 18821fc6d13cSAlexander This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not. 18831fc6d13cSAlexander 18841fc6d13cSAlexander In debug mode this also checks the memory for corruption at the end of the function call. 18851fc6d13cSAlexander 18861fc6d13cSAlexander Assumes the error return code of the function is an integer and that a value of 0 indicates success 18871fc6d13cSAlexander 18881fc6d13cSAlexander Developer Note: 18891fc6d13cSAlexander This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc. 18901fc6d13cSAlexander 18911fc6d13cSAlexander .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()`, `PetscCallExternal()` 18921fc6d13cSAlexander M*/ 18931fc6d13cSAlexander #define PetscCallExternalAbort(func, ...) \ 18941fc6d13cSAlexander do { \ 18951fc6d13cSAlexander PetscStackUpdateLine; \ 18961fc6d13cSAlexander int ierr_petsc_call_external_ = func(__VA_ARGS__); \ 18971fc6d13cSAlexander if (PetscUnlikely(ierr_petsc_call_external_ != 0)) { \ 18981fc6d13cSAlexander (void)PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_LIB, PETSC_ERROR_INITIAL, "Error in %s(): error code %d", PetscStringize(func), ierr_petsc_call_external_); \ 18991fc6d13cSAlexander PETSCABORTWITHIERR_Private(PETSC_COMM_SELF, PETSC_ERR_LIB); \ 19001fc6d13cSAlexander } \ 19011fc6d13cSAlexander } while (0) 19026d210af2SJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */ 1903