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 112f388eb8bSPatrick Sanan - 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 137fa190f98SMatthew G. Knepley - 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__); \ 15798921bdaSJacob Faibussowitsch 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 1732c71b3e2SJacob Faibussowitsch - message - Error message in 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 2084be741a6SBarry Smith - message - Error message in 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 2700e6b6b59SJacob Faibussowitsch - message - Error message in 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)) { \ 5273ba16761SJacob Faibussowitsch ierr_petsc_call_void_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_void_, PETSC_ERROR_REPEAT, " "); \ 5283ba16761SJacob Faibussowitsch (void)ierr_petsc_call_void_; \ 5299371c9d4SSatish Balay return; \ 5309371c9d4SSatish Balay } \ 5319566063dSJacob Faibussowitsch } while (0) 5329566063dSJacob Faibussowitsch #endif 5339566063dSJacob Faibussowitsch 5349566063dSJacob Faibussowitsch /*MC 5359566063dSJacob Faibussowitsch CHKERRQ - Checks error code returned from PETSc function 53630de9b25SBarry Smith 53730de9b25SBarry Smith Synopsis: 538aaa7dc30SBarry Smith #include <petscsys.h> 5399566063dSJacob Faibussowitsch void CHKERRQ(PetscErrorCode ierr) 5409566063dSJacob Faibussowitsch 5419566063dSJacob Faibussowitsch Not Collective 5429566063dSJacob Faibussowitsch 5432fe279fdSBarry Smith Input Parameter: 5449566063dSJacob Faibussowitsch . ierr - nonzero error code 5459566063dSJacob Faibussowitsch 5469566063dSJacob Faibussowitsch Level: deprecated 5479566063dSJacob Faibussowitsch 548667f096bSBarry Smith Note: 549667f096bSBarry Smith Deprecated in favor of `PetscCall()`. This routine behaves identically to it. 550667f096bSBarry Smith 551db781477SPatrick Sanan .seealso: `PetscCall()` 5529566063dSJacob Faibussowitsch M*/ 5539566063dSJacob Faibussowitsch #define CHKERRQ(...) PetscCall(__VA_ARGS__) 5549566063dSJacob Faibussowitsch #define CHKERRV(...) PetscCallVoid(__VA_ARGS__) 5559566063dSJacob Faibussowitsch 556a1fd7ae3SBarry Smith PETSC_EXTERN void PetscMPIErrorString(PetscMPIInt, size_t, char *); 557db9cea48SBarry Smith 5589566063dSJacob Faibussowitsch /*MC 5599566063dSJacob Faibussowitsch PetscCallMPI - Checks error code returned from MPI calls, if non-zero it calls the error 560648c30bcSBarry Smith handler and then returns a `PetscErrorCode` 5619566063dSJacob Faibussowitsch 5629566063dSJacob Faibussowitsch Synopsis: 5639566063dSJacob Faibussowitsch #include <petscerror.h> 56449c86fc7SBarry Smith void PetscCallMPI(MPI_Function(args)) 56530de9b25SBarry Smith 566eca87e8dSBarry Smith Not Collective 56730de9b25SBarry Smith 5682fe279fdSBarry Smith Input Parameter: 56949c86fc7SBarry Smith . MPI_Function - an MPI function that returns an MPI error code 57030de9b25SBarry Smith 571667f096bSBarry Smith Level: beginner 572667f096bSBarry Smith 5739566063dSJacob Faibussowitsch Notes: 57487497f52SBarry Smith Always returns the error code `PETSC_ERR_MPI`; the MPI error code and string are embedded in 5759566063dSJacob Faibussowitsch the string error message. Do not use this to call any other routines (for example PETSc 5763ba16761SJacob Faibussowitsch routines), it should only be used for direct MPI calls. The user may configure PETSc with the 5773ba16761SJacob Faibussowitsch `--with-strict-petscerrorcode` option to check this at compile-time, otherwise they must 5789566063dSJacob Faibussowitsch check this themselves. 5799566063dSJacob Faibussowitsch 580aaa8cc7dSPierre Jolivet This routine can only be used in functions returning `PetscErrorCode` themselves. If the 5813ba16761SJacob Faibussowitsch calling function returns a different type, use `PetscCallMPIAbort()` instead. 5823ba16761SJacob Faibussowitsch 5839566063dSJacob Faibussowitsch Example Usage: 5849566063dSJacob Faibussowitsch .vb 5859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(...)); // OK, calling MPI function 5869566063dSJacob Faibussowitsch 5879566063dSJacob Faibussowitsch PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead! 5889566063dSJacob Faibussowitsch .ve 5899566063dSJacob Faibussowitsch 59049c86fc7SBarry Smith Fortran Notes: 59187497f52SBarry Smith The Fortran function from which this is used must declare a variable `PetscErrorCode` ierr and ierr must be 59249c86fc7SBarry Smith the final argument to the MPI function being called. 59349c86fc7SBarry Smith 59449c86fc7SBarry Smith In the main program and in Fortran subroutines that do not have ierr as the final return parameter one 59587497f52SBarry Smith should use `PetscCallMPIA()` 59649c86fc7SBarry Smith 59749c86fc7SBarry Smith Fortran Usage: 59849c86fc7SBarry Smith .vb 59949c86fc7SBarry Smith PetscErrorCode ierr or integer ierr 60049c86fc7SBarry Smith ... 60149c86fc7SBarry Smith PetscCallMPI(MPI_Comm_size(...,ierr)) 60249c86fc7SBarry Smith PetscCallMPIA(MPI_Comm_size(...,ierr)) ! Will abort after calling error handler 60349c86fc7SBarry Smith 60449c86fc7SBarry Smith PetscCallMPI(MPI_Comm_size(...,eflag)) ! ERROR, final argument must be ierr 60549c86fc7SBarry Smith .ve 60649c86fc7SBarry Smith 607db781477SPatrick Sanan .seealso: `SETERRMPI()`, `PetscCall()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`, 6083ba16761SJacob Faibussowitsch `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 609648c30bcSBarry Smith `PetscError()`, `CHKMEMQ`, `PetscCallMPINull()` 610648c30bcSBarry Smith M*/ 611648c30bcSBarry Smith 612648c30bcSBarry Smith /*MC 613648c30bcSBarry Smith PetscCallMPINull - Checks error code returned from MPI calls, if non-zero it calls the error 614648c30bcSBarry Smith handler and then returns a `NULL` 615648c30bcSBarry Smith 616648c30bcSBarry Smith Synopsis: 617648c30bcSBarry Smith #include <petscerror.h> 618648c30bcSBarry Smith void PetscCallMPINull(MPI_Function(args)) 619648c30bcSBarry Smith 620648c30bcSBarry Smith Not Collective; No Fortran Support 621648c30bcSBarry Smith 622648c30bcSBarry Smith Input Parameter: 623648c30bcSBarry Smith . MPI_Function - an MPI function that returns an MPI error code 624648c30bcSBarry Smith 625648c30bcSBarry Smith Level: beginner 626648c30bcSBarry Smith 627648c30bcSBarry Smith Notes: 628648c30bcSBarry Smith Always passes the error code `PETSC_ERR_MPI` to the error handler `PetscError()`; the MPI error code and string are embedded in 629648c30bcSBarry Smith the string error message. Do not use this to call any other routines (for example PETSc 630648c30bcSBarry Smith routines), it should only be used for direct MPI calls. 631648c30bcSBarry Smith 632648c30bcSBarry Smith This routine can only be used in functions returning anything that can be returned as a `NULL` themselves. If the 633648c30bcSBarry Smith calling function returns a different type, use `PetscCallMPIAbort()` instead. 634648c30bcSBarry Smith 635648c30bcSBarry Smith Example Usage: 636648c30bcSBarry Smith .vb 637648c30bcSBarry Smith PetscCallMPINull(MPI_Comm_size(...)); // OK, calling MPI function 638648c30bcSBarry Smith 639648c30bcSBarry Smith PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead! 640648c30bcSBarry Smith .ve 641648c30bcSBarry Smith 642648c30bcSBarry Smith .seealso: `SETERRMPI()`, `PetscCall()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`, 643648c30bcSBarry Smith `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 644648c30bcSBarry Smith `PetscError()`, `CHKMEMQ`, `PetscCallMPI()` 6453ba16761SJacob Faibussowitsch M*/ 6463ba16761SJacob Faibussowitsch 6473ba16761SJacob Faibussowitsch /*MC 6483ba16761SJacob Faibussowitsch PetscCallMPIAbort - Like `PetscCallMPI()` but calls `MPI_Abort()` on error 6493ba16761SJacob Faibussowitsch 6503ba16761SJacob Faibussowitsch Synopsis: 6513ba16761SJacob Faibussowitsch #include <petscerror.h> 6523ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm comm, MPI_Function(args)) 6533ba16761SJacob Faibussowitsch 6543ba16761SJacob Faibussowitsch Not Collective 6553ba16761SJacob Faibussowitsch 6563ba16761SJacob Faibussowitsch Input Parameters: 6573ba16761SJacob Faibussowitsch + comm - the MPI communicator to abort on 6583ba16761SJacob Faibussowitsch - MPI_Function - an MPI function that returns an MPI error code 6593ba16761SJacob Faibussowitsch 660667f096bSBarry Smith Level: beginner 661667f096bSBarry Smith 6623ba16761SJacob Faibussowitsch Notes: 6633ba16761SJacob Faibussowitsch Usage is identical to `PetscCallMPI()`. See `PetscCallMPI()` for detailed discussion. 6643ba16761SJacob Faibussowitsch 6653ba16761SJacob Faibussowitsch This routine may be used in functions returning `void` or other non-`PetscErrorCode` types. 6663ba16761SJacob Faibussowitsch 667989c49a2SBarry Smith Fortran Note: 668989c49a2SBarry Smith In Fortran this is called `PetscCallMPIA()` and is intended to be used in the main program while `PetscCallMPI()` is 669989c49a2SBarry Smith used in Fortran subroutines. 670989c49a2SBarry Smith 671989c49a2SBarry Smith Developer Note: 672989c49a2SBarry Smith This should have the same name in Fortran. 673989c49a2SBarry Smith 6743ba16761SJacob Faibussowitsch .seealso: `PetscCallMPI()`, `PetscCallAbort()`, `SETERRABORT()` 67530de9b25SBarry Smith M*/ 6763fcd9f07SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 67763a3b9bcSJacob Faibussowitsch void PetscCallMPI(PetscMPIInt); 6783ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm, PetscMPIInt); 679648c30bcSBarry Smith void PetscCallMPINull(PetscMPIInt); 680064a246eSJacob Faibussowitsch #else 6813ba16761SJacob Faibussowitsch #define PetscCallMPI_Private(__PETSC_STACK_POP_FUNC__, __SETERR_FUNC__, __COMM__, ...) \ 6829371c9d4SSatish Balay do { \ 6833ba16761SJacob Faibussowitsch PetscMPIInt ierr_petsc_call_mpi_; \ 684ef1023bdSBarry Smith PetscStackUpdateLine; \ 685792fecdfSBarry Smith PetscStackPushExternal("MPI function"); \ 686d71ae5a4SJacob Faibussowitsch { \ 6873ba16761SJacob Faibussowitsch ierr_petsc_call_mpi_ = __VA_ARGS__; \ 688d71ae5a4SJacob Faibussowitsch } \ 6893ba16761SJacob Faibussowitsch __PETSC_STACK_POP_FUNC__; \ 6903ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_mpi_ != MPI_SUCCESS)) { \ 6913ba16761SJacob Faibussowitsch char petsc_mpi_7_errorstring[2 * MPI_MAX_ERROR_STRING]; \ 692a1fd7ae3SBarry Smith PetscMPIErrorString(ierr_petsc_call_mpi_, 2 * MPI_MAX_ERROR_STRING, (char *)petsc_mpi_7_errorstring); \ 693a1fd7ae3SBarry Smith __SETERR_FUNC__(__COMM__, PETSC_ERR_MPI, "MPI error %d %s", ierr_petsc_call_mpi_, petsc_mpi_7_errorstring); \ 6943fcd9f07SJacob Faibussowitsch } \ 6953fcd9f07SJacob Faibussowitsch } while (0) 6963ba16761SJacob Faibussowitsch 6973ba16761SJacob Faibussowitsch #define PetscCallMPI(...) PetscCallMPI_Private(PetscStackPop, SETERRQ, PETSC_COMM_SELF, __VA_ARGS__) 6983ba16761SJacob Faibussowitsch #define PetscCallMPIAbort(comm, ...) PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRABORT, comm, __VA_ARGS__) 699648c30bcSBarry Smith #define PetscCallMPINull(...) PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRQNULL, PETSC_COMM_SELF, __VA_ARGS__) 7009566063dSJacob Faibussowitsch #endif 701064a246eSJacob Faibussowitsch 7027037b0edSPatrick Sanan /*MC 7039566063dSJacob Faibussowitsch CHKERRMPI - Checks error code returned from MPI calls, if non-zero it calls the error 7049566063dSJacob Faibussowitsch handler and then returns 7059566063dSJacob Faibussowitsch 7069566063dSJacob Faibussowitsch Synopsis: 7079566063dSJacob Faibussowitsch #include <petscerror.h> 7089566063dSJacob Faibussowitsch void CHKERRMPI(PetscErrorCode ierr) 7099566063dSJacob Faibussowitsch 7109566063dSJacob Faibussowitsch Not Collective 7119566063dSJacob Faibussowitsch 7129566063dSJacob Faibussowitsch Input Parameter: 7139566063dSJacob Faibussowitsch . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 7149566063dSJacob Faibussowitsch 7159566063dSJacob Faibussowitsch Level: deprecated 7169566063dSJacob Faibussowitsch 717667f096bSBarry Smith Note: 718667f096bSBarry Smith Deprecated in favor of `PetscCallMPI()`. This routine behaves identically to it. 719667f096bSBarry Smith 720db781477SPatrick Sanan .seealso: `PetscCallMPI()` 7219566063dSJacob Faibussowitsch M*/ 7229566063dSJacob Faibussowitsch #define CHKERRMPI(...) PetscCallMPI(__VA_ARGS__) 7239566063dSJacob Faibussowitsch 7249566063dSJacob Faibussowitsch /*MC 725989c49a2SBarry Smith PetscCallAbort - Checks error code returned from PETSc function, if non-zero it aborts immediately by calling `MPI_Abort()` 7269566063dSJacob Faibussowitsch 7279566063dSJacob Faibussowitsch Synopsis: 7289566063dSJacob Faibussowitsch #include <petscerror.h> 7299566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm comm, PetscErrorCode ierr) 7309566063dSJacob Faibussowitsch 731c3339decSBarry Smith Collective 7329566063dSJacob Faibussowitsch 7339566063dSJacob Faibussowitsch Input Parameters: 7349566063dSJacob Faibussowitsch + comm - the MPI communicator on which to abort 7359566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 7369566063dSJacob Faibussowitsch 737667f096bSBarry Smith Level: intermediate 738667f096bSBarry Smith 7399566063dSJacob Faibussowitsch Notes: 74087497f52SBarry Smith This macro has identical type and usage semantics to `PetscCall()` with the important caveat 7419566063dSJacob Faibussowitsch that this macro does not return. Instead, if ierr is nonzero it calls the PETSc error handler 74287497f52SBarry Smith and then immediately calls `MPI_Abort()`. It can therefore be used anywhere. 7439566063dSJacob Faibussowitsch 74487497f52SBarry Smith As per `MPI_Abort()` semantics the communicator passed must be valid, although there is currently 74587497f52SBarry Smith no attempt made at handling any potential errors from `MPI_Abort()`. Note that while 74687497f52SBarry Smith `MPI_Abort()` is required to terminate only those processes which reside on comm, it is often 74787497f52SBarry Smith the case that `MPI_Abort()` terminates *all* processes. 7489566063dSJacob Faibussowitsch 7499566063dSJacob Faibussowitsch Example Usage: 7509566063dSJacob Faibussowitsch .vb 7519566063dSJacob Faibussowitsch PetscErrorCode boom(void) { return PETSC_ERR_MEM; } 7529566063dSJacob Faibussowitsch 7539566063dSJacob Faibussowitsch void foo(void) 7549566063dSJacob Faibussowitsch { 7559566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type 7569566063dSJacob Faibussowitsch } 7579566063dSJacob Faibussowitsch 7589566063dSJacob Faibussowitsch double bar(void) 7599566063dSJacob Faibussowitsch { 7609566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type 7619566063dSJacob Faibussowitsch } 7629566063dSJacob Faibussowitsch 7639566063dSJacob Faibussowitsch PetscCallAbort(MPI_COMM_NULL,boom()); // ERROR, communicator should be valid 7649566063dSJacob Faibussowitsch 7659566063dSJacob Faibussowitsch struct baz 7669566063dSJacob Faibussowitsch { 7679566063dSJacob Faibussowitsch baz() 7689566063dSJacob Faibussowitsch { 7699566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,boom()); // OK 7709566063dSJacob Faibussowitsch } 7719566063dSJacob Faibussowitsch 7729566063dSJacob Faibussowitsch ~baz() 7739566063dSJacob Faibussowitsch { 7749566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,boom()); // OK (in fact the only way to handle PETSc errors) 7759566063dSJacob Faibussowitsch } 7769566063dSJacob Faibussowitsch }; 7779566063dSJacob Faibussowitsch .ve 7789566063dSJacob Faibussowitsch 779989c49a2SBarry Smith Fortran Note: 780989c49a2SBarry Smith Use `PetscCallA()`. 781989c49a2SBarry Smith 782989c49a2SBarry Smith Developer Note: 783989c49a2SBarry Smith This should have the same name in Fortran as in C. 784989c49a2SBarry Smith 785db781477SPatrick Sanan .seealso: `SETERRABORT()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, 7865687f061SJacob Faibussowitsch `SETERRQ()`, `CHKMEMQ`, `PetscCallMPI()`, `PetscCallCXXAbort()` 7879566063dSJacob Faibussowitsch M*/ 7889566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 7899566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm, PetscErrorCode); 7909566063dSJacob Faibussowitsch void PetscCallContinue(PetscErrorCode); 7919566063dSJacob Faibussowitsch #else 7929371c9d4SSatish Balay #define PetscCallAbort(comm, ...) \ 7939371c9d4SSatish Balay do { \ 7947da6b3ccSLisandro Dalcin PetscErrorCode ierr_petsc_call_abort_; \ 7953ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 7967da6b3ccSLisandro Dalcin ierr_petsc_call_abort_ = __VA_ARGS__; \ 7973ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_abort_ != PETSC_SUCCESS)) { \ 7983ba16761SJacob Faibussowitsch ierr_petsc_call_abort_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_abort_, PETSC_ERROR_REPEAT, " "); \ 7993ba16761SJacob Faibussowitsch (void)MPI_Abort(comm, (PetscMPIInt)ierr_petsc_call_abort_); \ 8009566063dSJacob Faibussowitsch } \ 8019566063dSJacob Faibussowitsch } while (0) 8029371c9d4SSatish Balay #define PetscCallContinue(...) \ 8039371c9d4SSatish Balay do { \ 8047da6b3ccSLisandro Dalcin PetscErrorCode ierr_petsc_call_continue_; \ 8053ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 8067da6b3ccSLisandro Dalcin ierr_petsc_call_continue_ = __VA_ARGS__; \ 8073ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_continue_ != PETSC_SUCCESS)) { \ 8083ba16761SJacob Faibussowitsch ierr_petsc_call_continue_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_continue_, PETSC_ERROR_REPEAT, " "); \ 8093ba16761SJacob Faibussowitsch (void)ierr_petsc_call_continue_; \ 8103ba16761SJacob Faibussowitsch } \ 8119566063dSJacob Faibussowitsch } while (0) 8129566063dSJacob Faibussowitsch #endif 8139566063dSJacob Faibussowitsch 8149566063dSJacob Faibussowitsch /*MC 8159566063dSJacob Faibussowitsch CHKERRABORT - Checks error code returned from PETSc function. If non-zero it aborts immediately. 8169566063dSJacob Faibussowitsch 8179566063dSJacob Faibussowitsch Synopsis: 8189566063dSJacob Faibussowitsch #include <petscerror.h> 8199566063dSJacob Faibussowitsch void CHKERRABORT(MPI_Comm comm, PetscErrorCode ierr) 8209566063dSJacob Faibussowitsch 8219566063dSJacob Faibussowitsch Not Collective 8229566063dSJacob Faibussowitsch 8239566063dSJacob Faibussowitsch Input Parameters: 8249566063dSJacob Faibussowitsch + comm - the MPI communicator 8259566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 8269566063dSJacob Faibussowitsch 8279566063dSJacob Faibussowitsch Level: deprecated 8289566063dSJacob Faibussowitsch 829667f096bSBarry Smith Note: 830667f096bSBarry Smith Deprecated in favor of `PetscCallAbort()`. This routine behaves identically to it. 831667f096bSBarry Smith 832af27ebaaSBarry Smith .seealso: `PetscCallAbort()`, `PetscErrorCode` 8339566063dSJacob Faibussowitsch M*/ 8349566063dSJacob Faibussowitsch #define CHKERRABORT(comm, ...) PetscCallAbort(comm, __VA_ARGS__) 8359566063dSJacob Faibussowitsch #define CHKERRCONTINUE(...) PetscCallContinue(__VA_ARGS__) 8369566063dSJacob Faibussowitsch 8379566063dSJacob Faibussowitsch /*MC 83887497f52SBarry Smith CHKERRA - Fortran-only replacement for use of `CHKERRQ()` in the main program, which aborts immediately 839f388eb8bSPatrick Sanan 840f388eb8bSPatrick Sanan Synopsis: 841f388eb8bSPatrick Sanan #include <petscsys.h> 842f388eb8bSPatrick Sanan PetscErrorCode CHKERRA(PetscErrorCode ierr) 843f388eb8bSPatrick Sanan 844f388eb8bSPatrick Sanan Not Collective 845f388eb8bSPatrick Sanan 8462fe279fdSBarry Smith Input Parameter: 847f388eb8bSPatrick Sanan . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 848f388eb8bSPatrick Sanan 84987497f52SBarry Smith Level: deprecated 850f388eb8bSPatrick Sanan 85187497f52SBarry Smith Note: 85287497f52SBarry Smith This macro is rarely needed, normal usage is `PetscCallA()` in the main Fortran program. 853f388eb8bSPatrick Sanan 854989c49a2SBarry Smith Developer Note: 855989c49a2SBarry Smith Why isn't this named `CHKERRABORT()` in Fortran? 856989c49a2SBarry Smith 85787497f52SBarry Smith .seealso: `PetscCall()`, `PetscCallA()`, `PetscCallAbort()`, `CHKERRQ()`, `SETERRA()`, `SETERRQ()`, `SETERRABORT()` 858f388eb8bSPatrick Sanan M*/ 859f388eb8bSPatrick Sanan 8601e4f893aSSatish Balay PETSC_EXTERN PetscBool petscwaitonerrorflg; 8611e4f893aSSatish Balay PETSC_EXTERN PetscBool petscindebugger; 86235f00c14SToby Isaac PETSC_EXTERN PetscBool petscabortmpifinalize; 8637c66cc67SJunchao Zhang 864*1fc6d13cSAlexander #if defined(PETSC_CLANG_STATIC_ANALYZER) 865*1fc6d13cSAlexander void PETSCABORTWITHERR_Private(MPI_Comm, PetscErrorCode); 866*1fc6d13cSAlexander #else 867*1fc6d13cSAlexander #define PETSCABORTWITHIERR_Private(comm, ierr) \ 868*1fc6d13cSAlexander do { \ 869*1fc6d13cSAlexander PetscMPIInt size_; \ 870*1fc6d13cSAlexander MPI_Comm_size(comm, &size_); \ 871*1fc6d13cSAlexander if (PetscCIEnabledPortableErrorOutput && (size_ == PetscGlobalSize || petscabortmpifinalize) && ierr != PETSC_ERR_SIG) { \ 872*1fc6d13cSAlexander MPI_Finalize(); \ 873*1fc6d13cSAlexander exit(0); \ 874*1fc6d13cSAlexander } else if (PetscCIEnabledPortableErrorOutput && PetscGlobalSize == 1) { \ 875*1fc6d13cSAlexander exit(0); \ 876*1fc6d13cSAlexander } else { \ 877*1fc6d13cSAlexander MPI_Abort(comm, (PetscMPIInt)ierr); \ 878*1fc6d13cSAlexander } \ 879*1fc6d13cSAlexander } while (0) 880*1fc6d13cSAlexander #endif 881*1fc6d13cSAlexander 8827c66cc67SJunchao Zhang /*MC 883667f096bSBarry Smith PETSCABORT - Call `MPI_Abort()` with an informative error code 8847c66cc67SJunchao Zhang 8857c66cc67SJunchao Zhang Synopsis: 8867c66cc67SJunchao Zhang #include <petscsys.h> 8877c66cc67SJunchao Zhang PETSCABORT(MPI_Comm comm, PetscErrorCode ierr) 8887c66cc67SJunchao Zhang 8897cdbe19fSJose E. Roman Collective; No Fortran Support 8907c66cc67SJunchao Zhang 8917c66cc67SJunchao Zhang Input Parameters: 89295bd0b28SBarry Smith + comm - An MPI communicator, so that the error can be collective 8937c66cc67SJunchao Zhang - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 8947c66cc67SJunchao Zhang 895bf4d2887SBarry Smith Level: advanced 8967c66cc67SJunchao Zhang 897bf4d2887SBarry Smith Notes: 898989c49a2SBarry Smith If the option `-start_in_debugger` was used then this calls `abort()` to stop the program in the debugger. 899bf4d2887SBarry Smith 900989c49a2SBarry Smith if `PetscCIEnabledPortableErrorOutput` is set, which means the code is running in the PETSc test harness (make test), 901989c49a2SBarry Smith and `comm` is `MPI_COMM_WORLD` it strives to exit cleanly without calling `MPI_Abort()` and instead calling `MPI_Finalize()`. 902660278c0SBarry Smith 903989c49a2SBarry Smith This is currently only used when an error propagates up to the C `main()` program and is detected by a `PetscCall()`, `PetscCallMPI()`, 904989c49a2SBarry Smith or is set in `main()` with `SETERRQ()`. Abort calls such as `SETERRABORT()`, 905989c49a2SBarry Smith `PetscCheckAbort()`, `PetscCallMPIAbort()`, and `PetscCallAbort()` always call `MPI_Abort()` and do not have any special 906989c49a2SBarry Smith handling for the test harness. 907989c49a2SBarry Smith 908989c49a2SBarry Smith Developer Note: 909989c49a2SBarry Smith Should the other abort calls also pass through this call instead of calling `MPI_Abort()` directly? 910989c49a2SBarry Smith 911989c49a2SBarry Smith .seealso: `PetscError()`, `PetscCall()`, `SETERRABORT()`, `PetscCheckAbort()`, `PetscCallMPIAbort()`, `PetscCall()`, `PetscCallMPI()`, 912af27ebaaSBarry Smith `PetscCallAbort()`, `MPI_Abort()`, `PetscErrorCode` 9137c66cc67SJunchao Zhang M*/ 91429f88feaSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 91529f88feaSJacob Faibussowitsch void PETSCABORT(MPI_Comm, PetscErrorCode); 91629f88feaSJacob Faibussowitsch #else 9179371c9d4SSatish Balay #define PETSCABORT(comm, ...) \ 9189371c9d4SSatish Balay do { \ 9193ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_abort_; \ 9203ba16761SJacob Faibussowitsch if (petscwaitonerrorflg) { ierr_petsc_abort_ = PetscSleep(1000); } \ 9213ba16761SJacob Faibussowitsch if (petscindebugger) { \ 9223ba16761SJacob Faibussowitsch abort(); \ 9233ba16761SJacob Faibussowitsch } else { \ 9243ba16761SJacob Faibussowitsch ierr_petsc_abort_ = __VA_ARGS__; \ 925*1fc6d13cSAlexander PETSCABORTWITHIERR_Private(comm, ierr_petsc_abort_); \ 9263fcd9f07SJacob Faibussowitsch } \ 9277c66cc67SJunchao Zhang } while (0) 92829f88feaSJacob Faibussowitsch #endif 929986eef2eSBarry Smith 9309566063dSJacob Faibussowitsch #ifdef PETSC_CLANGUAGE_CXX 931986eef2eSBarry Smith /*MC 9329566063dSJacob Faibussowitsch PetscCallThrow - Checks error code, if non-zero it calls the C++ error handler which throws 9339566063dSJacob Faibussowitsch an exception 934986eef2eSBarry Smith 935986eef2eSBarry Smith Synopsis: 9369566063dSJacob Faibussowitsch #include <petscerror.h> 9379566063dSJacob Faibussowitsch void PetscCallThrow(PetscErrorCode ierr) 938986eef2eSBarry Smith 939986eef2eSBarry Smith Not Collective 940986eef2eSBarry Smith 9419566063dSJacob Faibussowitsch Input Parameter: 942986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 943986eef2eSBarry Smith 944667f096bSBarry Smith Level: beginner 945667f096bSBarry Smith 946986eef2eSBarry Smith Notes: 947af27ebaaSBarry Smith Requires PETSc to be configured with clanguage of c++. Throws a std::runtime_error() on error. 948986eef2eSBarry Smith 94987497f52SBarry Smith Once the error handler throws the exception you can use `PetscCallVoid()` which returns without 95087497f52SBarry Smith an error code (bad idea since the error is ignored) or `PetscCallAbort()` to have `MPI_Abort()` 9519566063dSJacob Faibussowitsch called immediately. 9529566063dSJacob Faibussowitsch 953db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscCall()`, `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, 954db781477SPatrick Sanan `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ` 955986eef2eSBarry Smith M*/ 9569371c9d4SSatish Balay #define PetscCallThrow(...) \ 9579371c9d4SSatish Balay do { \ 9583ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 9593ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_throw_ = __VA_ARGS__; \ 9603ba16761SJacob 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); \ 961ffc4695bSBarry Smith } while (0) 96285614651SBarry Smith 963cc26af49SMatthew Knepley /*MC 964cc26af49SMatthew Knepley CHKERRXX - Checks error code, if non-zero it calls the C++ error handler which throws an exception 965cc26af49SMatthew Knepley 966cc26af49SMatthew Knepley Synopsis: 9679566063dSJacob Faibussowitsch #include <petscerror.h> 9683af045c5SBarry Smith void CHKERRXX(PetscErrorCode ierr) 969cc26af49SMatthew Knepley 970eca87e8dSBarry Smith Not Collective 971cc26af49SMatthew Knepley 9729566063dSJacob Faibussowitsch Input Parameter: 9733af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 974cc26af49SMatthew Knepley 9759566063dSJacob Faibussowitsch Level: deprecated 976cc26af49SMatthew Knepley 977667f096bSBarry Smith Note: 978667f096bSBarry Smith Deprecated in favor of `PetscCallThrow()`. This routine behaves identically to it. 979667f096bSBarry Smith 980db781477SPatrick Sanan .seealso: `PetscCallThrow()` 981cc26af49SMatthew Knepley M*/ 9829566063dSJacob Faibussowitsch #define CHKERRXX(...) PetscCallThrow(__VA_ARGS__) 983fd705b32SMatthew Knepley #endif 984fd705b32SMatthew Knepley 9853ba16761SJacob Faibussowitsch #define PetscCallCXX_Private(__SETERR_FUNC__, __COMM__, ...) \ 9863ba16761SJacob Faibussowitsch do { \ 9873ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 9883ba16761SJacob Faibussowitsch try { \ 9893ba16761SJacob Faibussowitsch __VA_ARGS__; \ 9903ba16761SJacob Faibussowitsch } catch (const std::exception &e) { \ 9913ba16761SJacob Faibussowitsch __SETERR_FUNC__(__COMM__, PETSC_ERR_LIB, "%s", e.what()); \ 9923ba16761SJacob Faibussowitsch } \ 9933ba16761SJacob Faibussowitsch } while (0) 9943ba16761SJacob Faibussowitsch 9953f520e80SJunchao Zhang /*MC 9969566063dSJacob Faibussowitsch PetscCallCXX - Checks C++ function calls and if they throw an exception, catch it and then 9979566063dSJacob Faibussowitsch return a PETSc error code 9983f520e80SJunchao Zhang 9993f520e80SJunchao Zhang Synopsis: 10009566063dSJacob Faibussowitsch #include <petscerror.h> 10015687f061SJacob Faibussowitsch void PetscCallCXX(...) noexcept; 10023f520e80SJunchao Zhang 10033f520e80SJunchao Zhang Not Collective 10043f520e80SJunchao Zhang 10059566063dSJacob Faibussowitsch Input Parameter: 10065687f061SJacob Faibussowitsch . __VA_ARGS__ - An arbitrary expression 10075687f061SJacob Faibussowitsch 10085687f061SJacob Faibussowitsch Level: beginner 10093f520e80SJunchao Zhang 10103f520e80SJunchao Zhang Notes: 10115687f061SJacob Faibussowitsch `PetscCallCXX(...)` is a macro replacement for 10129566063dSJacob Faibussowitsch .vb 10139566063dSJacob Faibussowitsch try { 10145687f061SJacob Faibussowitsch __VA_ARGS__; 10159566063dSJacob Faibussowitsch } catch (const std::exception& e) { 10169566063dSJacob Faibussowitsch return ConvertToPetscErrorCode(e); 10179566063dSJacob Faibussowitsch } 10189566063dSJacob Faibussowitsch .ve 10199566063dSJacob Faibussowitsch Due to the fact that it catches any (reasonable) exception, it is essentially noexcept. 10203f520e80SJunchao Zhang 10215687f061SJacob Faibussowitsch If you cannot return a `PetscErrorCode` use `PetscCallCXXAbort()` instead. 10225687f061SJacob Faibussowitsch 10239566063dSJacob Faibussowitsch Example Usage: 10249566063dSJacob Faibussowitsch .vb 10259566063dSJacob Faibussowitsch void foo(void) { throw std::runtime_error("error"); } 10263f520e80SJunchao Zhang 10279566063dSJacob Faibussowitsch void bar() 10289566063dSJacob Faibussowitsch { 1029e8952933SJacob Faibussowitsch PetscCallCXX(foo()); // ERROR bar() does not return PetscErrorCode 10309566063dSJacob Faibussowitsch } 10319566063dSJacob Faibussowitsch 10329566063dSJacob Faibussowitsch PetscErrorCode baz() 10339566063dSJacob Faibussowitsch { 10349566063dSJacob Faibussowitsch PetscCallCXX(foo()); // OK 10359566063dSJacob Faibussowitsch 10369566063dSJacob Faibussowitsch PetscCallCXX( 10379566063dSJacob Faibussowitsch bar(); 1038d5b43468SJose E. Roman foo(); // OK multiple statements allowed 10399566063dSJacob Faibussowitsch ); 10409566063dSJacob Faibussowitsch } 10419566063dSJacob Faibussowitsch 10429566063dSJacob Faibussowitsch struct bop 10439566063dSJacob Faibussowitsch { 10449566063dSJacob Faibussowitsch bop() 10459566063dSJacob Faibussowitsch { 1046e8952933SJacob Faibussowitsch PetscCallCXX(foo()); // ERROR returns PetscErrorCode, cannot be used in constructors 10479566063dSJacob Faibussowitsch } 10489566063dSJacob Faibussowitsch }; 10499566063dSJacob Faibussowitsch 1050e8952933SJacob Faibussowitsch // ERROR contains do-while, cannot be used as function-try block 10519566063dSJacob Faibussowitsch PetscErrorCode qux() PetscCallCXX( 10529566063dSJacob Faibussowitsch bar(); 10539566063dSJacob Faibussowitsch baz(); 10549566063dSJacob Faibussowitsch foo(); 10559566063dSJacob Faibussowitsch return 0; 10569566063dSJacob Faibussowitsch ) 10579566063dSJacob Faibussowitsch .ve 10589566063dSJacob Faibussowitsch 10595687f061SJacob Faibussowitsch .seealso: `PetscCallCXXAbort()`, `PetscCallThrow()`, `SETERRQ()`, `PetscCall()`, 10605687f061SJacob Faibussowitsch `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 10615687f061SJacob Faibussowitsch `PetscError()`, `CHKMEMQ` 10623f520e80SJunchao Zhang M*/ 10633ba16761SJacob Faibussowitsch #define PetscCallCXX(...) PetscCallCXX_Private(SETERRQ, PETSC_COMM_SELF, __VA_ARGS__) 10645687f061SJacob Faibussowitsch 10655687f061SJacob Faibussowitsch /*MC 10665687f061SJacob Faibussowitsch PetscCallCXXAbort - Like `PetscCallCXX()` but calls `MPI_Abort()` instead of returning an 10675687f061SJacob Faibussowitsch error-code 10685687f061SJacob Faibussowitsch 10695687f061SJacob Faibussowitsch Synopsis: 10705687f061SJacob Faibussowitsch #include <petscerror.h> 10715687f061SJacob Faibussowitsch void PetscCallCXXAbort(MPI_Comm comm, ...) noexcept; 10725687f061SJacob Faibussowitsch 107320f4b53cSBarry Smith Collective; No Fortran Support 10745687f061SJacob Faibussowitsch 10755687f061SJacob Faibussowitsch Input Parameters: 10765687f061SJacob Faibussowitsch + comm - The MPI communicator to abort on 10775687f061SJacob Faibussowitsch - __VA_ARGS__ - An arbitrary expression 10785687f061SJacob Faibussowitsch 10795687f061SJacob Faibussowitsch Level: beginner 10805687f061SJacob Faibussowitsch 10815687f061SJacob Faibussowitsch Notes: 10825687f061SJacob Faibussowitsch This macro may be used to check C++ expressions for exceptions in cases where you cannot 10835687f061SJacob Faibussowitsch return an error code. This includes constructors, destructors, copy/move assignment functions 10845687f061SJacob Faibussowitsch or constructors among others. 10855687f061SJacob Faibussowitsch 10865687f061SJacob Faibussowitsch If an exception is caught, the macro calls `SETERRABORT()` on `comm`. The exception must 10875687f061SJacob Faibussowitsch derive from `std::exception` in order to be caught. 10885687f061SJacob Faibussowitsch 10895687f061SJacob Faibussowitsch If the routine _can_ return an error-code it is highly advised to use `PetscCallCXX()` 10905687f061SJacob Faibussowitsch instead. 10915687f061SJacob Faibussowitsch 10925687f061SJacob Faibussowitsch See `PetscCallCXX()` for additional discussion. 10935687f061SJacob Faibussowitsch 10945687f061SJacob Faibussowitsch Example Usage: 10955687f061SJacob Faibussowitsch .vb 10965687f061SJacob Faibussowitsch class Foo 10975687f061SJacob Faibussowitsch { 10985687f061SJacob Faibussowitsch std::vector<int> data_; 10995687f061SJacob Faibussowitsch 11005687f061SJacob Faibussowitsch public: 11015687f061SJacob Faibussowitsch // normally std::vector::reserve() may raise an exception, but since we handle it with 11025687f061SJacob Faibussowitsch // PetscCallCXXAbort() we may mark this routine as noexcept! 11035687f061SJacob Faibussowitsch Foo() noexcept 11045687f061SJacob Faibussowitsch { 11055687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, data_.reserve(10)); 11065687f061SJacob Faibussowitsch } 11075687f061SJacob Faibussowitsch }; 11085687f061SJacob Faibussowitsch 11095687f061SJacob Faibussowitsch std::vector<int> bar() 11105687f061SJacob Faibussowitsch { 11115687f061SJacob Faibussowitsch std::vector<int> v; 11125687f061SJacob Faibussowitsch 11135687f061SJacob Faibussowitsch PetscFunctionBegin; 11145687f061SJacob Faibussowitsch // OK! 11155687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1)); 11165687f061SJacob Faibussowitsch PetscFunctionReturn(v); 11175687f061SJacob Faibussowitsch } 11185687f061SJacob Faibussowitsch 11195687f061SJacob Faibussowitsch PetscErrorCode baz() 11205687f061SJacob Faibussowitsch { 11215687f061SJacob Faibussowitsch std::vector<int> v; 11225687f061SJacob Faibussowitsch 11235687f061SJacob Faibussowitsch PetscFunctionBegin; 11245687f061SJacob Faibussowitsch // WRONG! baz() returns a PetscErrorCode, prefer PetscCallCXX() instead 11255687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1)); 11263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11275687f061SJacob Faibussowitsch } 11285687f061SJacob Faibussowitsch .ve 11295687f061SJacob Faibussowitsch 11305687f061SJacob Faibussowitsch .seealso: `PetscCallCXX()`, `SETERRABORT()`, `PetscCallAbort()` 11315687f061SJacob Faibussowitsch M*/ 11323ba16761SJacob Faibussowitsch #define PetscCallCXXAbort(comm, ...) PetscCallCXX_Private(SETERRABORT, comm, __VA_ARGS__) 11333f520e80SJunchao Zhang 113430de9b25SBarry Smith /*MC 11359566063dSJacob Faibussowitsch CHKERRCXX - Checks C++ function calls and if they throw an exception, catch it and then 11369566063dSJacob Faibussowitsch return a PETSc error code 11379566063dSJacob Faibussowitsch 11389566063dSJacob Faibussowitsch Synopsis: 11399566063dSJacob Faibussowitsch #include <petscerror.h> 11409566063dSJacob Faibussowitsch void CHKERRCXX(func) noexcept; 11419566063dSJacob Faibussowitsch 11429566063dSJacob Faibussowitsch Not Collective 11439566063dSJacob Faibussowitsch 11449566063dSJacob Faibussowitsch Input Parameter: 11459566063dSJacob Faibussowitsch . func - C++ function calls 11469566063dSJacob Faibussowitsch 11479566063dSJacob Faibussowitsch Level: deprecated 11489566063dSJacob Faibussowitsch 1149667f096bSBarry Smith Note: 1150667f096bSBarry Smith Deprecated in favor of `PetscCallCXX()`. This routine behaves identically to it. 1151667f096bSBarry Smith 1152db781477SPatrick Sanan .seealso: `PetscCallCXX()` 11539566063dSJacob Faibussowitsch M*/ 11549566063dSJacob Faibussowitsch #define CHKERRCXX(...) PetscCallCXX(__VA_ARGS__) 11559566063dSJacob Faibussowitsch 11569566063dSJacob Faibussowitsch /*MC 115730de9b25SBarry Smith CHKMEMQ - Checks the memory for corruption, calls error handler if any is detected 115830de9b25SBarry Smith 115930de9b25SBarry Smith Synopsis: 1160aaa7dc30SBarry Smith #include <petscsys.h> 116191d3bdf4SKris Buschelman CHKMEMQ; 116230de9b25SBarry Smith 1163eca87e8dSBarry Smith Not Collective 1164eca87e8dSBarry Smith 116530de9b25SBarry Smith Level: beginner 116630de9b25SBarry Smith 116730de9b25SBarry Smith Notes: 1168af27ebaaSBarry Smith We recommend using Valgrind <https://petsc.org/release/faq/#valgrind> or for NVIDIA CUDA systems 1169af27ebaaSBarry Smith <https://docs.nvidia.com/cuda/cuda-memcheck/index.html> for finding memory problems. The ``CHKMEMQ`` macro is useful on systems that 11705ed36255SBarry Smith do not have valgrind, but is not as good as valgrind or cuda-memcheck. 11711957e957SBarry Smith 1172667f096bSBarry Smith Must run with the option `-malloc_debug` (`-malloc_test` in debug mode; or if `PetscMallocSetDebug()` called) to enable this option 117330de9b25SBarry Smith 117430de9b25SBarry Smith Once the error handler is called the calling function is then returned from with the given error code. 117530de9b25SBarry Smith 117630de9b25SBarry Smith By defaults prints location where memory that is corrupted was allocated. 117730de9b25SBarry Smith 11788af9f190SBarry Smith Use `CHKMEMA` for functions that return `void` 1179f621e05eSBarry Smith 1180db781477SPatrick Sanan .seealso: `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `SETERRQ()`, `PetscMallocValidate()` 118130de9b25SBarry Smith M*/ 11826d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 1183064a246eSJacob Faibussowitsch #define CHKMEMQ 1184064a246eSJacob Faibussowitsch #define CHKMEMA 11856d210af2SJacob Faibussowitsch #else 11869371c9d4SSatish Balay #define CHKMEMQ \ 11879371c9d4SSatish Balay do { \ 11883ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_memq_ = PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__); \ 11893ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_memq_ != PETSC_SUCCESS)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_memq_, PETSC_ERROR_REPEAT, " "); \ 119086d09637SLisandro Dalcin } while (0) 11916d210af2SJacob Faibussowitsch #define CHKMEMA PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__) 1192064a246eSJacob Faibussowitsch #endif 11939566063dSJacob Faibussowitsch 1194668f157eSBarry Smith /*E 1195668f157eSBarry Smith PetscErrorType - passed to the PETSc error handling routines indicating if this is the first or a later call to the error handlers 1196668f157eSBarry Smith 1197668f157eSBarry Smith Level: advanced 1198668f157eSBarry Smith 1199667f096bSBarry Smith Note: 120087497f52SBarry Smith `PETSC_ERROR_IN_CXX` indicates the error was detected in C++ and an exception should be generated 1201d736bfebSBarry Smith 1202667f096bSBarry Smith Developer Note: 1203667f096bSBarry Smith This is currently used to decide when to print the detailed information about the run in `PetscTraceBackErrorHandler()` 1204668f157eSBarry Smith 120587497f52SBarry Smith .seealso: `PetscError()`, `SETERRQ()` 1206668f157eSBarry Smith E*/ 12079371c9d4SSatish Balay typedef enum { 12089371c9d4SSatish Balay PETSC_ERROR_INITIAL = 0, 12099371c9d4SSatish Balay PETSC_ERROR_REPEAT = 1, 12109371c9d4SSatish Balay PETSC_ERROR_IN_CXX = 2 12119371c9d4SSatish Balay } PetscErrorType; 12124b209cf6SBarry Smith 1213eb9e708aSLisandro Dalcin #if defined(__clang_analyzer__) 1214eb9e708aSLisandro Dalcin __attribute__((analyzer_noreturn)) 1215eb9e708aSLisandro Dalcin #endif 12169371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode 12179371c9d4SSatish Balay PetscError(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, ...) PETSC_ATTRIBUTE_COLD PETSC_ATTRIBUTE_FORMAT(7, 8); 1218eb9e708aSLisandro Dalcin 1219014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscErrorPrintfInitialize(void); 12203ba16761SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscErrorMessage(PetscErrorCode, const char *[], char **); 1221d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscTraceBackErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1222d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscIgnoreErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1223d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscEmacsClientErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1224d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscMPIAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1225d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1226d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAttachDebuggerErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1227d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscReturnErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1228efca3c55SSatish Balay PETSC_EXTERN PetscErrorCode PetscPushErrorHandler(PetscErrorCode (*handler)(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *), void *); 1229014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopErrorHandler(void); 12308d359177SBarry Smith PETSC_EXTERN PetscErrorCode PetscSignalHandlerDefault(int, void *); 1231014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPushSignalHandler(PetscErrorCode (*)(int, void *), void *); 1232014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopSignalHandler(void); 123328559dc8SJed Brown PETSC_EXTERN PetscErrorCode PetscCheckPointerSetIntensity(PetscInt); 1234c2a741eeSJunchao Zhang PETSC_EXTERN void PetscSignalSegvCheckPointerOrMpi(void); 1235edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 13, 0, "PetscSignalSegvCheckPointerOrMpi()", ) static inline void PetscSignalSegvCheckPointer(void) 1236d71ae5a4SJacob Faibussowitsch { 12379371c9d4SSatish Balay PetscSignalSegvCheckPointerOrMpi(); 12389371c9d4SSatish Balay } 1239329f5518SBarry Smith 1240639ff905SBarry Smith /*MC 1241639ff905SBarry Smith PetscErrorPrintf - Prints error messages. 1242639ff905SBarry Smith 1243639ff905SBarry Smith Synopsis: 1244aaa7dc30SBarry Smith #include <petscsys.h> 1245639ff905SBarry Smith PetscErrorCode (*PetscErrorPrintf)(const char format[], ...); 1246639ff905SBarry Smith 12477cdbe19fSJose E. Roman Not Collective; No Fortran Support 12487cdbe19fSJose E. Roman 1249f899ff85SJose E. Roman Input Parameter: 1250cd05f99aSJacob Faibussowitsch . format - the usual `printf()` format string 1251639ff905SBarry Smith 1252639ff905SBarry Smith Options Database Keys: 12531957e957SBarry Smith + -error_output_stdout - cause error messages to be printed to stdout instead of the (default) stderr 1254e1bc860dSBarry Smith - -error_output_none - to turn off all printing of error messages (does not change the way the error is handled.) 1255639ff905SBarry Smith 1256cf53795eSBarry Smith Level: developer 1257cf53795eSBarry Smith 125895452b02SPatrick Sanan Notes: 125995452b02SPatrick Sanan Use 1260667f096bSBarry Smith .vb 126195bd0b28SBarry Smith PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the error is handled) and 1262667f096bSBarry Smith PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on or you can use your own function 1263667f096bSBarry Smith .ve 1264639ff905SBarry Smith Use 1265667f096bSBarry Smith .vb 126687497f52SBarry Smith `PETSC_STDERR` = FILE* obtained from a file open etc. to have stderr printed to the file. 126787497f52SBarry Smith `PETSC_STDOUT` = FILE* obtained from a file open etc. to have stdout printed to the file. 1268667f096bSBarry Smith .ve 1269639ff905SBarry Smith Use 1270af27ebaaSBarry Smith .vb 127187497f52SBarry Smith `PetscPushErrorHandler()` to provide your own error handler that determines what kind of messages to print 1272af27ebaaSBarry Smith .ve 1273639ff905SBarry Smith 1274db781477SPatrick Sanan .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscHelpPrintf()`, `PetscPrintf()`, `PetscPushErrorHandler()`, `PetscVFPrintf()`, `PetscHelpPrintf()` 1275639ff905SBarry Smith M*/ 12763ca90d2dSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode (*PetscErrorPrintf)(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2); 1277639ff905SBarry Smith 1278cf0818bdSBarry Smith /*E 1279cf0818bdSBarry Smith PetscFPTrap - types of floating point exceptions that may be trapped 1280cf0818bdSBarry Smith 128187497f52SBarry Smith Currently only `PETSC_FP_TRAP_OFF` and `PETSC_FP_TRAP_ON` are handled. All others are treated as `PETSC_FP_TRAP_ON`. 1282cf0818bdSBarry Smith 1283cf0818bdSBarry Smith Level: intermediate 1284cf0818bdSBarry Smith 12857de69702SBarry Smith .seealso: `PetscSetFPTrap()`, `PetscFPTrapPush()` 1286cf0818bdSBarry Smith E*/ 12879371c9d4SSatish Balay typedef enum { 12889371c9d4SSatish Balay PETSC_FP_TRAP_OFF = 0, 12899371c9d4SSatish Balay PETSC_FP_TRAP_INDIV = 1, 12909371c9d4SSatish Balay PETSC_FP_TRAP_FLTOPERR = 2, 12919371c9d4SSatish Balay PETSC_FP_TRAP_FLTOVF = 4, 12929371c9d4SSatish Balay PETSC_FP_TRAP_FLTUND = 8, 12939371c9d4SSatish Balay PETSC_FP_TRAP_FLTDIV = 16, 12949371c9d4SSatish Balay PETSC_FP_TRAP_FLTINEX = 32 12959371c9d4SSatish Balay } PetscFPTrap; 1296bd2b07b1SBarry Smith #define PETSC_FP_TRAP_ON (PetscFPTrap)(PETSC_FP_TRAP_INDIV | PETSC_FP_TRAP_FLTOPERR | PETSC_FP_TRAP_FLTOVF | PETSC_FP_TRAP_FLTDIV | PETSC_FP_TRAP_FLTINEX) 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) 13446d210af2SJacob Faibussowitsch #define PetscStackPopNoCheck 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) 1353dfb7d7afSStefano Zampini #elif defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY) 1354660278c0SBarry Smith 13559371c9d4SSatish Balay #define PetscStackPush_Private(stack__, file__, func__, line__, petsc_routine__, hot__) \ 13569371c9d4SSatish Balay do { \ 13575a96b57dSJacob Faibussowitsch if (stack__.currentsize < PETSCSTACKSIZE) { \ 13585a96b57dSJacob Faibussowitsch stack__.function[stack__.currentsize] = func__; \ 1359ef1023bdSBarry Smith if (petsc_routine__) { \ 1360ef1023bdSBarry Smith stack__.file[stack__.currentsize] = file__; \ 13615a96b57dSJacob Faibussowitsch stack__.line[stack__.currentsize] = line__; \ 1362ef1023bdSBarry Smith } else { \ 1363648bc8c4SBarry Smith stack__.file[stack__.currentsize] = PETSC_NULLPTR; \ 1364ef1023bdSBarry Smith stack__.line[stack__.currentsize] = 0; \ 1365ef1023bdSBarry Smith } \ 13665a96b57dSJacob Faibussowitsch stack__.petscroutine[stack__.currentsize] = petsc_routine__; \ 13675a96b57dSJacob Faibussowitsch } \ 13685a96b57dSJacob Faibussowitsch ++stack__.currentsize; \ 13695a96b57dSJacob Faibussowitsch stack__.hotdepth += (hot__ || stack__.hotdepth); \ 13705a96b57dSJacob Faibussowitsch } while (0) 13715a96b57dSJacob Faibussowitsch 13724be741a6SBarry Smith /* uses PetscCheckAbort() because may be used in a function that does not return an error code */ 13739371c9d4SSatish Balay #define PetscStackPop_Private(stack__, func__) \ 13749371c9d4SSatish Balay do { \ 13754be741a6SBarry 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__); \ 13765a96b57dSJacob Faibussowitsch if (--stack__.currentsize < PETSCSTACKSIZE) { \ 13779371c9d4SSatish 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", \ 13789371c9d4SSatish Balay stack__.function[stack__.currentsize], stack__.file[stack__.currentsize], stack__.line[stack__.currentsize], func__, __FILE__, __LINE__); \ 13795a96b57dSJacob Faibussowitsch stack__.function[stack__.currentsize] = PETSC_NULLPTR; \ 13805a96b57dSJacob Faibussowitsch stack__.file[stack__.currentsize] = PETSC_NULLPTR; \ 13815a96b57dSJacob Faibussowitsch stack__.line[stack__.currentsize] = 0; \ 13825a96b57dSJacob Faibussowitsch stack__.petscroutine[stack__.currentsize] = 0; \ 13835a96b57dSJacob Faibussowitsch } \ 13845a96b57dSJacob Faibussowitsch stack__.hotdepth = PetscMax(stack__.hotdepth - 1, 0); \ 13855a96b57dSJacob Faibussowitsch } while (0) 13865a96b57dSJacob Faibussowitsch 1387586f9135SBarry Smith /*MC 1388586f9135SBarry Smith PetscStackPushNoCheck - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is 1389586f9135SBarry Smith currently in the source code. 1390586f9135SBarry Smith 1391586f9135SBarry Smith Synopsis: 1392586f9135SBarry Smith #include <petscsys.h> 1393586f9135SBarry Smith void PetscStackPushNoCheck(char *funct,int petsc_routine,PetscBool hot); 1394586f9135SBarry Smith 13957cdbe19fSJose E. Roman Not Collective 13967cdbe19fSJose E. Roman 1397586f9135SBarry Smith Input Parameters: 1398586f9135SBarry Smith + funct - the function name 1399586f9135SBarry Smith . petsc_routine - 2 user function, 1 PETSc function, 0 some other function 1400586f9135SBarry Smith - hot - indicates that the function may be called often so expensive error checking should be turned off inside the function 1401586f9135SBarry Smith 1402586f9135SBarry Smith Level: developer 1403586f9135SBarry Smith 1404586f9135SBarry Smith Notes: 1405586f9135SBarry 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 140687497f52SBarry 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 1407586f9135SBarry Smith help debug the problem. 1408586f9135SBarry Smith 1409ef1023bdSBarry Smith This version does not check the memory corruption (an expensive operation), use `PetscStackPush()` to check the memory. 1410ef1023bdSBarry Smith 1411792fecdfSBarry Smith Use `PetscStackPushExternal()` for a function call that is about to be made to a non-PETSc or user function (such as BLAS etc). 1412ef1023bdSBarry Smith 1413586f9135SBarry Smith The default stack is a global variable called `petscstack`. 1414586f9135SBarry Smith 1415586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1416ef1023bdSBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPush()`, `PetscStackPop`, 1417792fecdfSBarry Smith `PetscStackPushExternal()` 1418586f9135SBarry Smith M*/ 14199371c9d4SSatish Balay #define PetscStackPushNoCheck(funct, petsc_routine, hot) \ 14209371c9d4SSatish Balay do { \ 1421e04113cfSBarry Smith PetscStackSAWsTakeAccess(); \ 14225a96b57dSJacob Faibussowitsch PetscStackPush_Private(petscstack, __FILE__, funct, __LINE__, petsc_routine, hot); \ 1423e04113cfSBarry Smith PetscStackSAWsGrantAccess(); \ 1424441dd030SJed Brown } while (0) 1425441dd030SJed Brown 1426586f9135SBarry Smith /*MC 142787497f52SBarry Smith PetscStackUpdateLine - in a function that has a `PetscFunctionBegin` or `PetscFunctionBeginUser` updates the stack line number to the 142837154d10SBarry Smith current line number. 142937154d10SBarry Smith 143037154d10SBarry Smith Synopsis: 143137154d10SBarry Smith #include <petscsys.h> 143237154d10SBarry Smith void PetscStackUpdateLine 143337154d10SBarry Smith 14347cdbe19fSJose E. Roman Not Collective 14357cdbe19fSJose E. Roman 143637154d10SBarry Smith Level: developer 143737154d10SBarry Smith 143837154d10SBarry Smith Notes: 143987497f52SBarry Smith Using `PetscCall()` and friends automatically handles this process 144087497f52SBarry Smith 144137154d10SBarry 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 144237154d10SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 144337154d10SBarry Smith help debug the problem. 144437154d10SBarry Smith 144595bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 144637154d10SBarry Smith 144737154d10SBarry Smith This is used by `PetscCall()` and is otherwise not like to be needed 144837154d10SBarry Smith 144937154d10SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`, `PetscCall()` 145037154d10SBarry Smith M*/ 145137154d10SBarry Smith #define PetscStackUpdateLine \ 14523ba16761SJacob Faibussowitsch do { \ 1453f1e99121SPierre Jolivet if (petscstack.currentsize > 0 && petscstack.currentsize < PETSCSTACKSIZE && petscstack.function[petscstack.currentsize - 1] == PETSC_FUNCTION_NAME) { petscstack.line[petscstack.currentsize - 1] = __LINE__; } \ 14543ba16761SJacob Faibussowitsch } while (0) 145537154d10SBarry Smith 145637154d10SBarry Smith /*MC 1457792fecdfSBarry Smith PetscStackPushExternal - Pushes a new function name onto the PETSc default stack that tracks where the running program is 1458ef1023bdSBarry Smith currently in the source code. Does not include the filename or line number since this is called by the calling routine 1459ef1023bdSBarry Smith for non-PETSc or user functions. 1460ef1023bdSBarry Smith 1461ef1023bdSBarry Smith Synopsis: 1462ef1023bdSBarry Smith #include <petscsys.h> 1463792fecdfSBarry Smith void PetscStackPushExternal(char *funct); 1464ef1023bdSBarry Smith 14657cdbe19fSJose E. Roman Not Collective 14667cdbe19fSJose E. Roman 14672fe279fdSBarry Smith Input Parameter: 1468ef1023bdSBarry Smith . funct - the function name 1469ef1023bdSBarry Smith 1470ef1023bdSBarry Smith Level: developer 1471ef1023bdSBarry Smith 1472ef1023bdSBarry Smith Notes: 147387497f52SBarry Smith Using `PetscCallExternal()` and friends automatically handles this process 147487497f52SBarry Smith 1475ef1023bdSBarry 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 1476ef1023bdSBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1477ef1023bdSBarry Smith help debug the problem. 1478ef1023bdSBarry Smith 1479ef1023bdSBarry Smith The default stack is a global variable called `petscstack`. 1480ef1023bdSBarry Smith 1481ef1023bdSBarry Smith This is to be used when calling an external package function such as a BLAS function. 1482ef1023bdSBarry Smith 1483ef1023bdSBarry Smith This also updates the stack line number for the current stack function. 1484ef1023bdSBarry Smith 1485ef1023bdSBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1486ef1023bdSBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1487ef1023bdSBarry Smith M*/ 14889371c9d4SSatish Balay #define PetscStackPushExternal(funct) \ 14899371c9d4SSatish Balay do { \ 14909371c9d4SSatish Balay PetscStackUpdateLine; \ 14919371c9d4SSatish Balay PetscStackPushNoCheck(funct, 0, PETSC_TRUE); \ 1492a8f51744SPierre Jolivet } while (0) 1493ef1023bdSBarry Smith 1494ef1023bdSBarry Smith /*MC 1495586f9135SBarry Smith PetscStackPopNoCheck - Pops a function name from the PETSc default stack that tracks where the running program is 1496586f9135SBarry Smith currently in the source code. 1497586f9135SBarry Smith 1498586f9135SBarry Smith Synopsis: 1499586f9135SBarry Smith #include <petscsys.h> 1500586f9135SBarry Smith void PetscStackPopNoCheck(char *funct); 1501586f9135SBarry Smith 15027cdbe19fSJose E. Roman Not Collective 15037cdbe19fSJose E. Roman 1504586f9135SBarry Smith Input Parameter: 1505586f9135SBarry Smith . funct - the function name 1506586f9135SBarry Smith 1507586f9135SBarry Smith Level: developer 1508586f9135SBarry Smith 1509586f9135SBarry Smith Notes: 151087497f52SBarry Smith Using `PetscCall()`, `PetscCallExternal()`, `PetscCallBack()` and friends negates the need to call this 151187497f52SBarry Smith 1512586f9135SBarry 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 1513586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1514586f9135SBarry Smith help debug the problem. 1515586f9135SBarry Smith 151695bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 1517586f9135SBarry Smith 1518586f9135SBarry Smith Developer Note: 1519586f9135SBarry Smith `PetscStackPopNoCheck()` takes a function argument while `PetscStackPop` does not, this difference is likely just historical. 1520586f9135SBarry Smith 1521586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1522586f9135SBarry Smith M*/ 15239371c9d4SSatish Balay #define PetscStackPopNoCheck(funct) \ 15249371c9d4SSatish Balay do { \ 1525362febeeSStefano Zampini PetscStackSAWsTakeAccess(); \ 15265a96b57dSJacob Faibussowitsch PetscStackPop_Private(petscstack, funct); \ 1527362febeeSStefano Zampini PetscStackSAWsGrantAccess(); \ 1528362febeeSStefano Zampini } while (0) 1529362febeeSStefano Zampini 15309371c9d4SSatish Balay #define PetscStackClearTop \ 15319371c9d4SSatish Balay do { \ 1532e04113cfSBarry Smith PetscStackSAWsTakeAccess(); \ 15339371c9d4SSatish Balay if (petscstack.currentsize > 0 && --petscstack.currentsize < PETSCSTACKSIZE) { \ 153427104ee2SJacob Faibussowitsch petscstack.function[petscstack.currentsize] = PETSC_NULLPTR; \ 153527104ee2SJacob Faibussowitsch petscstack.file[petscstack.currentsize] = PETSC_NULLPTR; \ 153627104ee2SJacob Faibussowitsch petscstack.line[petscstack.currentsize] = 0; \ 153727104ee2SJacob Faibussowitsch petscstack.petscroutine[petscstack.currentsize] = 0; \ 1538441dd030SJed Brown } \ 153927104ee2SJacob Faibussowitsch petscstack.hotdepth = PetscMax(petscstack.hotdepth - 1, 0); \ 1540e04113cfSBarry Smith PetscStackSAWsGrantAccess(); \ 1541441dd030SJed Brown } while (0) 1542441dd030SJed Brown 154330de9b25SBarry Smith /*MC 15441957e957SBarry Smith PetscFunctionBegin - First executable line of each PETSc function, used for error handling. Final 154587497f52SBarry Smith line of PETSc functions should be `PetscFunctionReturn`(0); 154630de9b25SBarry Smith 154730de9b25SBarry Smith Synopsis: 1548aaa7dc30SBarry Smith #include <petscsys.h> 154930de9b25SBarry Smith void PetscFunctionBegin; 155030de9b25SBarry Smith 155120f4b53cSBarry Smith Not Collective; No Fortran Support 1552eca87e8dSBarry Smith 155330de9b25SBarry Smith Usage: 155430de9b25SBarry Smith .vb 155530de9b25SBarry Smith int something; 155630de9b25SBarry Smith 155730de9b25SBarry Smith PetscFunctionBegin; 155830de9b25SBarry Smith .ve 155930de9b25SBarry Smith 156030de9b25SBarry Smith Level: developer 156130de9b25SBarry Smith 156220f4b53cSBarry Smith Note: 156320f4b53cSBarry Smith Use `PetscFunctionBeginUser` for application codes. 156420f4b53cSBarry Smith 1565586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()` 156630de9b25SBarry Smith 156730de9b25SBarry Smith M*/ 15689371c9d4SSatish Balay #define PetscFunctionBegin \ 15699371c9d4SSatish Balay do { \ 1570362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_FALSE); \ 1571a2f94806SJed Brown PetscRegister__FUNCT__(); \ 1572a2f94806SJed Brown } while (0) 1573a2f94806SJed Brown 1574a2f94806SJed Brown /*MC 157587497f52SBarry Smith PetscFunctionBeginHot - Substitute for `PetscFunctionBegin` to be used in functions that are called in 1576a2f94806SJed Brown performance-critical circumstances. Use of this function allows for lighter profiling by default. 1577a2f94806SJed Brown 1578a2f94806SJed Brown Synopsis: 1579aaa7dc30SBarry Smith #include <petscsys.h> 1580a2f94806SJed Brown void PetscFunctionBeginHot; 1581a2f94806SJed Brown 158220f4b53cSBarry Smith Not Collective; No Fortran Support 1583a2f94806SJed Brown 1584a2f94806SJed Brown Usage: 1585a2f94806SJed Brown .vb 1586a2f94806SJed Brown int something; 1587a2f94806SJed Brown 1588a2f94806SJed Brown PetscFunctionBeginHot; 1589a2f94806SJed Brown .ve 1590a2f94806SJed Brown 1591a2f94806SJed Brown Level: developer 1592a2f94806SJed Brown 1593586f9135SBarry Smith .seealso: `PetscFunctionBegin`, `PetscFunctionReturn()`, `PetscStackPushNoCheck()` 1594a2f94806SJed Brown 1595a2f94806SJed Brown M*/ 15969371c9d4SSatish Balay #define PetscFunctionBeginHot \ 15979371c9d4SSatish Balay do { \ 1598362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_TRUE); \ 15992d53ad75SBarry Smith PetscRegister__FUNCT__(); \ 160053c77d0aSJed Brown } while (0) 160153c77d0aSJed Brown 1602a8d2bbe5SBarry Smith /*MC 1603530d85c1SBarry Smith PetscFunctionBeginUser - First executable line of user provided routines 1604a8d2bbe5SBarry Smith 1605a8d2bbe5SBarry Smith Synopsis: 1606aaa7dc30SBarry Smith #include <petscsys.h> 1607a8d2bbe5SBarry Smith void PetscFunctionBeginUser; 1608a8d2bbe5SBarry Smith 160920f4b53cSBarry Smith Not Collective; No Fortran Support 1610a8d2bbe5SBarry Smith 1611a8d2bbe5SBarry Smith Usage: 1612a8d2bbe5SBarry Smith .vb 1613a8d2bbe5SBarry Smith int something; 1614a8d2bbe5SBarry Smith 1615ac285190SBarry Smith PetscFunctionBeginUser; 1616a8d2bbe5SBarry Smith .ve 1617a8d2bbe5SBarry Smith 161820f4b53cSBarry Smith Level: intermediate 161920f4b53cSBarry Smith 1620a8d2bbe5SBarry Smith Notes: 1621530d85c1SBarry Smith Functions that incorporate this must call `PetscFunctionReturn()` instead of return except for main(). 1622530d85c1SBarry Smith 1623530d85c1SBarry Smith May be used before `PetscInitialize()` 16241957e957SBarry Smith 1625530d85c1SBarry Smith This is identical to `PetscFunctionBegin` except it labels the routine as a user 1626ac285190SBarry Smith routine instead of as a PETSc library routine. 1627ac285190SBarry Smith 1628586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, `PetscFunctionBeginHot`, `PetscStackPushNoCheck()` 1629a8d2bbe5SBarry Smith M*/ 16309371c9d4SSatish Balay #define PetscFunctionBeginUser \ 16319371c9d4SSatish Balay do { \ 1632362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 2, PETSC_FALSE); \ 1633a8d2bbe5SBarry Smith PetscRegister__FUNCT__(); \ 1634a8d2bbe5SBarry Smith } while (0) 1635a8d2bbe5SBarry Smith 1636586f9135SBarry Smith /*MC 1637586f9135SBarry Smith PetscStackPush - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is 1638586f9135SBarry Smith currently in the source code and verifies the memory is not corrupted. 1639586f9135SBarry Smith 1640586f9135SBarry Smith Synopsis: 1641586f9135SBarry Smith #include <petscsys.h> 1642586f9135SBarry Smith void PetscStackPush(char *funct) 1643586f9135SBarry Smith 16447cdbe19fSJose E. Roman Not Collective 16457cdbe19fSJose E. Roman 1646586f9135SBarry Smith Input Parameter: 1647586f9135SBarry Smith . funct - the function name 1648586f9135SBarry Smith 1649586f9135SBarry Smith Level: developer 1650586f9135SBarry Smith 1651586f9135SBarry Smith Notes: 1652586f9135SBarry 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 1653586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1654586f9135SBarry Smith help debug the problem. 1655586f9135SBarry Smith 165695bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 1657586f9135SBarry Smith 1658586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1659586f9135SBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1660586f9135SBarry Smith M*/ 16619371c9d4SSatish Balay #define PetscStackPush(n) \ 16629371c9d4SSatish Balay do { \ 1663362febeeSStefano Zampini PetscStackPushNoCheck(n, 0, PETSC_FALSE); \ 166415681b3cSBarry Smith CHKMEMQ; \ 166515681b3cSBarry Smith } while (0) 16663a40ed3dSBarry Smith 1667586f9135SBarry Smith /*MC 1668586f9135SBarry Smith PetscStackPop - Pops a function name from the PETSc default stack that tracks where the running program is 1669586f9135SBarry Smith currently in the source code and verifies the memory is not corrupted. 1670586f9135SBarry Smith 1671586f9135SBarry Smith Synopsis: 1672586f9135SBarry Smith #include <petscsys.h> 1673586f9135SBarry Smith void PetscStackPop 1674586f9135SBarry Smith 16757cdbe19fSJose E. Roman Not Collective 16767cdbe19fSJose E. Roman 1677586f9135SBarry Smith Level: developer 1678586f9135SBarry Smith 1679586f9135SBarry Smith Notes: 1680586f9135SBarry 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 1681586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1682586f9135SBarry Smith help debug the problem. 1683586f9135SBarry Smith 168495bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 1685586f9135SBarry Smith 1686586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPopNoCheck()`, `PetscStackPush()` 1687586f9135SBarry Smith M*/ 16889371c9d4SSatish Balay #define PetscStackPop \ 16899371c9d4SSatish Balay do { \ 1690441dd030SJed Brown CHKMEMQ; \ 1691362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 169215681b3cSBarry Smith } while (0) 1693d64ed03dSBarry Smith 169430de9b25SBarry Smith /*MC 16950a57284eSJacob Faibussowitsch PetscFunctionReturn - Last executable line of each PETSc function used for error 16960a57284eSJacob Faibussowitsch handling. Replaces `return()`. 169730de9b25SBarry Smith 169830de9b25SBarry Smith Synopsis: 16990a57284eSJacob Faibussowitsch #include <petscerror.h> 17000a57284eSJacob Faibussowitsch void PetscFunctionReturn(...) 170130de9b25SBarry Smith 1702667f096bSBarry Smith Not Collective; No Fortran Support 1703eca87e8dSBarry Smith 17045687f061SJacob Faibussowitsch Level: beginner 17055687f061SJacob Faibussowitsch 17060a57284eSJacob Faibussowitsch Notes: 17070a57284eSJacob Faibussowitsch This routine is a macro, so while it does not "return" anything itself, it does return from 17080a57284eSJacob Faibussowitsch the function in the literal sense. 17090a57284eSJacob Faibussowitsch 17100a57284eSJacob Faibussowitsch Usually the return value is the integer literal `0` (for example in any function returning 17110a57284eSJacob Faibussowitsch `PetscErrorCode`), however it is possible to return any arbitrary type. The arguments of 17120a57284eSJacob Faibussowitsch this macro are placed before the `return` statement as-is. 17130a57284eSJacob Faibussowitsch 17140a57284eSJacob Faibussowitsch Any routine which returns via `PetscFunctionReturn()` must begin with a corresponding 17150a57284eSJacob Faibussowitsch `PetscFunctionBegin`. 17160a57284eSJacob Faibussowitsch 17170a57284eSJacob Faibussowitsch For routines which return `void` use `PetscFunctionReturnVoid()` instead. 17180a57284eSJacob Faibussowitsch 17190a57284eSJacob Faibussowitsch Example Usage: 172030de9b25SBarry Smith .vb 17210a57284eSJacob Faibussowitsch PetscErrorCode foo(int *x) 17220a57284eSJacob Faibussowitsch { 17230a57284eSJacob Faibussowitsch PetscFunctionBegin; // don't forget the begin! 17240a57284eSJacob Faibussowitsch *x = 10; 17253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 172630de9b25SBarry Smith } 172730de9b25SBarry Smith .ve 172830de9b25SBarry Smith 17290a57284eSJacob Faibussowitsch May return any arbitrary type\: 17300a57284eSJacob Faibussowitsch .vb 17310a57284eSJacob Faibussowitsch struct Foo 17320a57284eSJacob Faibussowitsch { 17330a57284eSJacob Faibussowitsch int x; 17340a57284eSJacob Faibussowitsch }; 17350a57284eSJacob Faibussowitsch 17360a57284eSJacob Faibussowitsch struct Foo make_foo(int value) 17370a57284eSJacob Faibussowitsch { 17380a57284eSJacob Faibussowitsch struct Foo f; 17390a57284eSJacob Faibussowitsch 17400a57284eSJacob Faibussowitsch PetscFunctionBegin; 17410a57284eSJacob Faibussowitsch f.x = value; 17420a57284eSJacob Faibussowitsch PetscFunctionReturn(f); 17430a57284eSJacob Faibussowitsch } 17440a57284eSJacob Faibussowitsch .ve 17450a57284eSJacob Faibussowitsch 17460a57284eSJacob Faibussowitsch .seealso: `PetscFunctionBegin`, `PetscFunctionBeginUser`, `PetscFunctionReturnVoid()`, 17470a57284eSJacob Faibussowitsch `PetscStackPopNoCheck()` 174830de9b25SBarry Smith M*/ 17495687f061SJacob Faibussowitsch #define PetscFunctionReturn(...) \ 17509371c9d4SSatish Balay do { \ 1751362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 17525687f061SJacob Faibussowitsch return __VA_ARGS__; \ 175327104ee2SJacob Faibussowitsch } while (0) 1754d64ed03dSBarry Smith 17550a57284eSJacob Faibussowitsch /*MC 17560a57284eSJacob Faibussowitsch PetscFunctionReturnVoid - Like `PetscFunctionReturn()` but returns `void` 17570a57284eSJacob Faibussowitsch 17580a57284eSJacob Faibussowitsch Synopsis: 17590a57284eSJacob Faibussowitsch #include <petscerror.h> 17600a57284eSJacob Faibussowitsch void PetscFunctionReturnVoid() 17610a57284eSJacob Faibussowitsch 17620a57284eSJacob Faibussowitsch Not Collective 17630a57284eSJacob Faibussowitsch 17645687f061SJacob Faibussowitsch Level: beginner 17655687f061SJacob Faibussowitsch 17665687f061SJacob Faibussowitsch Note: 17670a57284eSJacob Faibussowitsch Behaves identically to `PetscFunctionReturn()` except that it returns `void`. That is, this 17685687f061SJacob Faibussowitsch macro culminates with `return`. 17690a57284eSJacob Faibussowitsch 17700a57284eSJacob Faibussowitsch Example Usage: 17710a57284eSJacob Faibussowitsch .vb 17720a57284eSJacob Faibussowitsch void foo() 17730a57284eSJacob Faibussowitsch { 17740a57284eSJacob Faibussowitsch PetscFunctionBegin; // must start with PetscFunctionBegin! 17750a57284eSJacob Faibussowitsch bar(); 17760a57284eSJacob Faibussowitsch baz(); 17770a57284eSJacob Faibussowitsch PetscFunctionReturnVoid(); 17780a57284eSJacob Faibussowitsch } 17790a57284eSJacob Faibussowitsch .ve 17800a57284eSJacob Faibussowitsch 17810a57284eSJacob Faibussowitsch .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, PetscFunctionBeginUser` 17820a57284eSJacob Faibussowitsch M*/ 17839371c9d4SSatish Balay #define PetscFunctionReturnVoid() \ 17849371c9d4SSatish Balay do { \ 1785362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 178627104ee2SJacob Faibussowitsch return; \ 178727104ee2SJacob Faibussowitsch } while (0) 178827104ee2SJacob Faibussowitsch #else /* PETSC_USE_DEBUG */ 178927104ee2SJacob Faibussowitsch #define PetscStackPushNoCheck(funct, petsc_routine, hot) 179037154d10SBarry Smith #define PetscStackUpdateLine 1791792fecdfSBarry Smith #define PetscStackPushExternal(funct) 17923ba16761SJacob Faibussowitsch #define PetscStackPopNoCheck(...) 179327104ee2SJacob Faibussowitsch #define PetscStackClearTop 17943a40ed3dSBarry Smith #define PetscFunctionBegin 17950bdf7c52SPeter Brune #define PetscFunctionBeginUser 1796a2f94806SJed Brown #define PetscFunctionBeginHot 17970a57284eSJacob Faibussowitsch #define PetscFunctionReturn(...) return __VA_ARGS__ 17985665465eSBarry Smith #define PetscFunctionReturnVoid() return 1799812af9f3SBarry Smith #define PetscStackPop CHKMEMQ 1800812af9f3SBarry Smith #define PetscStackPush(f) CHKMEMQ 180127104ee2SJacob Faibussowitsch #endif /* PETSC_USE_DEBUG */ 18023a40ed3dSBarry Smith 18036d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 18043ba16761SJacob Faibussowitsch #define PetscStackCallExternalVoid(...) 18053ba16761SJacob Faibussowitsch template <typename F, typename... Args> 18063ba16761SJacob Faibussowitsch void PetscCallExternal(F, Args...); 1807*1fc6d13cSAlexander template <typename F, typename... Args> 1808*1fc6d13cSAlexander void PetscCallExternalAbort(F, Args...); 18096d210af2SJacob Faibussowitsch #else 1810586f9135SBarry Smith /*MC 1811e77caa6dSBarry Smith PetscStackCallExternalVoid - Calls an external library routine or user function after pushing the name of the routine on the stack. 1812eb6b5d47SBarry Smith 1813eb6b5d47SBarry Smith Input Parameters: 1814eb6b5d47SBarry Smith + name - string that gives the name of the function being called 1815586f9135SBarry Smith - routine - actual call to the routine, for example, functionname(a,b) 1816fd3f9acdSBarry Smith 1817586f9135SBarry Smith Level: developer 1818eb6b5d47SBarry Smith 181995bd0b28SBarry Smith Notes: 1820792fecdfSBarry Smith Often one should use `PetscCallExternal()` instead. This routine is intended for external library routines that DO NOT return error codes 1821eb6b5d47SBarry Smith 1822586f9135SBarry Smith In debug mode this also checks the memory for corruption at the end of the function call. 1823586f9135SBarry Smith 182495bd0b28SBarry Smith Certain external packages, such as BLAS/LAPACK may have their own macros, `PetscCallBLAS()` for managing the call, error checking, etc. 1825586f9135SBarry Smith 1826586f9135SBarry Smith Developer Note: 1827586f9135SBarry 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. 1828586f9135SBarry Smith 1829792fecdfSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscCallExternal()`, `PetscCallBLAS()` 1830586f9135SBarry Smith @*/ 18313ba16761SJacob Faibussowitsch #define PetscStackCallExternalVoid(name, ...) \ 18329371c9d4SSatish Balay do { \ 183328770333SStefano Zampini PetscStackPushExternal(name); \ 18343ba16761SJacob Faibussowitsch __VA_ARGS__; \ 18359371c9d4SSatish Balay PetscStackPop; \ 18369371c9d4SSatish Balay } while (0) 1837eb6b5d47SBarry Smith 1838586f9135SBarry Smith /*MC 1839792fecdfSBarry Smith PetscCallExternal - Calls an external library routine that returns an error code after pushing the name of the routine on the stack. 1840fd3f9acdSBarry Smith 1841fd3f9acdSBarry Smith Input Parameters: 1842fd3f9acdSBarry Smith + func - name of the routine 1843586f9135SBarry Smith - args - arguments to the routine 1844586f9135SBarry Smith 1845586f9135SBarry Smith Level: developer 1846fd3f9acdSBarry Smith 184795452b02SPatrick Sanan Notes: 1848e77caa6dSBarry Smith This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not. 1849dbf62e16SBarry Smith 1850586f9135SBarry Smith In debug mode this also checks the memory for corruption at the end of the function call. 1851fd3f9acdSBarry Smith 185287497f52SBarry Smith Assumes the error return code of the function is an integer and that a value of 0 indicates success 185387497f52SBarry Smith 1854586f9135SBarry Smith Developer Note: 1855d5b43468SJose E. Roman This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc. 1856586f9135SBarry Smith 1857*1fc6d13cSAlexander .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()`, `PetscCallExternalAbort()` 1858586f9135SBarry Smith M*/ 18599371c9d4SSatish Balay #define PetscCallExternal(func, ...) \ 18609371c9d4SSatish Balay do { \ 1861a74df02fSJacob Faibussowitsch PetscStackPush(PetscStringize(func)); \ 18623ba16761SJacob Faibussowitsch int ierr_petsc_call_external_ = func(__VA_ARGS__); \ 18631d4906efSStefano Zampini PetscStackPop; \ 18643ba16761SJacob 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_); \ 1865fd3f9acdSBarry Smith } while (0) 1866*1fc6d13cSAlexander 1867*1fc6d13cSAlexander /*MC 1868*1fc6d13cSAlexander 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 1869*1fc6d13cSAlexander 1870*1fc6d13cSAlexander Input Parameters: 1871*1fc6d13cSAlexander + func - name of the routine 1872*1fc6d13cSAlexander - args - arguments to the routine 1873*1fc6d13cSAlexander 1874*1fc6d13cSAlexander Level: developer 1875*1fc6d13cSAlexander 1876*1fc6d13cSAlexander Notes: 1877*1fc6d13cSAlexander This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not. 1878*1fc6d13cSAlexander 1879*1fc6d13cSAlexander In debug mode this also checks the memory for corruption at the end of the function call. 1880*1fc6d13cSAlexander 1881*1fc6d13cSAlexander Assumes the error return code of the function is an integer and that a value of 0 indicates success 1882*1fc6d13cSAlexander 1883*1fc6d13cSAlexander Developer Note: 1884*1fc6d13cSAlexander This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc. 1885*1fc6d13cSAlexander 1886*1fc6d13cSAlexander .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()`, `PetscCallExternal()` 1887*1fc6d13cSAlexander M*/ 1888*1fc6d13cSAlexander #define PetscCallExternalAbort(func, ...) \ 1889*1fc6d13cSAlexander do { \ 1890*1fc6d13cSAlexander PetscStackUpdateLine; \ 1891*1fc6d13cSAlexander int ierr_petsc_call_external_ = func(__VA_ARGS__); \ 1892*1fc6d13cSAlexander if (PetscUnlikely(ierr_petsc_call_external_ != 0)) { \ 1893*1fc6d13cSAlexander (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_); \ 1894*1fc6d13cSAlexander PETSCABORTWITHIERR_Private(PETSC_COMM_SELF, PETSC_ERR_LIB); \ 1895*1fc6d13cSAlexander } \ 1896*1fc6d13cSAlexander } while (0) 18976d210af2SJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */ 1898