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 60ffc4695bSBarry Smith /* 61ffc4695bSBarry Smith Returned from PETSc functions that are called from MPI, such as related to attributes 62ffc4695bSBarry Smith Do not confuse PETSC_MPI_ERROR_CODE and PETSC_ERR_MPI, the first is registered with MPI and returned to MPI as 63888a1fb3SBarry 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. 64ffc4695bSBarry Smith */ 65ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CLASS; 66ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CODE; 67ffc4695bSBarry Smith 68986eef2eSBarry Smith /*MC 69986eef2eSBarry Smith SETERRMPI - Macro to be called when an error has been detected within an MPI callback function 70986eef2eSBarry Smith 71989c49a2SBarry Smith No Fortran Support 72989c49a2SBarry Smith 73986eef2eSBarry Smith Synopsis: 74986eef2eSBarry Smith #include <petscsys.h> 7598921bdaSJacob Faibussowitsch PetscErrorCode SETERRMPI(MPI_Comm comm, PetscErrorCode ierr, char *message, ...) 76986eef2eSBarry Smith 77d083f849SBarry Smith Collective 78986eef2eSBarry Smith 79986eef2eSBarry Smith Input Parameters: 8095bd0b28SBarry Smith + comm - An MPI communicator, use `PETSC_COMM_SELF` unless you know all ranks of another communicator will detect the error 81986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 82986eef2eSBarry Smith - message - error message 83986eef2eSBarry Smith 84986eef2eSBarry Smith Level: developer 85986eef2eSBarry Smith 8695bd0b28SBarry Smith Note: 8787497f52SBarry 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` 8887497f52SBarry Smith which is registered with `MPI_Add_error_code()` when PETSc is initialized. 89986eef2eSBarry Smith 90af27ebaaSBarry Smith .seealso: `SETERRQ()`, `PetscCall()`, `PetscCallMPI()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `PetscErrorCode` 91986eef2eSBarry Smith M*/ 923ba16761SJacob Faibussowitsch #define SETERRMPI(comm, ierr, ...) return ((void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__), PETSC_MPI_ERROR_CODE) 93ee8199e6SMatthew G. Knepley 94ee8199e6SMatthew G. Knepley /*MC 95f388eb8bSPatrick Sanan SETERRA - Fortran-only macro that can be called when an error has been detected from the main program 96f388eb8bSPatrick Sanan 97f388eb8bSPatrick Sanan Synopsis: 98f388eb8bSPatrick Sanan #include <petscsys.h> 99f388eb8bSPatrick Sanan PetscErrorCode SETERRA(MPI_Comm comm, PetscErrorCode ierr, char *message) 100f388eb8bSPatrick Sanan 101f388eb8bSPatrick Sanan Collective 102f388eb8bSPatrick Sanan 103f388eb8bSPatrick Sanan Input Parameters: 10495bd0b28SBarry Smith + comm - An MPI communicator, so that the error can be collective 105f388eb8bSPatrick Sanan . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 106f388eb8bSPatrick Sanan - message - error message in the printf format 107f388eb8bSPatrick Sanan 108f388eb8bSPatrick Sanan Level: beginner 109f388eb8bSPatrick Sanan 110f388eb8bSPatrick Sanan Notes: 11187497f52SBarry Smith This should only be used with Fortran. With C/C++, use `SETERRQ()`. 112f388eb8bSPatrick Sanan 11387497f52SBarry Smith `SETERRQ()` may be called from Fortran subroutines but `SETERRA()` must be called from the 114f388eb8bSPatrick Sanan Fortran main program. 115f388eb8bSPatrick Sanan 116af27ebaaSBarry Smith .seealso: `SETERRQ()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()`, `PetscErrorCode` 117f388eb8bSPatrick Sanan M*/ 118f388eb8bSPatrick Sanan 119f388eb8bSPatrick Sanan /*MC 120fa190f98SMatthew G. Knepley SETERRABORT - Macro that can be called when an error has been detected, 121fa190f98SMatthew G. Knepley 122fa190f98SMatthew G. Knepley Synopsis: 123fa190f98SMatthew G. Knepley #include <petscsys.h> 12498921bdaSJacob Faibussowitsch PetscErrorCode SETERRABORT(MPI_Comm comm, PetscErrorCode ierr, char *message, ...) 125fa190f98SMatthew G. Knepley 126d083f849SBarry Smith Collective 127fa190f98SMatthew G. Knepley 128fa190f98SMatthew G. Knepley Input Parameters: 12995bd0b28SBarry Smith + comm - An MPI communicator, so that the error can be collective 1303af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 131fa190f98SMatthew G. Knepley - message - error message in the printf format 132fa190f98SMatthew G. Knepley 133fa190f98SMatthew G. Knepley Level: beginner 134fa190f98SMatthew G. Knepley 135fa190f98SMatthew G. Knepley Notes: 13687497f52SBarry Smith This function just calls `MPI_Abort()`. 13787497f52SBarry Smith 13887497f52SBarry Smith This should only be called in routines that cannot return an error code, such as in C++ constructors. 139fa190f98SMatthew G. Knepley 140989c49a2SBarry Smith Fortran Note: 141989c49a2SBarry Smith Use `SETERRA()` in Fortran main program and `SETERRQ()` in Fortran subroutines 142989c49a2SBarry Smith 143989c49a2SBarry Smith Developer Note: 144989c49a2SBarry Smith In Fortran `SETERRA()` could be called `SETERRABORT()` since they serve the same purpose 145989c49a2SBarry Smith 146af27ebaaSBarry Smith .seealso: `SETERRQ()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `PetscCall()`, `CHKMEMQ`, `PetscErrorCode` 147fa190f98SMatthew G. Knepley M*/ 1489371c9d4SSatish Balay #define SETERRABORT(comm, ierr, ...) \ 1499371c9d4SSatish Balay do { \ 1503ba16761SJacob Faibussowitsch (void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \ 15198921bdaSJacob Faibussowitsch MPI_Abort(comm, ierr); \ 15298921bdaSJacob Faibussowitsch } while (0) 1539a00fa46SSatish Balay 15430de9b25SBarry Smith /*MC 155139c8099SSatish Balay PetscCheck - Checks that a particular condition is true; if not true, then returns the provided error code 1562c71b3e2SJacob Faibussowitsch 1572c71b3e2SJacob Faibussowitsch Synopsis: 1582c71b3e2SJacob Faibussowitsch #include <petscerror.h> 1592c71b3e2SJacob Faibussowitsch void PetscCheck(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 1602c71b3e2SJacob Faibussowitsch 1617cdbe19fSJose E. Roman Collective; No Fortran Support 1622c71b3e2SJacob Faibussowitsch 1632c71b3e2SJacob Faibussowitsch Input Parameters: 1642c71b3e2SJacob Faibussowitsch + cond - The boolean condition 1652c71b3e2SJacob Faibussowitsch . comm - The communicator on which the check can be collective on 1662c71b3e2SJacob Faibussowitsch . ierr - A nonzero error code, see include/petscerror.h for the complete list 1672c71b3e2SJacob Faibussowitsch - message - Error message in printf format 1682c71b3e2SJacob Faibussowitsch 169667f096bSBarry Smith Level: beginner 170667f096bSBarry Smith 1712c71b3e2SJacob Faibussowitsch Notes: 1722c71b3e2SJacob Faibussowitsch Enabled in both optimized and debug builds. 1732c71b3e2SJacob Faibussowitsch 174a2454668SBarry Smith As a general rule, `PetscCheck()` is used to check "usage error" (for example, passing an incorrect value as a function argument), 175a2454668SBarry Smith `PetscAssert()` is used to "check for bugs in PETSc" (for example, is a value in a PETSc data structure nonsensical). 176a2454668SBarry Smith However, for functions that are called in a "hot spot", for example, thousands of times in a loop, `PetscAssert()` should be used instead 177a2454668SBarry Smith of `PetscCheck()` since the former is compiled out in PETSc's optimization code. 178a2454668SBarry Smith 17987497f52SBarry Smith Calls `SETERRQ()` if the assertion fails, so can only be called from functions returning a 18087497f52SBarry Smith `PetscErrorCode` (or equivalent type after conversion). 1812c71b3e2SJacob Faibussowitsch 182af27ebaaSBarry Smith .seealso: `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheckAbort()`, `PetscErrorCode` 1839566063dSJacob Faibussowitsch M*/ 1849371c9d4SSatish Balay #define PetscCheck(cond, comm, ierr, ...) \ 1859371c9d4SSatish Balay do { \ 1869371c9d4SSatish Balay if (PetscUnlikely(!(cond))) SETERRQ(comm, ierr, __VA_ARGS__); \ 1879371c9d4SSatish Balay } while (0) 1882c71b3e2SJacob Faibussowitsch 1892c71b3e2SJacob Faibussowitsch /*MC 1904be741a6SBarry Smith PetscCheckAbort - Check that a particular condition is true, otherwise prints error and aborts 1914be741a6SBarry Smith 1924be741a6SBarry Smith Synopsis: 1934be741a6SBarry Smith #include <petscerror.h> 1944be741a6SBarry Smith void PetscCheckAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 1954be741a6SBarry Smith 1967cdbe19fSJose E. Roman Collective; No Fortran Support 1974be741a6SBarry Smith 1984be741a6SBarry Smith Input Parameters: 1994be741a6SBarry Smith + cond - The boolean condition 2004be741a6SBarry Smith . comm - The communicator on which the check can be collective on 2014be741a6SBarry Smith . ierr - A nonzero error code, see include/petscerror.h for the complete list 2024be741a6SBarry Smith - message - Error message in printf format 2034be741a6SBarry Smith 204667f096bSBarry Smith Level: developer 205667f096bSBarry Smith 2064be741a6SBarry Smith Notes: 2074be741a6SBarry Smith Enabled in both optimized and debug builds. 2084be741a6SBarry Smith 20987497f52SBarry Smith Calls `SETERRABORT()` if the assertion fails, can be called from a function that does not return an 21087497f52SBarry Smith error code, such as a C++ constructor. usually `PetscCheck()` should be used. 2114be741a6SBarry Smith 212af27ebaaSBarry Smith .seealso: `PetscAssertAbort()`, `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheck()`, `SETERRABORT()`, `PetscErrorCode` 2134be741a6SBarry Smith M*/ 2149371c9d4SSatish Balay #define PetscCheckAbort(cond, comm, ierr, ...) \ 2150e6b6b59SJacob Faibussowitsch do { \ 2160e6b6b59SJacob Faibussowitsch if (PetscUnlikely(!(cond))) SETERRABORT(comm, ierr, __VA_ARGS__); \ 217c7481402SJacob Faibussowitsch } while (0) 2184be741a6SBarry Smith 2194be741a6SBarry Smith /*MC 2209ace16cdSJacob Faibussowitsch PetscAssert - Assert that a particular condition is true 2219ace16cdSJacob Faibussowitsch 2229ace16cdSJacob Faibussowitsch Synopsis: 2239ace16cdSJacob Faibussowitsch #include <petscerror.h> 2249ace16cdSJacob Faibussowitsch void PetscAssert(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 2259ace16cdSJacob Faibussowitsch 2267cdbe19fSJose E. Roman Collective; No Fortran Support 2279ace16cdSJacob Faibussowitsch 2289ace16cdSJacob Faibussowitsch Input Parameters: 2299ace16cdSJacob Faibussowitsch + cond - The boolean condition 2309ace16cdSJacob Faibussowitsch . comm - The communicator on which the check can be collective on 2319ace16cdSJacob Faibussowitsch . ierr - A nonzero error code, see include/petscerror.h for the complete list 232af27ebaaSBarry Smith - message - Error message in `printf()` format 2339ace16cdSJacob Faibussowitsch 234667f096bSBarry Smith Level: beginner 235667f096bSBarry Smith 2369ace16cdSJacob Faibussowitsch Notes: 237c7481402SJacob Faibussowitsch Equivalent to `PetscCheck()` if debugging is enabled, and `PetscAssume(cond)` otherwise. 2382c71b3e2SJacob Faibussowitsch 23987497f52SBarry Smith See `PetscCheck()` for usage and behaviour. 24087497f52SBarry Smith 24187497f52SBarry Smith This is needed instead of simply using `assert()` because this correctly handles the collective nature of errors under MPI 2429ace16cdSJacob Faibussowitsch 243af27ebaaSBarry Smith .seealso: `PetscCheck()`, `SETERRQ()`, `PetscError()`, `PetscAssertAbort()`, `PetscErrorCode` 2449566063dSJacob Faibussowitsch M*/ 245c7481402SJacob Faibussowitsch #if PetscDefined(USE_DEBUG) 246c7481402SJacob Faibussowitsch #define PetscAssert(cond, comm, ierr, ...) PetscCheck(cond, comm, ierr, __VA_ARGS__) 247c7481402SJacob Faibussowitsch #else 248c7481402SJacob Faibussowitsch #define PetscAssert(cond, ...) PetscAssume(cond) 249c7481402SJacob Faibussowitsch #endif 2509ace16cdSJacob Faibussowitsch 2519ace16cdSJacob Faibussowitsch /*MC 2520e6b6b59SJacob Faibussowitsch PetscAssertAbort - Assert that a particular condition is true, otherwise prints error and aborts 2530e6b6b59SJacob Faibussowitsch 2540e6b6b59SJacob Faibussowitsch Synopsis: 2550e6b6b59SJacob Faibussowitsch #include <petscerror.h> 2560e6b6b59SJacob Faibussowitsch void PetscAssertAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 2570e6b6b59SJacob Faibussowitsch 2587cdbe19fSJose E. Roman Collective; No Fortran Support 2590e6b6b59SJacob Faibussowitsch 2600e6b6b59SJacob Faibussowitsch Input Parameters: 2610e6b6b59SJacob Faibussowitsch + cond - The boolean condition 2620e6b6b59SJacob Faibussowitsch . comm - The communicator on which the check can be collective on 2630e6b6b59SJacob Faibussowitsch . ierr - A nonzero error code, see include/petscerror.h for the complete list 2640e6b6b59SJacob Faibussowitsch - message - Error message in printf format 2650e6b6b59SJacob Faibussowitsch 266667f096bSBarry Smith Level: beginner 267667f096bSBarry Smith 26895bd0b28SBarry Smith Note: 2690e6b6b59SJacob Faibussowitsch Enabled only in debug builds. See `PetscCheckAbort()` for usage. 2700e6b6b59SJacob Faibussowitsch 2710e6b6b59SJacob Faibussowitsch .seealso: `PetscCheckAbort()`, `PetscAssert()`, `PetscCheck()`, `SETERRABORT()`, `PetscError()` 2720e6b6b59SJacob Faibussowitsch M*/ 2733ba16761SJacob Faibussowitsch #if PetscDefined(USE_DEBUG) 2743ba16761SJacob Faibussowitsch #define PetscAssertAbort(cond, comm, ierr, ...) PetscCheckAbort(cond, comm, ierr, __VA_ARGS__) 2753ba16761SJacob Faibussowitsch #else 2763ba16761SJacob Faibussowitsch #define PetscAssertAbort(cond, comm, ierr, ...) PetscAssume(cond) 2773ba16761SJacob Faibussowitsch #endif 2780e6b6b59SJacob Faibussowitsch 2790e6b6b59SJacob Faibussowitsch /*MC 2803ba16761SJacob Faibussowitsch PetscCall - Calls a PETSc function and then checks the resulting error code, if it is 2813ba16761SJacob Faibussowitsch non-zero it calls the error handler and returns from the current function with the error 2823ba16761SJacob Faibussowitsch code. 2839566063dSJacob Faibussowitsch 2849566063dSJacob Faibussowitsch Synopsis: 2859566063dSJacob Faibussowitsch #include <petscerror.h> 28649c86fc7SBarry Smith void PetscCall(PetscFunction(args)) 2879566063dSJacob Faibussowitsch 2889566063dSJacob Faibussowitsch Not Collective 2899566063dSJacob Faibussowitsch 2909566063dSJacob Faibussowitsch Input Parameter: 29149c86fc7SBarry Smith . PetscFunction - any PETSc function that returns an error code 2929566063dSJacob Faibussowitsch 293667f096bSBarry Smith Level: beginner 294667f096bSBarry Smith 2959566063dSJacob Faibussowitsch Notes: 2969566063dSJacob Faibussowitsch Once the error handler is called the calling function is then returned from with the given 29787497f52SBarry Smith error code. Experienced users can set the error handler with `PetscPushErrorHandler()`. 2989566063dSJacob Faibussowitsch 29987497f52SBarry Smith `PetscCall()` cannot be used in functions returning a datatype not convertible to 3008af9f190SBarry Smith `PetscErrorCode`. For example, `PetscCall()` may not be used in functions returning `void`, use 3013ba16761SJacob Faibussowitsch `PetscCallAbort()` or `PetscCallVoid()` in this case. 3029566063dSJacob Faibussowitsch 3039566063dSJacob Faibussowitsch Example Usage: 3049566063dSJacob Faibussowitsch .vb 3059566063dSJacob Faibussowitsch PetscCall(PetscInitiailize(...)); // OK to call even when PETSc is not yet initialized! 3069566063dSJacob Faibussowitsch 3079566063dSJacob Faibussowitsch struct my_struct 3089566063dSJacob Faibussowitsch { 3099566063dSJacob Faibussowitsch void *data; 3109566063dSJacob Faibussowitsch } my_complex_type; 3119566063dSJacob Faibussowitsch 3129566063dSJacob Faibussowitsch struct my_struct bar(void) 3139566063dSJacob Faibussowitsch { 3146aad120cSJose E. Roman PetscCall(foo(15)); // ERROR PetscErrorCode not convertible to struct my_struct! 3159566063dSJacob Faibussowitsch } 3169566063dSJacob Faibussowitsch 3176aad120cSJose E. Roman PetscCall(bar()) // ERROR input not convertible to PetscErrorCode 3189566063dSJacob Faibussowitsch .ve 3199566063dSJacob Faibussowitsch 32087497f52SBarry Smith It is also possible to call this directly on a `PetscErrorCode` variable 32149c86fc7SBarry Smith .vb 32249c86fc7SBarry Smith PetscCall(ierr); // check if ierr is nonzero 32349c86fc7SBarry Smith .ve 32449c86fc7SBarry Smith 325792fecdfSBarry Smith Should not be used to call callback functions provided by users, `PetscCallBack()` should be used in that situation. 326ef1023bdSBarry Smith 3276a8be23eSBarry Smith `PetscUseTypeMethod()` or `PetscTryTypeMethod()` should be used when calling functions pointers contained in a PETSc object's `ops` array 3286a8be23eSBarry Smith 32949c86fc7SBarry Smith Fortran Notes: 33098d50953SPierre Jolivet The Fortran function in which this is used must declare a `PetscErrorCode` variable necessarily named `ierr`, and `ierr` must be 33187497f52SBarry Smith the final argument to the PETSc function being called. 33249c86fc7SBarry Smith 33398d50953SPierre Jolivet In the main program and in Fortran subroutines that do not have `ierr` as the final return parameter, one 33487497f52SBarry Smith should use `PetscCallA()` 33549c86fc7SBarry Smith 33649c86fc7SBarry Smith Example Fortran Usage: 33749c86fc7SBarry Smith .vb 33849c86fc7SBarry Smith PetscErrorCode ierr 33949c86fc7SBarry Smith Vec v 34049c86fc7SBarry Smith 34149c86fc7SBarry Smith ... 34249c86fc7SBarry Smith PetscCall(VecShift(v, 1.0, ierr)) 34349c86fc7SBarry Smith PetscCallA(VecShift(v, 1.0, ierr)) 34449c86fc7SBarry Smith .ve 34549c86fc7SBarry Smith 346667f096bSBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`, 347667f096bSBarry Smith `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, 3483ba16761SJacob Faibussowitsch `CHKERRMPI()`, `PetscCallBack()`, `PetscCallAbort()`, `PetscCallVoid()` 3499566063dSJacob Faibussowitsch M*/ 350ef1023bdSBarry Smith 351ef1023bdSBarry Smith /*MC 35298d50953SPierre 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 35398d50953SPierre Jolivet `PetscCall()` which should be used in other Fortran subroutines 354989c49a2SBarry Smith 355989c49a2SBarry Smith Synopsis: 356989c49a2SBarry Smith #include <petscsys.h> 357989c49a2SBarry Smith PetscErrorCode PetscCallA(PetscFunction(arguments, ierr)) 358989c49a2SBarry Smith 359989c49a2SBarry Smith Collective 360989c49a2SBarry Smith 361989c49a2SBarry Smith Input Parameter: 362989c49a2SBarry Smith . PetscFunction(arguments,ierr) - the call to the function 363989c49a2SBarry Smith 364989c49a2SBarry Smith Level: beginner 365989c49a2SBarry Smith 366989c49a2SBarry Smith Notes: 367989c49a2SBarry Smith This should only be used with Fortran. With C/C++, use `PetscCall()` always. 368989c49a2SBarry Smith 36998d50953SPierre Jolivet The Fortran function in which this is used must declare a `PetscErrorCode` variable necessarily named `ierr` 370989c49a2SBarry Smith Use `SETERRA()` to set an error in a Fortran main program and `SETERRQ()` in Fortran subroutines 371989c49a2SBarry Smith 372989c49a2SBarry Smith .seealso: `SETERRQ()`, `SETERRA()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()` 373989c49a2SBarry Smith M*/ 374989c49a2SBarry Smith 375989c49a2SBarry Smith /*MC 376792fecdfSBarry 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 377ef1023bdSBarry Smith handler and returns from the current function with the error code. 378ef1023bdSBarry Smith 379ef1023bdSBarry Smith Synopsis: 380ef1023bdSBarry Smith #include <petscerror.h> 381792fecdfSBarry Smith void PetscCallBack(const char *functionname, PetscFunction(args)) 382ef1023bdSBarry Smith 3837cdbe19fSJose E. Roman Not Collective; No Fortran Support 384ef1023bdSBarry Smith 385ef1023bdSBarry Smith Input Parameters: 386ef1023bdSBarry Smith + functionname - the name of the function being called, this can be a string with spaces that describes the meaning of the callback 38787497f52SBarry Smith - PetscFunction - user provided callback function that returns an error code 388ef1023bdSBarry Smith 389ef1023bdSBarry Smith Example Usage: 390ef1023bdSBarry Smith .vb 391792fecdfSBarry Smith PetscCallBack("XXX callback to do something", a->callback(...)); 392ef1023bdSBarry Smith .ve 393ef1023bdSBarry Smith 394ef1023bdSBarry Smith Level: developer 395ef1023bdSBarry Smith 396667f096bSBarry Smith Notes: 397*7e6f8dd6SBarry Smith `PetscUseTypeMethod()` and ` PetscTryTypeMethod()` are the preferred API for this functionality. But when the callback functions are associated with a 398*7e6f8dd6SBarry Smith `DMSNES` or `DMTS` this API must be used. 399*7e6f8dd6SBarry Smith 400667f096bSBarry Smith Once the error handler is called the calling function is then returned from with the given 401667f096bSBarry Smith error code. Experienced users can set the error handler with `PetscPushErrorHandler()`. 402667f096bSBarry Smith 403667f096bSBarry Smith `PetscCallBack()` should only be called in PETSc when a call is being made to a user provided call-back routine. 404667f096bSBarry Smith 405*7e6f8dd6SBarry Smith Developer Note: 406*7e6f8dd6SBarry 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 407*7e6f8dd6SBarry Smith 40887497f52SBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscCall()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()` 409*7e6f8dd6SBarry Smith `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, `CHKERRMPI()`, `PetscCall()`, `PetscUseTypeMethod()`, `PetscTryTypeMethod()` 410ef1023bdSBarry Smith M*/ 411ef1023bdSBarry Smith 4123ba16761SJacob Faibussowitsch /*MC 4138af9f190SBarry Smith PetscCallVoid - Like `PetscCall()` but for use in functions that return `void` 4143ba16761SJacob Faibussowitsch 4153ba16761SJacob Faibussowitsch Synopsis: 4163ba16761SJacob Faibussowitsch #include <petscerror.h> 4178af9f190SBarry Smith void PetscCallVoid(PetscFunction(args)) 4183ba16761SJacob Faibussowitsch 4197cdbe19fSJose E. Roman Not Collective; No Fortran Support 4203ba16761SJacob Faibussowitsch 4213ba16761SJacob Faibussowitsch Input Parameter: 4223ba16761SJacob Faibussowitsch . PetscFunction - any PETSc function that returns an error code 4233ba16761SJacob Faibussowitsch 4243ba16761SJacob Faibussowitsch Example Usage: 4253ba16761SJacob Faibussowitsch .vb 4263ba16761SJacob Faibussowitsch void foo() 4273ba16761SJacob Faibussowitsch { 4283ba16761SJacob Faibussowitsch KSP ksp; 4293ba16761SJacob Faibussowitsch 4303ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 4313ba16761SJacob Faibussowitsch // OK, properly handles PETSc error codes 4323ba16761SJacob Faibussowitsch PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp)); 4338af9f190SBarry Smith PetscFunctionReturnVoid(); 4343ba16761SJacob Faibussowitsch } 4353ba16761SJacob Faibussowitsch 4363ba16761SJacob Faibussowitsch PetscErrorCode bar() 4373ba16761SJacob Faibussowitsch { 4383ba16761SJacob Faibussowitsch KSP ksp; 4393ba16761SJacob Faibussowitsch 4403ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 4413ba16761SJacob Faibussowitsch // ERROR, Non-void function 'bar' should return a value 4423ba16761SJacob Faibussowitsch PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp)); 4433ba16761SJacob Faibussowitsch // OK, returning PetscErrorCode 4443ba16761SJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp)); 4453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4463ba16761SJacob Faibussowitsch } 447667f096bSBarry Smith .ve 4483ba16761SJacob Faibussowitsch 4493ba16761SJacob Faibussowitsch Level: beginner 4503ba16761SJacob Faibussowitsch 451667f096bSBarry Smith Notes: 452667f096bSBarry Smith Has identical usage to `PetscCall()`, except that it returns `void` on error instead of a 453667f096bSBarry Smith `PetscErrorCode`. See `PetscCall()` for more detailed discussion. 454667f096bSBarry Smith 455667f096bSBarry Smith Note that users should prefer `PetscCallAbort()` to this routine. While this routine does 456667f096bSBarry Smith "handle" errors by returning from the enclosing function, it effectively gobbles the 457667f096bSBarry Smith error. Since the enclosing function itself returns `void`, its callers have no way of knowing 458667f096bSBarry Smith that the routine returned early due to an error. `PetscCallAbort()` at least ensures that the 459667f096bSBarry Smith program crashes gracefully. 460667f096bSBarry Smith 4618af9f190SBarry Smith .seealso: `PetscCall()`, `PetscErrorCode`, `PetscCallAbort()` 4623ba16761SJacob Faibussowitsch M*/ 4639566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 4649566063dSJacob Faibussowitsch void PetscCall(PetscErrorCode); 465792fecdfSBarry Smith void PetscCallBack(const char *, PetscErrorCode); 4669566063dSJacob Faibussowitsch void PetscCallVoid(PetscErrorCode); 4679566063dSJacob Faibussowitsch #else 4689371c9d4SSatish Balay #define PetscCall(...) \ 4699371c9d4SSatish Balay do { \ 4703ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_q_; \ 47137154d10SBarry Smith PetscStackUpdateLine; \ 4723ba16761SJacob Faibussowitsch ierr_petsc_call_q_ = __VA_ARGS__; \ 4733ba16761SJacob 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, " "); \ 4749566063dSJacob Faibussowitsch } while (0) 4759371c9d4SSatish Balay #define PetscCallBack(function, ...) \ 4769371c9d4SSatish Balay do { \ 4773ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_q_; \ 478e33ced7fSLisandro Dalcin PetscStackUpdateLine; \ 479e33ced7fSLisandro Dalcin PetscStackPushExternal(function); \ 4803ba16761SJacob Faibussowitsch ierr_petsc_call_q_ = __VA_ARGS__; \ 481ef1023bdSBarry Smith PetscStackPop; \ 4823ba16761SJacob 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, " "); \ 483ef1023bdSBarry Smith } while (0) 4849371c9d4SSatish Balay #define PetscCallVoid(...) \ 4859371c9d4SSatish Balay do { \ 4863ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_void_; \ 487e33ced7fSLisandro Dalcin PetscStackUpdateLine; \ 4883ba16761SJacob Faibussowitsch ierr_petsc_call_void_ = __VA_ARGS__; \ 4893ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_void_ != PETSC_SUCCESS)) { \ 4903ba16761SJacob Faibussowitsch ierr_petsc_call_void_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_void_, PETSC_ERROR_REPEAT, " "); \ 4913ba16761SJacob Faibussowitsch (void)ierr_petsc_call_void_; \ 4929371c9d4SSatish Balay return; \ 4939371c9d4SSatish Balay } \ 4949566063dSJacob Faibussowitsch } while (0) 4959566063dSJacob Faibussowitsch #endif 4969566063dSJacob Faibussowitsch 4979566063dSJacob Faibussowitsch /*MC 4989566063dSJacob Faibussowitsch CHKERRQ - Checks error code returned from PETSc function 49930de9b25SBarry Smith 50030de9b25SBarry Smith Synopsis: 501aaa7dc30SBarry Smith #include <petscsys.h> 5029566063dSJacob Faibussowitsch void CHKERRQ(PetscErrorCode ierr) 5039566063dSJacob Faibussowitsch 5049566063dSJacob Faibussowitsch Not Collective 5059566063dSJacob Faibussowitsch 5062fe279fdSBarry Smith Input Parameter: 5079566063dSJacob Faibussowitsch . ierr - nonzero error code 5089566063dSJacob Faibussowitsch 5099566063dSJacob Faibussowitsch Level: deprecated 5109566063dSJacob Faibussowitsch 511667f096bSBarry Smith Note: 512667f096bSBarry Smith Deprecated in favor of `PetscCall()`. This routine behaves identically to it. 513667f096bSBarry Smith 514db781477SPatrick Sanan .seealso: `PetscCall()` 5159566063dSJacob Faibussowitsch M*/ 5169566063dSJacob Faibussowitsch #define CHKERRQ(...) PetscCall(__VA_ARGS__) 5179566063dSJacob Faibussowitsch #define CHKERRV(...) PetscCallVoid(__VA_ARGS__) 5189566063dSJacob Faibussowitsch 519db9cea48SBarry Smith PETSC_EXTERN void PetscMPIErrorString(PetscMPIInt, char *); 520db9cea48SBarry Smith 5219566063dSJacob Faibussowitsch /*MC 5229566063dSJacob Faibussowitsch PetscCallMPI - Checks error code returned from MPI calls, if non-zero it calls the error 5239566063dSJacob Faibussowitsch handler and then returns 5249566063dSJacob Faibussowitsch 5259566063dSJacob Faibussowitsch Synopsis: 5269566063dSJacob Faibussowitsch #include <petscerror.h> 52749c86fc7SBarry Smith void PetscCallMPI(MPI_Function(args)) 52830de9b25SBarry Smith 529eca87e8dSBarry Smith Not Collective 53030de9b25SBarry Smith 5312fe279fdSBarry Smith Input Parameter: 53249c86fc7SBarry Smith . MPI_Function - an MPI function that returns an MPI error code 53330de9b25SBarry Smith 534667f096bSBarry Smith Level: beginner 535667f096bSBarry Smith 5369566063dSJacob Faibussowitsch Notes: 53787497f52SBarry Smith Always returns the error code `PETSC_ERR_MPI`; the MPI error code and string are embedded in 5389566063dSJacob Faibussowitsch the string error message. Do not use this to call any other routines (for example PETSc 5393ba16761SJacob Faibussowitsch routines), it should only be used for direct MPI calls. The user may configure PETSc with the 5403ba16761SJacob Faibussowitsch `--with-strict-petscerrorcode` option to check this at compile-time, otherwise they must 5419566063dSJacob Faibussowitsch check this themselves. 5429566063dSJacob Faibussowitsch 543aaa8cc7dSPierre Jolivet This routine can only be used in functions returning `PetscErrorCode` themselves. If the 5443ba16761SJacob Faibussowitsch calling function returns a different type, use `PetscCallMPIAbort()` instead. 5453ba16761SJacob Faibussowitsch 5469566063dSJacob Faibussowitsch Example Usage: 5479566063dSJacob Faibussowitsch .vb 5489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(...)); // OK, calling MPI function 5499566063dSJacob Faibussowitsch 5509566063dSJacob Faibussowitsch PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead! 5519566063dSJacob Faibussowitsch .ve 5529566063dSJacob Faibussowitsch 55349c86fc7SBarry Smith Fortran Notes: 55487497f52SBarry Smith The Fortran function from which this is used must declare a variable `PetscErrorCode` ierr and ierr must be 55549c86fc7SBarry Smith the final argument to the MPI function being called. 55649c86fc7SBarry Smith 55749c86fc7SBarry Smith In the main program and in Fortran subroutines that do not have ierr as the final return parameter one 55887497f52SBarry Smith should use `PetscCallMPIA()` 55949c86fc7SBarry Smith 56049c86fc7SBarry Smith Fortran Usage: 56149c86fc7SBarry Smith .vb 56249c86fc7SBarry Smith PetscErrorCode ierr or integer ierr 56349c86fc7SBarry Smith ... 56449c86fc7SBarry Smith PetscCallMPI(MPI_Comm_size(...,ierr)) 56549c86fc7SBarry Smith PetscCallMPIA(MPI_Comm_size(...,ierr)) ! Will abort after calling error handler 56649c86fc7SBarry Smith 56749c86fc7SBarry Smith PetscCallMPI(MPI_Comm_size(...,eflag)) ! ERROR, final argument must be ierr 56849c86fc7SBarry Smith .ve 56949c86fc7SBarry Smith 570db781477SPatrick Sanan .seealso: `SETERRMPI()`, `PetscCall()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`, 5713ba16761SJacob Faibussowitsch `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 5723ba16761SJacob Faibussowitsch `PetscError()`, `CHKMEMQ` 5733ba16761SJacob Faibussowitsch M*/ 5743ba16761SJacob Faibussowitsch 5753ba16761SJacob Faibussowitsch /*MC 5763ba16761SJacob Faibussowitsch PetscCallMPIAbort - Like `PetscCallMPI()` but calls `MPI_Abort()` on error 5773ba16761SJacob Faibussowitsch 5783ba16761SJacob Faibussowitsch Synopsis: 5793ba16761SJacob Faibussowitsch #include <petscerror.h> 5803ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm comm, MPI_Function(args)) 5813ba16761SJacob Faibussowitsch 5823ba16761SJacob Faibussowitsch Not Collective 5833ba16761SJacob Faibussowitsch 5843ba16761SJacob Faibussowitsch Input Parameters: 5853ba16761SJacob Faibussowitsch + comm - the MPI communicator to abort on 5863ba16761SJacob Faibussowitsch - MPI_Function - an MPI function that returns an MPI error code 5873ba16761SJacob Faibussowitsch 588667f096bSBarry Smith Level: beginner 589667f096bSBarry Smith 5903ba16761SJacob Faibussowitsch Notes: 5913ba16761SJacob Faibussowitsch Usage is identical to `PetscCallMPI()`. See `PetscCallMPI()` for detailed discussion. 5923ba16761SJacob Faibussowitsch 5933ba16761SJacob Faibussowitsch This routine may be used in functions returning `void` or other non-`PetscErrorCode` types. 5943ba16761SJacob Faibussowitsch 595989c49a2SBarry Smith Fortran Note: 596989c49a2SBarry Smith In Fortran this is called `PetscCallMPIA()` and is intended to be used in the main program while `PetscCallMPI()` is 597989c49a2SBarry Smith used in Fortran subroutines. 598989c49a2SBarry Smith 599989c49a2SBarry Smith Developer Note: 600989c49a2SBarry Smith This should have the same name in Fortran. 601989c49a2SBarry Smith 6023ba16761SJacob Faibussowitsch .seealso: `PetscCallMPI()`, `PetscCallAbort()`, `SETERRABORT()` 60330de9b25SBarry Smith M*/ 6043fcd9f07SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 60563a3b9bcSJacob Faibussowitsch void PetscCallMPI(PetscMPIInt); 6063ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm, PetscMPIInt); 607064a246eSJacob Faibussowitsch #else 6083ba16761SJacob Faibussowitsch #define PetscCallMPI_Private(__PETSC_STACK_POP_FUNC__, __SETERR_FUNC__, __COMM__, ...) \ 6099371c9d4SSatish Balay do { \ 6103ba16761SJacob Faibussowitsch PetscMPIInt ierr_petsc_call_mpi_; \ 611ef1023bdSBarry Smith PetscStackUpdateLine; \ 612792fecdfSBarry Smith PetscStackPushExternal("MPI function"); \ 613d71ae5a4SJacob Faibussowitsch { \ 6143ba16761SJacob Faibussowitsch ierr_petsc_call_mpi_ = __VA_ARGS__; \ 615d71ae5a4SJacob Faibussowitsch } \ 6163ba16761SJacob Faibussowitsch __PETSC_STACK_POP_FUNC__; \ 6173ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_mpi_ != MPI_SUCCESS)) { \ 6183ba16761SJacob Faibussowitsch char petsc_mpi_7_errorstring[2 * MPI_MAX_ERROR_STRING]; \ 6193ba16761SJacob Faibussowitsch PetscMPIErrorString(ierr_petsc_call_mpi_, (char *)petsc_mpi_7_errorstring); \ 6203ba16761SJacob Faibussowitsch __SETERR_FUNC__(__COMM__, PETSC_ERR_MPI, "MPI error %d %s", (int)ierr_petsc_call_mpi_, petsc_mpi_7_errorstring); \ 6213fcd9f07SJacob Faibussowitsch } \ 6223fcd9f07SJacob Faibussowitsch } while (0) 6233ba16761SJacob Faibussowitsch 6243ba16761SJacob Faibussowitsch #define PetscCallMPI(...) PetscCallMPI_Private(PetscStackPop, SETERRQ, PETSC_COMM_SELF, __VA_ARGS__) 6253ba16761SJacob Faibussowitsch #define PetscCallMPIAbort(comm, ...) PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRABORT, comm, __VA_ARGS__) 6269566063dSJacob Faibussowitsch #endif 627064a246eSJacob Faibussowitsch 6287037b0edSPatrick Sanan /*MC 6299566063dSJacob Faibussowitsch CHKERRMPI - Checks error code returned from MPI calls, if non-zero it calls the error 6309566063dSJacob Faibussowitsch handler and then returns 6319566063dSJacob Faibussowitsch 6329566063dSJacob Faibussowitsch Synopsis: 6339566063dSJacob Faibussowitsch #include <petscerror.h> 6349566063dSJacob Faibussowitsch void CHKERRMPI(PetscErrorCode ierr) 6359566063dSJacob Faibussowitsch 6369566063dSJacob Faibussowitsch Not Collective 6379566063dSJacob Faibussowitsch 6389566063dSJacob Faibussowitsch Input Parameter: 6399566063dSJacob Faibussowitsch . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 6409566063dSJacob Faibussowitsch 6419566063dSJacob Faibussowitsch Level: deprecated 6429566063dSJacob Faibussowitsch 643667f096bSBarry Smith Note: 644667f096bSBarry Smith Deprecated in favor of `PetscCallMPI()`. This routine behaves identically to it. 645667f096bSBarry Smith 646db781477SPatrick Sanan .seealso: `PetscCallMPI()` 6479566063dSJacob Faibussowitsch M*/ 6489566063dSJacob Faibussowitsch #define CHKERRMPI(...) PetscCallMPI(__VA_ARGS__) 6499566063dSJacob Faibussowitsch 6509566063dSJacob Faibussowitsch /*MC 651989c49a2SBarry Smith PetscCallAbort - Checks error code returned from PETSc function, if non-zero it aborts immediately by calling `MPI_Abort()` 6529566063dSJacob Faibussowitsch 6539566063dSJacob Faibussowitsch Synopsis: 6549566063dSJacob Faibussowitsch #include <petscerror.h> 6559566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm comm, PetscErrorCode ierr) 6569566063dSJacob Faibussowitsch 657c3339decSBarry Smith Collective 6589566063dSJacob Faibussowitsch 6599566063dSJacob Faibussowitsch Input Parameters: 6609566063dSJacob Faibussowitsch + comm - the MPI communicator on which to abort 6619566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 6629566063dSJacob Faibussowitsch 663667f096bSBarry Smith Level: intermediate 664667f096bSBarry Smith 6659566063dSJacob Faibussowitsch Notes: 66687497f52SBarry Smith This macro has identical type and usage semantics to `PetscCall()` with the important caveat 6679566063dSJacob Faibussowitsch that this macro does not return. Instead, if ierr is nonzero it calls the PETSc error handler 66887497f52SBarry Smith and then immediately calls `MPI_Abort()`. It can therefore be used anywhere. 6699566063dSJacob Faibussowitsch 67087497f52SBarry Smith As per `MPI_Abort()` semantics the communicator passed must be valid, although there is currently 67187497f52SBarry Smith no attempt made at handling any potential errors from `MPI_Abort()`. Note that while 67287497f52SBarry Smith `MPI_Abort()` is required to terminate only those processes which reside on comm, it is often 67387497f52SBarry Smith the case that `MPI_Abort()` terminates *all* processes. 6749566063dSJacob Faibussowitsch 6759566063dSJacob Faibussowitsch Example Usage: 6769566063dSJacob Faibussowitsch .vb 6779566063dSJacob Faibussowitsch PetscErrorCode boom(void) { return PETSC_ERR_MEM; } 6789566063dSJacob Faibussowitsch 6799566063dSJacob Faibussowitsch void foo(void) 6809566063dSJacob Faibussowitsch { 6819566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type 6829566063dSJacob Faibussowitsch } 6839566063dSJacob Faibussowitsch 6849566063dSJacob Faibussowitsch double bar(void) 6859566063dSJacob Faibussowitsch { 6869566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type 6879566063dSJacob Faibussowitsch } 6889566063dSJacob Faibussowitsch 6899566063dSJacob Faibussowitsch PetscCallAbort(MPI_COMM_NULL,boom()); // ERROR, communicator should be valid 6909566063dSJacob Faibussowitsch 6919566063dSJacob Faibussowitsch struct baz 6929566063dSJacob Faibussowitsch { 6939566063dSJacob Faibussowitsch baz() 6949566063dSJacob Faibussowitsch { 6959566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,boom()); // OK 6969566063dSJacob Faibussowitsch } 6979566063dSJacob Faibussowitsch 6989566063dSJacob Faibussowitsch ~baz() 6999566063dSJacob Faibussowitsch { 7009566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,boom()); // OK (in fact the only way to handle PETSc errors) 7019566063dSJacob Faibussowitsch } 7029566063dSJacob Faibussowitsch }; 7039566063dSJacob Faibussowitsch .ve 7049566063dSJacob Faibussowitsch 705989c49a2SBarry Smith Fortran Note: 706989c49a2SBarry Smith Use `PetscCallA()`. 707989c49a2SBarry Smith 708989c49a2SBarry Smith Developer Note: 709989c49a2SBarry Smith This should have the same name in Fortran as in C. 710989c49a2SBarry Smith 711db781477SPatrick Sanan .seealso: `SETERRABORT()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, 7125687f061SJacob Faibussowitsch `SETERRQ()`, `CHKMEMQ`, `PetscCallMPI()`, `PetscCallCXXAbort()` 7139566063dSJacob Faibussowitsch M*/ 7149566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 7159566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm, PetscErrorCode); 7169566063dSJacob Faibussowitsch void PetscCallContinue(PetscErrorCode); 7179566063dSJacob Faibussowitsch #else 7189371c9d4SSatish Balay #define PetscCallAbort(comm, ...) \ 7199371c9d4SSatish Balay do { \ 7207da6b3ccSLisandro Dalcin PetscErrorCode ierr_petsc_call_abort_; \ 7213ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 7227da6b3ccSLisandro Dalcin ierr_petsc_call_abort_ = __VA_ARGS__; \ 7233ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_abort_ != PETSC_SUCCESS)) { \ 7243ba16761SJacob Faibussowitsch ierr_petsc_call_abort_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_abort_, PETSC_ERROR_REPEAT, " "); \ 7253ba16761SJacob Faibussowitsch (void)MPI_Abort(comm, (PetscMPIInt)ierr_petsc_call_abort_); \ 7269566063dSJacob Faibussowitsch } \ 7279566063dSJacob Faibussowitsch } while (0) 7289371c9d4SSatish Balay #define PetscCallContinue(...) \ 7299371c9d4SSatish Balay do { \ 7307da6b3ccSLisandro Dalcin PetscErrorCode ierr_petsc_call_continue_; \ 7313ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 7327da6b3ccSLisandro Dalcin ierr_petsc_call_continue_ = __VA_ARGS__; \ 7333ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_continue_ != PETSC_SUCCESS)) { \ 7343ba16761SJacob Faibussowitsch ierr_petsc_call_continue_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_continue_, PETSC_ERROR_REPEAT, " "); \ 7353ba16761SJacob Faibussowitsch (void)ierr_petsc_call_continue_; \ 7363ba16761SJacob Faibussowitsch } \ 7379566063dSJacob Faibussowitsch } while (0) 7389566063dSJacob Faibussowitsch #endif 7399566063dSJacob Faibussowitsch 7409566063dSJacob Faibussowitsch /*MC 7419566063dSJacob Faibussowitsch CHKERRABORT - Checks error code returned from PETSc function. If non-zero it aborts immediately. 7429566063dSJacob Faibussowitsch 7439566063dSJacob Faibussowitsch Synopsis: 7449566063dSJacob Faibussowitsch #include <petscerror.h> 7459566063dSJacob Faibussowitsch void CHKERRABORT(MPI_Comm comm, PetscErrorCode ierr) 7469566063dSJacob Faibussowitsch 7479566063dSJacob Faibussowitsch Not Collective 7489566063dSJacob Faibussowitsch 7499566063dSJacob Faibussowitsch Input Parameters: 7509566063dSJacob Faibussowitsch + comm - the MPI communicator 7519566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 7529566063dSJacob Faibussowitsch 7539566063dSJacob Faibussowitsch Level: deprecated 7549566063dSJacob Faibussowitsch 755667f096bSBarry Smith Note: 756667f096bSBarry Smith Deprecated in favor of `PetscCallAbort()`. This routine behaves identically to it. 757667f096bSBarry Smith 758af27ebaaSBarry Smith .seealso: `PetscCallAbort()`, `PetscErrorCode` 7599566063dSJacob Faibussowitsch M*/ 7609566063dSJacob Faibussowitsch #define CHKERRABORT(comm, ...) PetscCallAbort(comm, __VA_ARGS__) 7619566063dSJacob Faibussowitsch #define CHKERRCONTINUE(...) PetscCallContinue(__VA_ARGS__) 7629566063dSJacob Faibussowitsch 7639566063dSJacob Faibussowitsch /*MC 76487497f52SBarry Smith CHKERRA - Fortran-only replacement for use of `CHKERRQ()` in the main program, which aborts immediately 765f388eb8bSPatrick Sanan 766f388eb8bSPatrick Sanan Synopsis: 767f388eb8bSPatrick Sanan #include <petscsys.h> 768f388eb8bSPatrick Sanan PetscErrorCode CHKERRA(PetscErrorCode ierr) 769f388eb8bSPatrick Sanan 770f388eb8bSPatrick Sanan Not Collective 771f388eb8bSPatrick Sanan 7722fe279fdSBarry Smith Input Parameter: 773f388eb8bSPatrick Sanan . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 774f388eb8bSPatrick Sanan 77587497f52SBarry Smith Level: deprecated 776f388eb8bSPatrick Sanan 77787497f52SBarry Smith Note: 77887497f52SBarry Smith This macro is rarely needed, normal usage is `PetscCallA()` in the main Fortran program. 779f388eb8bSPatrick Sanan 780989c49a2SBarry Smith Developer Note: 781989c49a2SBarry Smith Why isn't this named `CHKERRABORT()` in Fortran? 782989c49a2SBarry Smith 78387497f52SBarry Smith .seealso: `PetscCall()`, `PetscCallA()`, `PetscCallAbort()`, `CHKERRQ()`, `SETERRA()`, `SETERRQ()`, `SETERRABORT()` 784f388eb8bSPatrick Sanan M*/ 785f388eb8bSPatrick Sanan 7861e4f893aSSatish Balay PETSC_EXTERN PetscBool petscwaitonerrorflg; 7871e4f893aSSatish Balay PETSC_EXTERN PetscBool petscindebugger; 78835f00c14SToby Isaac PETSC_EXTERN PetscBool petscabortmpifinalize; 7897c66cc67SJunchao Zhang 7907c66cc67SJunchao Zhang /*MC 791667f096bSBarry Smith PETSCABORT - Call `MPI_Abort()` with an informative error code 7927c66cc67SJunchao Zhang 7937c66cc67SJunchao Zhang Synopsis: 7947c66cc67SJunchao Zhang #include <petscsys.h> 7957c66cc67SJunchao Zhang PETSCABORT(MPI_Comm comm, PetscErrorCode ierr) 7967c66cc67SJunchao Zhang 7977cdbe19fSJose E. Roman Collective; No Fortran Support 7987c66cc67SJunchao Zhang 7997c66cc67SJunchao Zhang Input Parameters: 80095bd0b28SBarry Smith + comm - An MPI communicator, so that the error can be collective 8017c66cc67SJunchao Zhang - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 8027c66cc67SJunchao Zhang 803bf4d2887SBarry Smith Level: advanced 8047c66cc67SJunchao Zhang 805bf4d2887SBarry Smith Notes: 806989c49a2SBarry Smith If the option `-start_in_debugger` was used then this calls `abort()` to stop the program in the debugger. 807bf4d2887SBarry Smith 808989c49a2SBarry Smith if `PetscCIEnabledPortableErrorOutput` is set, which means the code is running in the PETSc test harness (make test), 809989c49a2SBarry Smith and `comm` is `MPI_COMM_WORLD` it strives to exit cleanly without calling `MPI_Abort()` and instead calling `MPI_Finalize()`. 810660278c0SBarry Smith 811989c49a2SBarry Smith This is currently only used when an error propagates up to the C `main()` program and is detected by a `PetscCall()`, `PetscCallMPI()`, 812989c49a2SBarry Smith or is set in `main()` with `SETERRQ()`. Abort calls such as `SETERRABORT()`, 813989c49a2SBarry Smith `PetscCheckAbort()`, `PetscCallMPIAbort()`, and `PetscCallAbort()` always call `MPI_Abort()` and do not have any special 814989c49a2SBarry Smith handling for the test harness. 815989c49a2SBarry Smith 816989c49a2SBarry Smith Developer Note: 817989c49a2SBarry Smith Should the other abort calls also pass through this call instead of calling `MPI_Abort()` directly? 818989c49a2SBarry Smith 819989c49a2SBarry Smith .seealso: `PetscError()`, `PetscCall()`, `SETERRABORT()`, `PetscCheckAbort()`, `PetscCallMPIAbort()`, `PetscCall()`, `PetscCallMPI()`, 820af27ebaaSBarry Smith `PetscCallAbort()`, `MPI_Abort()`, `PetscErrorCode` 8217c66cc67SJunchao Zhang M*/ 82229f88feaSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 82329f88feaSJacob Faibussowitsch void PETSCABORT(MPI_Comm, PetscErrorCode); 82429f88feaSJacob Faibussowitsch #else 8259371c9d4SSatish Balay #define PETSCABORT(comm, ...) \ 8269371c9d4SSatish Balay do { \ 8273ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_abort_; \ 8283ba16761SJacob Faibussowitsch if (petscwaitonerrorflg) { ierr_petsc_abort_ = PetscSleep(1000); } \ 8293ba16761SJacob Faibussowitsch if (petscindebugger) { \ 8303ba16761SJacob Faibussowitsch abort(); \ 8313ba16761SJacob Faibussowitsch } else { \ 8327da6b3ccSLisandro Dalcin PetscMPIInt size_; \ 8333ba16761SJacob Faibussowitsch ierr_petsc_abort_ = __VA_ARGS__; \ 8347da6b3ccSLisandro Dalcin MPI_Comm_size(comm, &size_); \ 83535f00c14SToby Isaac if (PetscCIEnabledPortableErrorOutput && (size_ == PetscGlobalSize || petscabortmpifinalize) && ierr_petsc_abort_ != PETSC_ERR_SIG) { \ 8369371c9d4SSatish Balay MPI_Finalize(); \ 8379371c9d4SSatish Balay exit(0); \ 838660278c0SBarry Smith } else if (PetscCIEnabledPortableErrorOutput && PetscGlobalSize == 1) { \ 839660278c0SBarry Smith exit(0); \ 840660278c0SBarry Smith } else { \ 841660278c0SBarry Smith MPI_Abort(comm, (PetscMPIInt)ierr_petsc_abort_); \ 842660278c0SBarry Smith } \ 8433fcd9f07SJacob Faibussowitsch } \ 8447c66cc67SJunchao Zhang } while (0) 84529f88feaSJacob Faibussowitsch #endif 846986eef2eSBarry Smith 8479566063dSJacob Faibussowitsch #ifdef PETSC_CLANGUAGE_CXX 848986eef2eSBarry Smith /*MC 8499566063dSJacob Faibussowitsch PetscCallThrow - Checks error code, if non-zero it calls the C++ error handler which throws 8509566063dSJacob Faibussowitsch an exception 851986eef2eSBarry Smith 852986eef2eSBarry Smith Synopsis: 8539566063dSJacob Faibussowitsch #include <petscerror.h> 8549566063dSJacob Faibussowitsch void PetscCallThrow(PetscErrorCode ierr) 855986eef2eSBarry Smith 856986eef2eSBarry Smith Not Collective 857986eef2eSBarry Smith 8589566063dSJacob Faibussowitsch Input Parameter: 859986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 860986eef2eSBarry Smith 861667f096bSBarry Smith Level: beginner 862667f096bSBarry Smith 863986eef2eSBarry Smith Notes: 864af27ebaaSBarry Smith Requires PETSc to be configured with clanguage of c++. Throws a std::runtime_error() on error. 865986eef2eSBarry Smith 86687497f52SBarry Smith Once the error handler throws the exception you can use `PetscCallVoid()` which returns without 86787497f52SBarry Smith an error code (bad idea since the error is ignored) or `PetscCallAbort()` to have `MPI_Abort()` 8689566063dSJacob Faibussowitsch called immediately. 8699566063dSJacob Faibussowitsch 870db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscCall()`, `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, 871db781477SPatrick Sanan `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ` 872986eef2eSBarry Smith M*/ 8739371c9d4SSatish Balay #define PetscCallThrow(...) \ 8749371c9d4SSatish Balay do { \ 8753ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 8763ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_throw_ = __VA_ARGS__; \ 8773ba16761SJacob 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); \ 878ffc4695bSBarry Smith } while (0) 87985614651SBarry Smith 880cc26af49SMatthew Knepley /*MC 881cc26af49SMatthew Knepley CHKERRXX - Checks error code, if non-zero it calls the C++ error handler which throws an exception 882cc26af49SMatthew Knepley 883cc26af49SMatthew Knepley Synopsis: 8849566063dSJacob Faibussowitsch #include <petscerror.h> 8853af045c5SBarry Smith void CHKERRXX(PetscErrorCode ierr) 886cc26af49SMatthew Knepley 887eca87e8dSBarry Smith Not Collective 888cc26af49SMatthew Knepley 8899566063dSJacob Faibussowitsch Input Parameter: 8903af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 891cc26af49SMatthew Knepley 8929566063dSJacob Faibussowitsch Level: deprecated 893cc26af49SMatthew Knepley 894667f096bSBarry Smith Note: 895667f096bSBarry Smith Deprecated in favor of `PetscCallThrow()`. This routine behaves identically to it. 896667f096bSBarry Smith 897db781477SPatrick Sanan .seealso: `PetscCallThrow()` 898cc26af49SMatthew Knepley M*/ 8999566063dSJacob Faibussowitsch #define CHKERRXX(...) PetscCallThrow(__VA_ARGS__) 900fd705b32SMatthew Knepley #endif 901fd705b32SMatthew Knepley 9023ba16761SJacob Faibussowitsch #define PetscCallCXX_Private(__SETERR_FUNC__, __COMM__, ...) \ 9033ba16761SJacob Faibussowitsch do { \ 9043ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 9053ba16761SJacob Faibussowitsch try { \ 9063ba16761SJacob Faibussowitsch __VA_ARGS__; \ 9073ba16761SJacob Faibussowitsch } catch (const std::exception &e) { \ 9083ba16761SJacob Faibussowitsch __SETERR_FUNC__(__COMM__, PETSC_ERR_LIB, "%s", e.what()); \ 9093ba16761SJacob Faibussowitsch } \ 9103ba16761SJacob Faibussowitsch } while (0) 9113ba16761SJacob Faibussowitsch 9123f520e80SJunchao Zhang /*MC 9139566063dSJacob Faibussowitsch PetscCallCXX - Checks C++ function calls and if they throw an exception, catch it and then 9149566063dSJacob Faibussowitsch return a PETSc error code 9153f520e80SJunchao Zhang 9163f520e80SJunchao Zhang Synopsis: 9179566063dSJacob Faibussowitsch #include <petscerror.h> 9185687f061SJacob Faibussowitsch void PetscCallCXX(...) noexcept; 9193f520e80SJunchao Zhang 9203f520e80SJunchao Zhang Not Collective 9213f520e80SJunchao Zhang 9229566063dSJacob Faibussowitsch Input Parameter: 9235687f061SJacob Faibussowitsch . __VA_ARGS__ - An arbitrary expression 9245687f061SJacob Faibussowitsch 9255687f061SJacob Faibussowitsch Level: beginner 9263f520e80SJunchao Zhang 9273f520e80SJunchao Zhang Notes: 9285687f061SJacob Faibussowitsch `PetscCallCXX(...)` is a macro replacement for 9299566063dSJacob Faibussowitsch .vb 9309566063dSJacob Faibussowitsch try { 9315687f061SJacob Faibussowitsch __VA_ARGS__; 9329566063dSJacob Faibussowitsch } catch (const std::exception& e) { 9339566063dSJacob Faibussowitsch return ConvertToPetscErrorCode(e); 9349566063dSJacob Faibussowitsch } 9359566063dSJacob Faibussowitsch .ve 9369566063dSJacob Faibussowitsch Due to the fact that it catches any (reasonable) exception, it is essentially noexcept. 9373f520e80SJunchao Zhang 9385687f061SJacob Faibussowitsch If you cannot return a `PetscErrorCode` use `PetscCallCXXAbort()` instead. 9395687f061SJacob Faibussowitsch 9409566063dSJacob Faibussowitsch Example Usage: 9419566063dSJacob Faibussowitsch .vb 9429566063dSJacob Faibussowitsch void foo(void) { throw std::runtime_error("error"); } 9433f520e80SJunchao Zhang 9449566063dSJacob Faibussowitsch void bar() 9459566063dSJacob Faibussowitsch { 946e8952933SJacob Faibussowitsch PetscCallCXX(foo()); // ERROR bar() does not return PetscErrorCode 9479566063dSJacob Faibussowitsch } 9489566063dSJacob Faibussowitsch 9499566063dSJacob Faibussowitsch PetscErrorCode baz() 9509566063dSJacob Faibussowitsch { 9519566063dSJacob Faibussowitsch PetscCallCXX(foo()); // OK 9529566063dSJacob Faibussowitsch 9539566063dSJacob Faibussowitsch PetscCallCXX( 9549566063dSJacob Faibussowitsch bar(); 955d5b43468SJose E. Roman foo(); // OK multiple statements allowed 9569566063dSJacob Faibussowitsch ); 9579566063dSJacob Faibussowitsch } 9589566063dSJacob Faibussowitsch 9599566063dSJacob Faibussowitsch struct bop 9609566063dSJacob Faibussowitsch { 9619566063dSJacob Faibussowitsch bop() 9629566063dSJacob Faibussowitsch { 963e8952933SJacob Faibussowitsch PetscCallCXX(foo()); // ERROR returns PetscErrorCode, cannot be used in constructors 9649566063dSJacob Faibussowitsch } 9659566063dSJacob Faibussowitsch }; 9669566063dSJacob Faibussowitsch 967e8952933SJacob Faibussowitsch // ERROR contains do-while, cannot be used as function-try block 9689566063dSJacob Faibussowitsch PetscErrorCode qux() PetscCallCXX( 9699566063dSJacob Faibussowitsch bar(); 9709566063dSJacob Faibussowitsch baz(); 9719566063dSJacob Faibussowitsch foo(); 9729566063dSJacob Faibussowitsch return 0; 9739566063dSJacob Faibussowitsch ) 9749566063dSJacob Faibussowitsch .ve 9759566063dSJacob Faibussowitsch 9765687f061SJacob Faibussowitsch .seealso: `PetscCallCXXAbort()`, `PetscCallThrow()`, `SETERRQ()`, `PetscCall()`, 9775687f061SJacob Faibussowitsch `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 9785687f061SJacob Faibussowitsch `PetscError()`, `CHKMEMQ` 9793f520e80SJunchao Zhang M*/ 9803ba16761SJacob Faibussowitsch #define PetscCallCXX(...) PetscCallCXX_Private(SETERRQ, PETSC_COMM_SELF, __VA_ARGS__) 9815687f061SJacob Faibussowitsch 9825687f061SJacob Faibussowitsch /*MC 9835687f061SJacob Faibussowitsch PetscCallCXXAbort - Like `PetscCallCXX()` but calls `MPI_Abort()` instead of returning an 9845687f061SJacob Faibussowitsch error-code 9855687f061SJacob Faibussowitsch 9865687f061SJacob Faibussowitsch Synopsis: 9875687f061SJacob Faibussowitsch #include <petscerror.h> 9885687f061SJacob Faibussowitsch void PetscCallCXXAbort(MPI_Comm comm, ...) noexcept; 9895687f061SJacob Faibussowitsch 99020f4b53cSBarry Smith Collective; No Fortran Support 9915687f061SJacob Faibussowitsch 9925687f061SJacob Faibussowitsch Input Parameters: 9935687f061SJacob Faibussowitsch + comm - The MPI communicator to abort on 9945687f061SJacob Faibussowitsch - __VA_ARGS__ - An arbitrary expression 9955687f061SJacob Faibussowitsch 9965687f061SJacob Faibussowitsch Level: beginner 9975687f061SJacob Faibussowitsch 9985687f061SJacob Faibussowitsch Notes: 9995687f061SJacob Faibussowitsch This macro may be used to check C++ expressions for exceptions in cases where you cannot 10005687f061SJacob Faibussowitsch return an error code. This includes constructors, destructors, copy/move assignment functions 10015687f061SJacob Faibussowitsch or constructors among others. 10025687f061SJacob Faibussowitsch 10035687f061SJacob Faibussowitsch If an exception is caught, the macro calls `SETERRABORT()` on `comm`. The exception must 10045687f061SJacob Faibussowitsch derive from `std::exception` in order to be caught. 10055687f061SJacob Faibussowitsch 10065687f061SJacob Faibussowitsch If the routine _can_ return an error-code it is highly advised to use `PetscCallCXX()` 10075687f061SJacob Faibussowitsch instead. 10085687f061SJacob Faibussowitsch 10095687f061SJacob Faibussowitsch See `PetscCallCXX()` for additional discussion. 10105687f061SJacob Faibussowitsch 10115687f061SJacob Faibussowitsch Example Usage: 10125687f061SJacob Faibussowitsch .vb 10135687f061SJacob Faibussowitsch class Foo 10145687f061SJacob Faibussowitsch { 10155687f061SJacob Faibussowitsch std::vector<int> data_; 10165687f061SJacob Faibussowitsch 10175687f061SJacob Faibussowitsch public: 10185687f061SJacob Faibussowitsch // normally std::vector::reserve() may raise an exception, but since we handle it with 10195687f061SJacob Faibussowitsch // PetscCallCXXAbort() we may mark this routine as noexcept! 10205687f061SJacob Faibussowitsch Foo() noexcept 10215687f061SJacob Faibussowitsch { 10225687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, data_.reserve(10)); 10235687f061SJacob Faibussowitsch } 10245687f061SJacob Faibussowitsch }; 10255687f061SJacob Faibussowitsch 10265687f061SJacob Faibussowitsch std::vector<int> bar() 10275687f061SJacob Faibussowitsch { 10285687f061SJacob Faibussowitsch std::vector<int> v; 10295687f061SJacob Faibussowitsch 10305687f061SJacob Faibussowitsch PetscFunctionBegin; 10315687f061SJacob Faibussowitsch // OK! 10325687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1)); 10335687f061SJacob Faibussowitsch PetscFunctionReturn(v); 10345687f061SJacob Faibussowitsch } 10355687f061SJacob Faibussowitsch 10365687f061SJacob Faibussowitsch PetscErrorCode baz() 10375687f061SJacob Faibussowitsch { 10385687f061SJacob Faibussowitsch std::vector<int> v; 10395687f061SJacob Faibussowitsch 10405687f061SJacob Faibussowitsch PetscFunctionBegin; 10415687f061SJacob Faibussowitsch // WRONG! baz() returns a PetscErrorCode, prefer PetscCallCXX() instead 10425687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1)); 10433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10445687f061SJacob Faibussowitsch } 10455687f061SJacob Faibussowitsch .ve 10465687f061SJacob Faibussowitsch 10475687f061SJacob Faibussowitsch .seealso: `PetscCallCXX()`, `SETERRABORT()`, `PetscCallAbort()` 10485687f061SJacob Faibussowitsch M*/ 10493ba16761SJacob Faibussowitsch #define PetscCallCXXAbort(comm, ...) PetscCallCXX_Private(SETERRABORT, comm, __VA_ARGS__) 10503f520e80SJunchao Zhang 105130de9b25SBarry Smith /*MC 10529566063dSJacob Faibussowitsch CHKERRCXX - Checks C++ function calls and if they throw an exception, catch it and then 10539566063dSJacob Faibussowitsch return a PETSc error code 10549566063dSJacob Faibussowitsch 10559566063dSJacob Faibussowitsch Synopsis: 10569566063dSJacob Faibussowitsch #include <petscerror.h> 10579566063dSJacob Faibussowitsch void CHKERRCXX(func) noexcept; 10589566063dSJacob Faibussowitsch 10599566063dSJacob Faibussowitsch Not Collective 10609566063dSJacob Faibussowitsch 10619566063dSJacob Faibussowitsch Input Parameter: 10629566063dSJacob Faibussowitsch . func - C++ function calls 10639566063dSJacob Faibussowitsch 10649566063dSJacob Faibussowitsch Level: deprecated 10659566063dSJacob Faibussowitsch 1066667f096bSBarry Smith Note: 1067667f096bSBarry Smith Deprecated in favor of `PetscCallCXX()`. This routine behaves identically to it. 1068667f096bSBarry Smith 1069db781477SPatrick Sanan .seealso: `PetscCallCXX()` 10709566063dSJacob Faibussowitsch M*/ 10719566063dSJacob Faibussowitsch #define CHKERRCXX(...) PetscCallCXX(__VA_ARGS__) 10729566063dSJacob Faibussowitsch 10739566063dSJacob Faibussowitsch /*MC 107430de9b25SBarry Smith CHKMEMQ - Checks the memory for corruption, calls error handler if any is detected 107530de9b25SBarry Smith 107630de9b25SBarry Smith Synopsis: 1077aaa7dc30SBarry Smith #include <petscsys.h> 107891d3bdf4SKris Buschelman CHKMEMQ; 107930de9b25SBarry Smith 1080eca87e8dSBarry Smith Not Collective 1081eca87e8dSBarry Smith 108230de9b25SBarry Smith Level: beginner 108330de9b25SBarry Smith 108430de9b25SBarry Smith Notes: 1085af27ebaaSBarry Smith We recommend using Valgrind <https://petsc.org/release/faq/#valgrind> or for NVIDIA CUDA systems 1086af27ebaaSBarry Smith <https://docs.nvidia.com/cuda/cuda-memcheck/index.html> for finding memory problems. The ``CHKMEMQ`` macro is useful on systems that 10875ed36255SBarry Smith do not have valgrind, but is not as good as valgrind or cuda-memcheck. 10881957e957SBarry Smith 1089667f096bSBarry Smith Must run with the option `-malloc_debug` (`-malloc_test` in debug mode; or if `PetscMallocSetDebug()` called) to enable this option 109030de9b25SBarry Smith 109130de9b25SBarry Smith Once the error handler is called the calling function is then returned from with the given error code. 109230de9b25SBarry Smith 109330de9b25SBarry Smith By defaults prints location where memory that is corrupted was allocated. 109430de9b25SBarry Smith 10958af9f190SBarry Smith Use `CHKMEMA` for functions that return `void` 1096f621e05eSBarry Smith 1097db781477SPatrick Sanan .seealso: `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `SETERRQ()`, `PetscMallocValidate()` 109830de9b25SBarry Smith M*/ 10996d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 1100064a246eSJacob Faibussowitsch #define CHKMEMQ 1101064a246eSJacob Faibussowitsch #define CHKMEMA 11026d210af2SJacob Faibussowitsch #else 11039371c9d4SSatish Balay #define CHKMEMQ \ 11049371c9d4SSatish Balay do { \ 11053ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_memq_ = PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__); \ 11063ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_memq_ != PETSC_SUCCESS)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_memq_, PETSC_ERROR_REPEAT, " "); \ 110786d09637SLisandro Dalcin } while (0) 11086d210af2SJacob Faibussowitsch #define CHKMEMA PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__) 1109064a246eSJacob Faibussowitsch #endif 11109566063dSJacob Faibussowitsch 1111668f157eSBarry Smith /*E 1112668f157eSBarry Smith PetscErrorType - passed to the PETSc error handling routines indicating if this is the first or a later call to the error handlers 1113668f157eSBarry Smith 1114668f157eSBarry Smith Level: advanced 1115668f157eSBarry Smith 1116667f096bSBarry Smith Note: 111787497f52SBarry Smith `PETSC_ERROR_IN_CXX` indicates the error was detected in C++ and an exception should be generated 1118d736bfebSBarry Smith 1119667f096bSBarry Smith Developer Note: 1120667f096bSBarry Smith This is currently used to decide when to print the detailed information about the run in `PetscTraceBackErrorHandler()` 1121668f157eSBarry Smith 112287497f52SBarry Smith .seealso: `PetscError()`, `SETERRQ()` 1123668f157eSBarry Smith E*/ 11249371c9d4SSatish Balay typedef enum { 11259371c9d4SSatish Balay PETSC_ERROR_INITIAL = 0, 11269371c9d4SSatish Balay PETSC_ERROR_REPEAT = 1, 11279371c9d4SSatish Balay PETSC_ERROR_IN_CXX = 2 11289371c9d4SSatish Balay } PetscErrorType; 11294b209cf6SBarry Smith 1130eb9e708aSLisandro Dalcin #if defined(__clang_analyzer__) 1131eb9e708aSLisandro Dalcin __attribute__((analyzer_noreturn)) 1132eb9e708aSLisandro Dalcin #endif 11339371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode 11349371c9d4SSatish Balay PetscError(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, ...) PETSC_ATTRIBUTE_COLD PETSC_ATTRIBUTE_FORMAT(7, 8); 1135eb9e708aSLisandro Dalcin 1136014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscErrorPrintfInitialize(void); 11373ba16761SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscErrorMessage(PetscErrorCode, const char *[], char **); 1138d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscTraceBackErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1139d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscIgnoreErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1140d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscEmacsClientErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1141d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscMPIAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1142d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1143d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAttachDebuggerErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1144d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscReturnErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1145efca3c55SSatish Balay PETSC_EXTERN PetscErrorCode PetscPushErrorHandler(PetscErrorCode (*handler)(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *), void *); 1146014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopErrorHandler(void); 11478d359177SBarry Smith PETSC_EXTERN PetscErrorCode PetscSignalHandlerDefault(int, void *); 1148014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPushSignalHandler(PetscErrorCode (*)(int, void *), void *); 1149014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopSignalHandler(void); 115028559dc8SJed Brown PETSC_EXTERN PetscErrorCode PetscCheckPointerSetIntensity(PetscInt); 1151c2a741eeSJunchao Zhang PETSC_EXTERN void PetscSignalSegvCheckPointerOrMpi(void); 1152edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 13, 0, "PetscSignalSegvCheckPointerOrMpi()", ) static inline void PetscSignalSegvCheckPointer(void) 1153d71ae5a4SJacob Faibussowitsch { 11549371c9d4SSatish Balay PetscSignalSegvCheckPointerOrMpi(); 11559371c9d4SSatish Balay } 1156329f5518SBarry Smith 1157639ff905SBarry Smith /*MC 1158639ff905SBarry Smith PetscErrorPrintf - Prints error messages. 1159639ff905SBarry Smith 1160639ff905SBarry Smith Synopsis: 1161aaa7dc30SBarry Smith #include <petscsys.h> 1162639ff905SBarry Smith PetscErrorCode (*PetscErrorPrintf)(const char format[], ...); 1163639ff905SBarry Smith 11647cdbe19fSJose E. Roman Not Collective; No Fortran Support 11657cdbe19fSJose E. Roman 1166f899ff85SJose E. Roman Input Parameter: 1167cd05f99aSJacob Faibussowitsch . format - the usual `printf()` format string 1168639ff905SBarry Smith 1169639ff905SBarry Smith Options Database Keys: 11701957e957SBarry Smith + -error_output_stdout - cause error messages to be printed to stdout instead of the (default) stderr 1171e1bc860dSBarry Smith - -error_output_none - to turn off all printing of error messages (does not change the way the error is handled.) 1172639ff905SBarry Smith 1173cf53795eSBarry Smith Level: developer 1174cf53795eSBarry Smith 117595452b02SPatrick Sanan Notes: 117695452b02SPatrick Sanan Use 1177667f096bSBarry Smith .vb 117895bd0b28SBarry Smith PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the error is handled) and 1179667f096bSBarry Smith PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on or you can use your own function 1180667f096bSBarry Smith .ve 1181639ff905SBarry Smith Use 1182667f096bSBarry Smith .vb 118387497f52SBarry Smith `PETSC_STDERR` = FILE* obtained from a file open etc. to have stderr printed to the file. 118487497f52SBarry Smith `PETSC_STDOUT` = FILE* obtained from a file open etc. to have stdout printed to the file. 1185667f096bSBarry Smith .ve 1186639ff905SBarry Smith Use 1187af27ebaaSBarry Smith .vb 118887497f52SBarry Smith `PetscPushErrorHandler()` to provide your own error handler that determines what kind of messages to print 1189af27ebaaSBarry Smith .ve 1190639ff905SBarry Smith 1191db781477SPatrick Sanan .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscHelpPrintf()`, `PetscPrintf()`, `PetscPushErrorHandler()`, `PetscVFPrintf()`, `PetscHelpPrintf()` 1192639ff905SBarry Smith M*/ 11933ca90d2dSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode (*PetscErrorPrintf)(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2); 1194639ff905SBarry Smith 1195cf0818bdSBarry Smith /*E 1196cf0818bdSBarry Smith PetscFPTrap - types of floating point exceptions that may be trapped 1197cf0818bdSBarry Smith 119887497f52SBarry Smith Currently only `PETSC_FP_TRAP_OFF` and `PETSC_FP_TRAP_ON` are handled. All others are treated as `PETSC_FP_TRAP_ON`. 1199cf0818bdSBarry Smith 1200cf0818bdSBarry Smith Level: intermediate 1201cf0818bdSBarry Smith 12027de69702SBarry Smith .seealso: `PetscSetFPTrap()`, `PetscFPTrapPush()` 1203cf0818bdSBarry Smith E*/ 12049371c9d4SSatish Balay typedef enum { 12059371c9d4SSatish Balay PETSC_FP_TRAP_OFF = 0, 12069371c9d4SSatish Balay PETSC_FP_TRAP_INDIV = 1, 12079371c9d4SSatish Balay PETSC_FP_TRAP_FLTOPERR = 2, 12089371c9d4SSatish Balay PETSC_FP_TRAP_FLTOVF = 4, 12099371c9d4SSatish Balay PETSC_FP_TRAP_FLTUND = 8, 12109371c9d4SSatish Balay PETSC_FP_TRAP_FLTDIV = 16, 12119371c9d4SSatish Balay PETSC_FP_TRAP_FLTINEX = 32 12129371c9d4SSatish Balay } PetscFPTrap; 1213bd2b07b1SBarry 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) 1214014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSetFPTrap(PetscFPTrap); 1215014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPush(PetscFPTrap); 1216014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPop(void); 1217aba4c478SBarry Smith PETSC_EXTERN PetscErrorCode PetscDetermineInitialFPTrap(void); 121854a8ef01SBarry Smith 12193a40ed3dSBarry Smith /* 12203a40ed3dSBarry Smith Allows the code to build a stack frame as it runs 12213a40ed3dSBarry Smith */ 12223a40ed3dSBarry Smith 122399cd645aSJed Brown #define PETSCSTACKSIZE 64 12243a40ed3dSBarry Smith typedef struct { 12250e33f6ddSBarry Smith const char *function[PETSCSTACKSIZE]; 12260e33f6ddSBarry Smith const char *file[PETSCSTACKSIZE]; 1227184914b5SBarry Smith int line[PETSCSTACKSIZE]; 1228362febeeSStefano Zampini int petscroutine[PETSCSTACKSIZE]; /* 0 external called from petsc, 1 petsc functions, 2 petsc user functions */ 1229184914b5SBarry Smith int currentsize; 1230a2f94806SJed Brown int hotdepth; 12314be741a6SBarry Smith PetscBool check; /* option to check for correct Push/Pop semantics, true for default petscstack but not other stacks */ 12323a40ed3dSBarry Smith } PetscStack; 1233dfb7d7afSStefano Zampini #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY) 123427104ee2SJacob Faibussowitsch PETSC_EXTERN PetscStack petscstack; 123527104ee2SJacob Faibussowitsch #endif 12363a40ed3dSBarry Smith 12375d12eec7SSatish Balay #if defined(PETSC_SERIALIZE_FUNCTIONS) 12385d12eec7SSatish Balay #include <petsc/private/petscfptimpl.h> 12395d12eec7SSatish Balay /* 12405d12eec7SSatish Balay Registers the current function into the global function pointer to function name table 12415d12eec7SSatish Balay 12425d12eec7SSatish Balay Have to fix this to handle errors but cannot return error since used in PETSC_VIEWER_DRAW_() etc 12435d12eec7SSatish Balay */ 12449371c9d4SSatish Balay #define PetscRegister__FUNCT__() \ 12459371c9d4SSatish Balay do { \ 12465d12eec7SSatish Balay static PetscBool __chked = PETSC_FALSE; \ 12475d12eec7SSatish Balay if (!__chked) { \ 12489371c9d4SSatish Balay void *ptr; \ 12493ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscDLSym(NULL, PETSC_FUNCTION_NAME, &ptr)); \ 12505d12eec7SSatish Balay __chked = PETSC_TRUE; \ 12519371c9d4SSatish Balay } \ 12529371c9d4SSatish Balay } while (0) 12535d12eec7SSatish Balay #else 12545d12eec7SSatish Balay #define PetscRegister__FUNCT__() 12555d12eec7SSatish Balay #endif 12565d12eec7SSatish Balay 1257eae3dc7dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) || defined(__clang_analyzer__) 12586d210af2SJacob Faibussowitsch #define PetscStackPushNoCheck(funct, petsc_routine, hot) 125937154d10SBarry Smith #define PetscStackUpdateLine 1260792fecdfSBarry Smith #define PetscStackPushExternal(funct) 12616d210af2SJacob Faibussowitsch #define PetscStackPopNoCheck 12626d210af2SJacob Faibussowitsch #define PetscStackClearTop 12636d210af2SJacob Faibussowitsch #define PetscFunctionBegin 12646d210af2SJacob Faibussowitsch #define PetscFunctionBeginUser 12656d210af2SJacob Faibussowitsch #define PetscFunctionBeginHot 12660a57284eSJacob Faibussowitsch #define PetscFunctionReturn(...) return __VA_ARGS__ 12676d210af2SJacob Faibussowitsch #define PetscFunctionReturnVoid() return 12686d210af2SJacob Faibussowitsch #define PetscStackPop 12696d210af2SJacob Faibussowitsch #define PetscStackPush(f) 1270dfb7d7afSStefano Zampini #elif defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY) 1271660278c0SBarry Smith 12729371c9d4SSatish Balay #define PetscStackPush_Private(stack__, file__, func__, line__, petsc_routine__, hot__) \ 12739371c9d4SSatish Balay do { \ 12745a96b57dSJacob Faibussowitsch if (stack__.currentsize < PETSCSTACKSIZE) { \ 12755a96b57dSJacob Faibussowitsch stack__.function[stack__.currentsize] = func__; \ 1276ef1023bdSBarry Smith if (petsc_routine__) { \ 1277ef1023bdSBarry Smith stack__.file[stack__.currentsize] = file__; \ 12785a96b57dSJacob Faibussowitsch stack__.line[stack__.currentsize] = line__; \ 1279ef1023bdSBarry Smith } else { \ 1280648bc8c4SBarry Smith stack__.file[stack__.currentsize] = PETSC_NULLPTR; \ 1281ef1023bdSBarry Smith stack__.line[stack__.currentsize] = 0; \ 1282ef1023bdSBarry Smith } \ 12835a96b57dSJacob Faibussowitsch stack__.petscroutine[stack__.currentsize] = petsc_routine__; \ 12845a96b57dSJacob Faibussowitsch } \ 12855a96b57dSJacob Faibussowitsch ++stack__.currentsize; \ 12865a96b57dSJacob Faibussowitsch stack__.hotdepth += (hot__ || stack__.hotdepth); \ 12875a96b57dSJacob Faibussowitsch } while (0) 12885a96b57dSJacob Faibussowitsch 12894be741a6SBarry Smith /* uses PetscCheckAbort() because may be used in a function that does not return an error code */ 12909371c9d4SSatish Balay #define PetscStackPop_Private(stack__, func__) \ 12919371c9d4SSatish Balay do { \ 12924be741a6SBarry 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__); \ 12935a96b57dSJacob Faibussowitsch if (--stack__.currentsize < PETSCSTACKSIZE) { \ 12949371c9d4SSatish 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", \ 12959371c9d4SSatish Balay stack__.function[stack__.currentsize], stack__.file[stack__.currentsize], stack__.line[stack__.currentsize], func__, __FILE__, __LINE__); \ 12965a96b57dSJacob Faibussowitsch stack__.function[stack__.currentsize] = PETSC_NULLPTR; \ 12975a96b57dSJacob Faibussowitsch stack__.file[stack__.currentsize] = PETSC_NULLPTR; \ 12985a96b57dSJacob Faibussowitsch stack__.line[stack__.currentsize] = 0; \ 12995a96b57dSJacob Faibussowitsch stack__.petscroutine[stack__.currentsize] = 0; \ 13005a96b57dSJacob Faibussowitsch } \ 13015a96b57dSJacob Faibussowitsch stack__.hotdepth = PetscMax(stack__.hotdepth - 1, 0); \ 13025a96b57dSJacob Faibussowitsch } while (0) 13035a96b57dSJacob Faibussowitsch 1304586f9135SBarry Smith /*MC 1305586f9135SBarry Smith PetscStackPushNoCheck - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is 1306586f9135SBarry Smith currently in the source code. 1307586f9135SBarry Smith 1308586f9135SBarry Smith Synopsis: 1309586f9135SBarry Smith #include <petscsys.h> 1310586f9135SBarry Smith void PetscStackPushNoCheck(char *funct,int petsc_routine,PetscBool hot); 1311586f9135SBarry Smith 13127cdbe19fSJose E. Roman Not Collective 13137cdbe19fSJose E. Roman 1314586f9135SBarry Smith Input Parameters: 1315586f9135SBarry Smith + funct - the function name 1316586f9135SBarry Smith . petsc_routine - 2 user function, 1 PETSc function, 0 some other function 1317586f9135SBarry Smith - hot - indicates that the function may be called often so expensive error checking should be turned off inside the function 1318586f9135SBarry Smith 1319586f9135SBarry Smith Level: developer 1320586f9135SBarry Smith 1321586f9135SBarry Smith Notes: 1322586f9135SBarry 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 132387497f52SBarry 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 1324586f9135SBarry Smith help debug the problem. 1325586f9135SBarry Smith 1326ef1023bdSBarry Smith This version does not check the memory corruption (an expensive operation), use `PetscStackPush()` to check the memory. 1327ef1023bdSBarry Smith 1328792fecdfSBarry Smith Use `PetscStackPushExternal()` for a function call that is about to be made to a non-PETSc or user function (such as BLAS etc). 1329ef1023bdSBarry Smith 1330586f9135SBarry Smith The default stack is a global variable called `petscstack`. 1331586f9135SBarry Smith 1332586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1333ef1023bdSBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPush()`, `PetscStackPop`, 1334792fecdfSBarry Smith `PetscStackPushExternal()` 1335586f9135SBarry Smith M*/ 13369371c9d4SSatish Balay #define PetscStackPushNoCheck(funct, petsc_routine, hot) \ 13379371c9d4SSatish Balay do { \ 1338e04113cfSBarry Smith PetscStackSAWsTakeAccess(); \ 13395a96b57dSJacob Faibussowitsch PetscStackPush_Private(petscstack, __FILE__, funct, __LINE__, petsc_routine, hot); \ 1340e04113cfSBarry Smith PetscStackSAWsGrantAccess(); \ 1341441dd030SJed Brown } while (0) 1342441dd030SJed Brown 1343586f9135SBarry Smith /*MC 134487497f52SBarry Smith PetscStackUpdateLine - in a function that has a `PetscFunctionBegin` or `PetscFunctionBeginUser` updates the stack line number to the 134537154d10SBarry Smith current line number. 134637154d10SBarry Smith 134737154d10SBarry Smith Synopsis: 134837154d10SBarry Smith #include <petscsys.h> 134937154d10SBarry Smith void PetscStackUpdateLine 135037154d10SBarry Smith 13517cdbe19fSJose E. Roman Not Collective 13527cdbe19fSJose E. Roman 135337154d10SBarry Smith Level: developer 135437154d10SBarry Smith 135537154d10SBarry Smith Notes: 135687497f52SBarry Smith Using `PetscCall()` and friends automatically handles this process 135787497f52SBarry Smith 135837154d10SBarry 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 135937154d10SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 136037154d10SBarry Smith help debug the problem. 136137154d10SBarry Smith 136295bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 136337154d10SBarry Smith 136437154d10SBarry Smith This is used by `PetscCall()` and is otherwise not like to be needed 136537154d10SBarry Smith 136637154d10SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`, `PetscCall()` 136737154d10SBarry Smith M*/ 136837154d10SBarry Smith #define PetscStackUpdateLine \ 13693ba16761SJacob Faibussowitsch do { \ 1370f1e99121SPierre Jolivet if (petscstack.currentsize > 0 && petscstack.currentsize < PETSCSTACKSIZE && petscstack.function[petscstack.currentsize - 1] == PETSC_FUNCTION_NAME) { petscstack.line[petscstack.currentsize - 1] = __LINE__; } \ 13713ba16761SJacob Faibussowitsch } while (0) 137237154d10SBarry Smith 137337154d10SBarry Smith /*MC 1374792fecdfSBarry Smith PetscStackPushExternal - Pushes a new function name onto the PETSc default stack that tracks where the running program is 1375ef1023bdSBarry Smith currently in the source code. Does not include the filename or line number since this is called by the calling routine 1376ef1023bdSBarry Smith for non-PETSc or user functions. 1377ef1023bdSBarry Smith 1378ef1023bdSBarry Smith Synopsis: 1379ef1023bdSBarry Smith #include <petscsys.h> 1380792fecdfSBarry Smith void PetscStackPushExternal(char *funct); 1381ef1023bdSBarry Smith 13827cdbe19fSJose E. Roman Not Collective 13837cdbe19fSJose E. Roman 13842fe279fdSBarry Smith Input Parameter: 1385ef1023bdSBarry Smith . funct - the function name 1386ef1023bdSBarry Smith 1387ef1023bdSBarry Smith Level: developer 1388ef1023bdSBarry Smith 1389ef1023bdSBarry Smith Notes: 139087497f52SBarry Smith Using `PetscCallExternal()` and friends automatically handles this process 139187497f52SBarry Smith 1392ef1023bdSBarry 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 1393ef1023bdSBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1394ef1023bdSBarry Smith help debug the problem. 1395ef1023bdSBarry Smith 1396ef1023bdSBarry Smith The default stack is a global variable called `petscstack`. 1397ef1023bdSBarry Smith 1398ef1023bdSBarry Smith This is to be used when calling an external package function such as a BLAS function. 1399ef1023bdSBarry Smith 1400ef1023bdSBarry Smith This also updates the stack line number for the current stack function. 1401ef1023bdSBarry Smith 1402ef1023bdSBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1403ef1023bdSBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1404ef1023bdSBarry Smith M*/ 14059371c9d4SSatish Balay #define PetscStackPushExternal(funct) \ 14069371c9d4SSatish Balay do { \ 14079371c9d4SSatish Balay PetscStackUpdateLine; \ 14089371c9d4SSatish Balay PetscStackPushNoCheck(funct, 0, PETSC_TRUE); \ 1409a8f51744SPierre Jolivet } while (0) 1410ef1023bdSBarry Smith 1411ef1023bdSBarry Smith /*MC 1412586f9135SBarry Smith PetscStackPopNoCheck - Pops a function name from the PETSc default stack that tracks where the running program is 1413586f9135SBarry Smith currently in the source code. 1414586f9135SBarry Smith 1415586f9135SBarry Smith Synopsis: 1416586f9135SBarry Smith #include <petscsys.h> 1417586f9135SBarry Smith void PetscStackPopNoCheck(char *funct); 1418586f9135SBarry Smith 14197cdbe19fSJose E. Roman Not Collective 14207cdbe19fSJose E. Roman 1421586f9135SBarry Smith Input Parameter: 1422586f9135SBarry Smith . funct - the function name 1423586f9135SBarry Smith 1424586f9135SBarry Smith Level: developer 1425586f9135SBarry Smith 1426586f9135SBarry Smith Notes: 142787497f52SBarry Smith Using `PetscCall()`, `PetscCallExternal()`, `PetscCallBack()` and friends negates the need to call this 142887497f52SBarry Smith 1429586f9135SBarry 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 1430586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1431586f9135SBarry Smith help debug the problem. 1432586f9135SBarry Smith 143395bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 1434586f9135SBarry Smith 1435586f9135SBarry Smith Developer Note: 1436586f9135SBarry Smith `PetscStackPopNoCheck()` takes a function argument while `PetscStackPop` does not, this difference is likely just historical. 1437586f9135SBarry Smith 1438586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1439586f9135SBarry Smith M*/ 14409371c9d4SSatish Balay #define PetscStackPopNoCheck(funct) \ 14419371c9d4SSatish Balay do { \ 1442362febeeSStefano Zampini PetscStackSAWsTakeAccess(); \ 14435a96b57dSJacob Faibussowitsch PetscStackPop_Private(petscstack, funct); \ 1444362febeeSStefano Zampini PetscStackSAWsGrantAccess(); \ 1445362febeeSStefano Zampini } while (0) 1446362febeeSStefano Zampini 14479371c9d4SSatish Balay #define PetscStackClearTop \ 14489371c9d4SSatish Balay do { \ 1449e04113cfSBarry Smith PetscStackSAWsTakeAccess(); \ 14509371c9d4SSatish Balay if (petscstack.currentsize > 0 && --petscstack.currentsize < PETSCSTACKSIZE) { \ 145127104ee2SJacob Faibussowitsch petscstack.function[petscstack.currentsize] = PETSC_NULLPTR; \ 145227104ee2SJacob Faibussowitsch petscstack.file[petscstack.currentsize] = PETSC_NULLPTR; \ 145327104ee2SJacob Faibussowitsch petscstack.line[petscstack.currentsize] = 0; \ 145427104ee2SJacob Faibussowitsch petscstack.petscroutine[petscstack.currentsize] = 0; \ 1455441dd030SJed Brown } \ 145627104ee2SJacob Faibussowitsch petscstack.hotdepth = PetscMax(petscstack.hotdepth - 1, 0); \ 1457e04113cfSBarry Smith PetscStackSAWsGrantAccess(); \ 1458441dd030SJed Brown } while (0) 1459441dd030SJed Brown 146030de9b25SBarry Smith /*MC 14611957e957SBarry Smith PetscFunctionBegin - First executable line of each PETSc function, used for error handling. Final 146287497f52SBarry Smith line of PETSc functions should be `PetscFunctionReturn`(0); 146330de9b25SBarry Smith 146430de9b25SBarry Smith Synopsis: 1465aaa7dc30SBarry Smith #include <petscsys.h> 146630de9b25SBarry Smith void PetscFunctionBegin; 146730de9b25SBarry Smith 146820f4b53cSBarry Smith Not Collective; No Fortran Support 1469eca87e8dSBarry Smith 147030de9b25SBarry Smith Usage: 147130de9b25SBarry Smith .vb 147230de9b25SBarry Smith int something; 147330de9b25SBarry Smith 147430de9b25SBarry Smith PetscFunctionBegin; 147530de9b25SBarry Smith .ve 147630de9b25SBarry Smith 147730de9b25SBarry Smith Level: developer 147830de9b25SBarry Smith 147920f4b53cSBarry Smith Note: 148020f4b53cSBarry Smith Use `PetscFunctionBeginUser` for application codes. 148120f4b53cSBarry Smith 1482586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()` 148330de9b25SBarry Smith 148430de9b25SBarry Smith M*/ 14859371c9d4SSatish Balay #define PetscFunctionBegin \ 14869371c9d4SSatish Balay do { \ 1487362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_FALSE); \ 1488a2f94806SJed Brown PetscRegister__FUNCT__(); \ 1489a2f94806SJed Brown } while (0) 1490a2f94806SJed Brown 1491a2f94806SJed Brown /*MC 149287497f52SBarry Smith PetscFunctionBeginHot - Substitute for `PetscFunctionBegin` to be used in functions that are called in 1493a2f94806SJed Brown performance-critical circumstances. Use of this function allows for lighter profiling by default. 1494a2f94806SJed Brown 1495a2f94806SJed Brown Synopsis: 1496aaa7dc30SBarry Smith #include <petscsys.h> 1497a2f94806SJed Brown void PetscFunctionBeginHot; 1498a2f94806SJed Brown 149920f4b53cSBarry Smith Not Collective; No Fortran Support 1500a2f94806SJed Brown 1501a2f94806SJed Brown Usage: 1502a2f94806SJed Brown .vb 1503a2f94806SJed Brown int something; 1504a2f94806SJed Brown 1505a2f94806SJed Brown PetscFunctionBeginHot; 1506a2f94806SJed Brown .ve 1507a2f94806SJed Brown 1508a2f94806SJed Brown Level: developer 1509a2f94806SJed Brown 1510586f9135SBarry Smith .seealso: `PetscFunctionBegin`, `PetscFunctionReturn()`, `PetscStackPushNoCheck()` 1511a2f94806SJed Brown 1512a2f94806SJed Brown M*/ 15139371c9d4SSatish Balay #define PetscFunctionBeginHot \ 15149371c9d4SSatish Balay do { \ 1515362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_TRUE); \ 15162d53ad75SBarry Smith PetscRegister__FUNCT__(); \ 151753c77d0aSJed Brown } while (0) 151853c77d0aSJed Brown 1519a8d2bbe5SBarry Smith /*MC 1520530d85c1SBarry Smith PetscFunctionBeginUser - First executable line of user provided routines 1521a8d2bbe5SBarry Smith 1522a8d2bbe5SBarry Smith Synopsis: 1523aaa7dc30SBarry Smith #include <petscsys.h> 1524a8d2bbe5SBarry Smith void PetscFunctionBeginUser; 1525a8d2bbe5SBarry Smith 152620f4b53cSBarry Smith Not Collective; No Fortran Support 1527a8d2bbe5SBarry Smith 1528a8d2bbe5SBarry Smith Usage: 1529a8d2bbe5SBarry Smith .vb 1530a8d2bbe5SBarry Smith int something; 1531a8d2bbe5SBarry Smith 1532ac285190SBarry Smith PetscFunctionBeginUser; 1533a8d2bbe5SBarry Smith .ve 1534a8d2bbe5SBarry Smith 153520f4b53cSBarry Smith Level: intermediate 153620f4b53cSBarry Smith 1537a8d2bbe5SBarry Smith Notes: 1538530d85c1SBarry Smith Functions that incorporate this must call `PetscFunctionReturn()` instead of return except for main(). 1539530d85c1SBarry Smith 1540530d85c1SBarry Smith May be used before `PetscInitialize()` 15411957e957SBarry Smith 1542530d85c1SBarry Smith This is identical to `PetscFunctionBegin` except it labels the routine as a user 1543ac285190SBarry Smith routine instead of as a PETSc library routine. 1544ac285190SBarry Smith 1545586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, `PetscFunctionBeginHot`, `PetscStackPushNoCheck()` 1546a8d2bbe5SBarry Smith M*/ 15479371c9d4SSatish Balay #define PetscFunctionBeginUser \ 15489371c9d4SSatish Balay do { \ 1549362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 2, PETSC_FALSE); \ 1550a8d2bbe5SBarry Smith PetscRegister__FUNCT__(); \ 1551a8d2bbe5SBarry Smith } while (0) 1552a8d2bbe5SBarry Smith 1553586f9135SBarry Smith /*MC 1554586f9135SBarry Smith PetscStackPush - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is 1555586f9135SBarry Smith currently in the source code and verifies the memory is not corrupted. 1556586f9135SBarry Smith 1557586f9135SBarry Smith Synopsis: 1558586f9135SBarry Smith #include <petscsys.h> 1559586f9135SBarry Smith void PetscStackPush(char *funct) 1560586f9135SBarry Smith 15617cdbe19fSJose E. Roman Not Collective 15627cdbe19fSJose E. Roman 1563586f9135SBarry Smith Input Parameter: 1564586f9135SBarry Smith . funct - the function name 1565586f9135SBarry Smith 1566586f9135SBarry Smith Level: developer 1567586f9135SBarry Smith 1568586f9135SBarry Smith Notes: 1569586f9135SBarry 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 1570586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1571586f9135SBarry Smith help debug the problem. 1572586f9135SBarry Smith 157395bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 1574586f9135SBarry Smith 1575586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1576586f9135SBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1577586f9135SBarry Smith M*/ 15789371c9d4SSatish Balay #define PetscStackPush(n) \ 15799371c9d4SSatish Balay do { \ 1580362febeeSStefano Zampini PetscStackPushNoCheck(n, 0, PETSC_FALSE); \ 158115681b3cSBarry Smith CHKMEMQ; \ 158215681b3cSBarry Smith } while (0) 15833a40ed3dSBarry Smith 1584586f9135SBarry Smith /*MC 1585586f9135SBarry Smith PetscStackPop - Pops a function name from the PETSc default stack that tracks where the running program is 1586586f9135SBarry Smith currently in the source code and verifies the memory is not corrupted. 1587586f9135SBarry Smith 1588586f9135SBarry Smith Synopsis: 1589586f9135SBarry Smith #include <petscsys.h> 1590586f9135SBarry Smith void PetscStackPop 1591586f9135SBarry Smith 15927cdbe19fSJose E. Roman Not Collective 15937cdbe19fSJose E. Roman 1594586f9135SBarry Smith Level: developer 1595586f9135SBarry Smith 1596586f9135SBarry Smith Notes: 1597586f9135SBarry 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 1598586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1599586f9135SBarry Smith help debug the problem. 1600586f9135SBarry Smith 160195bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 1602586f9135SBarry Smith 1603586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPopNoCheck()`, `PetscStackPush()` 1604586f9135SBarry Smith M*/ 16059371c9d4SSatish Balay #define PetscStackPop \ 16069371c9d4SSatish Balay do { \ 1607441dd030SJed Brown CHKMEMQ; \ 1608362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 160915681b3cSBarry Smith } while (0) 1610d64ed03dSBarry Smith 161130de9b25SBarry Smith /*MC 16120a57284eSJacob Faibussowitsch PetscFunctionReturn - Last executable line of each PETSc function used for error 16130a57284eSJacob Faibussowitsch handling. Replaces `return()`. 161430de9b25SBarry Smith 161530de9b25SBarry Smith Synopsis: 16160a57284eSJacob Faibussowitsch #include <petscerror.h> 16170a57284eSJacob Faibussowitsch void PetscFunctionReturn(...) 161830de9b25SBarry Smith 1619667f096bSBarry Smith Not Collective; No Fortran Support 1620eca87e8dSBarry Smith 16215687f061SJacob Faibussowitsch Level: beginner 16225687f061SJacob Faibussowitsch 16230a57284eSJacob Faibussowitsch Notes: 16240a57284eSJacob Faibussowitsch This routine is a macro, so while it does not "return" anything itself, it does return from 16250a57284eSJacob Faibussowitsch the function in the literal sense. 16260a57284eSJacob Faibussowitsch 16270a57284eSJacob Faibussowitsch Usually the return value is the integer literal `0` (for example in any function returning 16280a57284eSJacob Faibussowitsch `PetscErrorCode`), however it is possible to return any arbitrary type. The arguments of 16290a57284eSJacob Faibussowitsch this macro are placed before the `return` statement as-is. 16300a57284eSJacob Faibussowitsch 16310a57284eSJacob Faibussowitsch Any routine which returns via `PetscFunctionReturn()` must begin with a corresponding 16320a57284eSJacob Faibussowitsch `PetscFunctionBegin`. 16330a57284eSJacob Faibussowitsch 16340a57284eSJacob Faibussowitsch For routines which return `void` use `PetscFunctionReturnVoid()` instead. 16350a57284eSJacob Faibussowitsch 16360a57284eSJacob Faibussowitsch Example Usage: 163730de9b25SBarry Smith .vb 16380a57284eSJacob Faibussowitsch PetscErrorCode foo(int *x) 16390a57284eSJacob Faibussowitsch { 16400a57284eSJacob Faibussowitsch PetscFunctionBegin; // don't forget the begin! 16410a57284eSJacob Faibussowitsch *x = 10; 16423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 164330de9b25SBarry Smith } 164430de9b25SBarry Smith .ve 164530de9b25SBarry Smith 16460a57284eSJacob Faibussowitsch May return any arbitrary type\: 16470a57284eSJacob Faibussowitsch .vb 16480a57284eSJacob Faibussowitsch struct Foo 16490a57284eSJacob Faibussowitsch { 16500a57284eSJacob Faibussowitsch int x; 16510a57284eSJacob Faibussowitsch }; 16520a57284eSJacob Faibussowitsch 16530a57284eSJacob Faibussowitsch struct Foo make_foo(int value) 16540a57284eSJacob Faibussowitsch { 16550a57284eSJacob Faibussowitsch struct Foo f; 16560a57284eSJacob Faibussowitsch 16570a57284eSJacob Faibussowitsch PetscFunctionBegin; 16580a57284eSJacob Faibussowitsch f.x = value; 16590a57284eSJacob Faibussowitsch PetscFunctionReturn(f); 16600a57284eSJacob Faibussowitsch } 16610a57284eSJacob Faibussowitsch .ve 16620a57284eSJacob Faibussowitsch 16630a57284eSJacob Faibussowitsch .seealso: `PetscFunctionBegin`, `PetscFunctionBeginUser`, `PetscFunctionReturnVoid()`, 16640a57284eSJacob Faibussowitsch `PetscStackPopNoCheck()` 166530de9b25SBarry Smith M*/ 16665687f061SJacob Faibussowitsch #define PetscFunctionReturn(...) \ 16679371c9d4SSatish Balay do { \ 1668362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 16695687f061SJacob Faibussowitsch return __VA_ARGS__; \ 167027104ee2SJacob Faibussowitsch } while (0) 1671d64ed03dSBarry Smith 16720a57284eSJacob Faibussowitsch /*MC 16730a57284eSJacob Faibussowitsch PetscFunctionReturnVoid - Like `PetscFunctionReturn()` but returns `void` 16740a57284eSJacob Faibussowitsch 16750a57284eSJacob Faibussowitsch Synopsis: 16760a57284eSJacob Faibussowitsch #include <petscerror.h> 16770a57284eSJacob Faibussowitsch void PetscFunctionReturnVoid() 16780a57284eSJacob Faibussowitsch 16790a57284eSJacob Faibussowitsch Not Collective 16800a57284eSJacob Faibussowitsch 16815687f061SJacob Faibussowitsch Level: beginner 16825687f061SJacob Faibussowitsch 16835687f061SJacob Faibussowitsch Note: 16840a57284eSJacob Faibussowitsch Behaves identically to `PetscFunctionReturn()` except that it returns `void`. That is, this 16855687f061SJacob Faibussowitsch macro culminates with `return`. 16860a57284eSJacob Faibussowitsch 16870a57284eSJacob Faibussowitsch Example Usage: 16880a57284eSJacob Faibussowitsch .vb 16890a57284eSJacob Faibussowitsch void foo() 16900a57284eSJacob Faibussowitsch { 16910a57284eSJacob Faibussowitsch PetscFunctionBegin; // must start with PetscFunctionBegin! 16920a57284eSJacob Faibussowitsch bar(); 16930a57284eSJacob Faibussowitsch baz(); 16940a57284eSJacob Faibussowitsch PetscFunctionReturnVoid(); 16950a57284eSJacob Faibussowitsch } 16960a57284eSJacob Faibussowitsch .ve 16970a57284eSJacob Faibussowitsch 16980a57284eSJacob Faibussowitsch .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, PetscFunctionBeginUser` 16990a57284eSJacob Faibussowitsch M*/ 17009371c9d4SSatish Balay #define PetscFunctionReturnVoid() \ 17019371c9d4SSatish Balay do { \ 1702362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 170327104ee2SJacob Faibussowitsch return; \ 170427104ee2SJacob Faibussowitsch } while (0) 170527104ee2SJacob Faibussowitsch #else /* PETSC_USE_DEBUG */ 170627104ee2SJacob Faibussowitsch #define PetscStackPushNoCheck(funct, petsc_routine, hot) 170737154d10SBarry Smith #define PetscStackUpdateLine 1708792fecdfSBarry Smith #define PetscStackPushExternal(funct) 17093ba16761SJacob Faibussowitsch #define PetscStackPopNoCheck(...) 171027104ee2SJacob Faibussowitsch #define PetscStackClearTop 17113a40ed3dSBarry Smith #define PetscFunctionBegin 17120bdf7c52SPeter Brune #define PetscFunctionBeginUser 1713a2f94806SJed Brown #define PetscFunctionBeginHot 17140a57284eSJacob Faibussowitsch #define PetscFunctionReturn(...) return __VA_ARGS__ 17155665465eSBarry Smith #define PetscFunctionReturnVoid() return 1716812af9f3SBarry Smith #define PetscStackPop CHKMEMQ 1717812af9f3SBarry Smith #define PetscStackPush(f) CHKMEMQ 171827104ee2SJacob Faibussowitsch #endif /* PETSC_USE_DEBUG */ 17193a40ed3dSBarry Smith 17206d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 17213ba16761SJacob Faibussowitsch #define PetscStackCallExternalVoid(...) 17223ba16761SJacob Faibussowitsch template <typename F, typename... Args> 17233ba16761SJacob Faibussowitsch void PetscCallExternal(F, Args...); 17246d210af2SJacob Faibussowitsch #else 1725586f9135SBarry Smith /*MC 1726e77caa6dSBarry Smith PetscStackCallExternalVoid - Calls an external library routine or user function after pushing the name of the routine on the stack. 1727eb6b5d47SBarry Smith 1728eb6b5d47SBarry Smith Input Parameters: 1729eb6b5d47SBarry Smith + name - string that gives the name of the function being called 1730586f9135SBarry Smith - routine - actual call to the routine, for example, functionname(a,b) 1731fd3f9acdSBarry Smith 1732586f9135SBarry Smith Level: developer 1733eb6b5d47SBarry Smith 173495bd0b28SBarry Smith Notes: 1735792fecdfSBarry Smith Often one should use `PetscCallExternal()` instead. This routine is intended for external library routines that DO NOT return error codes 1736eb6b5d47SBarry Smith 1737586f9135SBarry Smith In debug mode this also checks the memory for corruption at the end of the function call. 1738586f9135SBarry Smith 173995bd0b28SBarry Smith Certain external packages, such as BLAS/LAPACK may have their own macros, `PetscCallBLAS()` for managing the call, error checking, etc. 1740586f9135SBarry Smith 1741586f9135SBarry Smith Developer Note: 1742586f9135SBarry 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. 1743586f9135SBarry Smith 1744792fecdfSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscCallExternal()`, `PetscCallBLAS()` 1745586f9135SBarry Smith @*/ 17463ba16761SJacob Faibussowitsch #define PetscStackCallExternalVoid(name, ...) \ 17479371c9d4SSatish Balay do { \ 174828770333SStefano Zampini PetscStackPushExternal(name); \ 17493ba16761SJacob Faibussowitsch __VA_ARGS__; \ 17509371c9d4SSatish Balay PetscStackPop; \ 17519371c9d4SSatish Balay } while (0) 1752eb6b5d47SBarry Smith 1753586f9135SBarry Smith /*MC 1754792fecdfSBarry Smith PetscCallExternal - Calls an external library routine that returns an error code after pushing the name of the routine on the stack. 1755fd3f9acdSBarry Smith 1756fd3f9acdSBarry Smith Input Parameters: 1757fd3f9acdSBarry Smith + func - name of the routine 1758586f9135SBarry Smith - args - arguments to the routine 1759586f9135SBarry Smith 1760586f9135SBarry Smith Level: developer 1761fd3f9acdSBarry Smith 176295452b02SPatrick Sanan Notes: 1763e77caa6dSBarry Smith This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not. 1764dbf62e16SBarry Smith 1765586f9135SBarry Smith In debug mode this also checks the memory for corruption at the end of the function call. 1766fd3f9acdSBarry Smith 176787497f52SBarry Smith Assumes the error return code of the function is an integer and that a value of 0 indicates success 176887497f52SBarry Smith 1769586f9135SBarry Smith Developer Note: 1770d5b43468SJose E. Roman This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc. 1771586f9135SBarry Smith 1772e77caa6dSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()` 1773586f9135SBarry Smith M*/ 17749371c9d4SSatish Balay #define PetscCallExternal(func, ...) \ 17759371c9d4SSatish Balay do { \ 1776a74df02fSJacob Faibussowitsch PetscStackPush(PetscStringize(func)); \ 17773ba16761SJacob Faibussowitsch int ierr_petsc_call_external_ = func(__VA_ARGS__); \ 17781d4906efSStefano Zampini PetscStackPop; \ 17793ba16761SJacob 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_); \ 1780fd3f9acdSBarry Smith } while (0) 17816d210af2SJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */ 1782