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()`, 52*af27ebaaSBarry 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 90*af27ebaaSBarry 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 116*af27ebaaSBarry 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 146*af27ebaaSBarry 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 1552c71b3e2SJacob Faibussowitsch PetscCheck - Check that a particular condition is true 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 182*af27ebaaSBarry 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 212*af27ebaaSBarry 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 232*af27ebaaSBarry 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 243*af27ebaaSBarry 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 30087497f52SBarry 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: 33049c86fc7SBarry Smith The Fortran function from which this is used must declare a variable PetscErrorCode ierr and ierr must be 33187497f52SBarry Smith the final argument to the PETSc function being called. 33249c86fc7SBarry Smith 33349c86fc7SBarry Smith 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 352989c49a2SBarry Smith PetscCallA - Fortran-only macro that should be used in the main program to call PETSc functions instead of using 353989c49a2SBarry Smith 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 369989c49a2SBarry Smith Use `SETERRA()` to set an error in a Fortran main program and `SETERRQ()` in Fortran subroutines 370989c49a2SBarry Smith 371989c49a2SBarry Smith .seealso: `SETERRQ()`, `SETERRA()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()` 372989c49a2SBarry Smith M*/ 373989c49a2SBarry Smith 374989c49a2SBarry Smith /*MC 375792fecdfSBarry 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 376ef1023bdSBarry Smith handler and returns from the current function with the error code. 377ef1023bdSBarry Smith 378ef1023bdSBarry Smith Synopsis: 379ef1023bdSBarry Smith #include <petscerror.h> 380792fecdfSBarry Smith void PetscCallBack(const char *functionname,PetscFunction(args)) 381ef1023bdSBarry Smith 3827cdbe19fSJose E. Roman Not Collective; No Fortran Support 383ef1023bdSBarry Smith 384ef1023bdSBarry Smith Input Parameters: 385ef1023bdSBarry Smith + functionname - the name of the function being called, this can be a string with spaces that describes the meaning of the callback 38687497f52SBarry Smith - PetscFunction - user provided callback function that returns an error code 387ef1023bdSBarry Smith 388ef1023bdSBarry Smith Example Usage: 389ef1023bdSBarry Smith .vb 390792fecdfSBarry Smith PetscCallBack("XXX callback to do something",a->callback(...)); 391ef1023bdSBarry Smith .ve 392ef1023bdSBarry Smith 393ef1023bdSBarry Smith Level: developer 394ef1023bdSBarry Smith 395667f096bSBarry Smith Notes: 396667f096bSBarry Smith Once the error handler is called the calling function is then returned from with the given 397667f096bSBarry Smith error code. Experienced users can set the error handler with `PetscPushErrorHandler()`. 398667f096bSBarry Smith 399667f096bSBarry Smith `PetscCallBack()` should only be called in PETSc when a call is being made to a user provided call-back routine. 400667f096bSBarry Smith 40187497f52SBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscCall()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()` 402ef1023bdSBarry Smith `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, `CHKERRMPI()`, `PetscCall()` 403ef1023bdSBarry Smith M*/ 404ef1023bdSBarry Smith 4053ba16761SJacob Faibussowitsch /*MC 4063ba16761SJacob Faibussowitsch PetscCallVoid - Like `PetscCall()` but for functions returning `void` 4073ba16761SJacob Faibussowitsch 4083ba16761SJacob Faibussowitsch Synopsis: 4093ba16761SJacob Faibussowitsch #include <petscerror.h> 4103ba16761SJacob Faibussowitsch void PetscCall(PetscFunction(args)) 4113ba16761SJacob Faibussowitsch 4127cdbe19fSJose E. Roman Not Collective; No Fortran Support 4133ba16761SJacob Faibussowitsch 4143ba16761SJacob Faibussowitsch Input Parameter: 4153ba16761SJacob Faibussowitsch . PetscFunction - any PETSc function that returns an error code 4163ba16761SJacob Faibussowitsch 4173ba16761SJacob Faibussowitsch Example Usage: 4183ba16761SJacob Faibussowitsch .vb 4193ba16761SJacob Faibussowitsch void foo() 4203ba16761SJacob Faibussowitsch { 4213ba16761SJacob Faibussowitsch KSP ksp; 4223ba16761SJacob Faibussowitsch 4233ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 4243ba16761SJacob Faibussowitsch // OK, properly handles PETSc error codes 4253ba16761SJacob Faibussowitsch PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp)); 4263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4273ba16761SJacob Faibussowitsch } 4283ba16761SJacob Faibussowitsch 4293ba16761SJacob Faibussowitsch PetscErrorCode bar() 4303ba16761SJacob Faibussowitsch { 4313ba16761SJacob Faibussowitsch KSP ksp; 4323ba16761SJacob Faibussowitsch 4333ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 4343ba16761SJacob Faibussowitsch // ERROR, Non-void function 'bar' should return a value 4353ba16761SJacob Faibussowitsch PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp)); 4363ba16761SJacob Faibussowitsch // OK, returning PetscErrorCode 4373ba16761SJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp)); 4383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4393ba16761SJacob Faibussowitsch } 440667f096bSBarry Smith .ve 4413ba16761SJacob Faibussowitsch 4423ba16761SJacob Faibussowitsch Level: beginner 4433ba16761SJacob Faibussowitsch 444667f096bSBarry Smith Notes: 445667f096bSBarry Smith Has identical usage to `PetscCall()`, except that it returns `void` on error instead of a 446667f096bSBarry Smith `PetscErrorCode`. See `PetscCall()` for more detailed discussion. 447667f096bSBarry Smith 448667f096bSBarry Smith Note that users should prefer `PetscCallAbort()` to this routine. While this routine does 449667f096bSBarry Smith "handle" errors by returning from the enclosing function, it effectively gobbles the 450667f096bSBarry Smith error. Since the enclosing function itself returns `void`, its callers have no way of knowing 451667f096bSBarry Smith that the routine returned early due to an error. `PetscCallAbort()` at least ensures that the 452667f096bSBarry Smith program crashes gracefully. 453667f096bSBarry Smith 4543ba16761SJacob Faibussowitsch .seealso: `PetscCall()`, `PetscErrorCode` 4553ba16761SJacob Faibussowitsch M*/ 4569566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 4579566063dSJacob Faibussowitsch void PetscCall(PetscErrorCode); 458792fecdfSBarry Smith void PetscCallBack(const char *, PetscErrorCode); 4599566063dSJacob Faibussowitsch void PetscCallVoid(PetscErrorCode); 4609566063dSJacob Faibussowitsch #else 4619371c9d4SSatish Balay #define PetscCall(...) \ 4629371c9d4SSatish Balay do { \ 4633ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_q_; \ 46437154d10SBarry Smith PetscStackUpdateLine; \ 4653ba16761SJacob Faibussowitsch ierr_petsc_call_q_ = __VA_ARGS__; \ 4663ba16761SJacob 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, " "); \ 4679566063dSJacob Faibussowitsch } while (0) 4689371c9d4SSatish Balay #define PetscCallBack(function, ...) \ 4699371c9d4SSatish Balay do { \ 4703ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_q_; \ 471e33ced7fSLisandro Dalcin PetscStackUpdateLine; \ 472e33ced7fSLisandro Dalcin PetscStackPushExternal(function); \ 4733ba16761SJacob Faibussowitsch ierr_petsc_call_q_ = __VA_ARGS__; \ 474ef1023bdSBarry Smith PetscStackPop; \ 4753ba16761SJacob 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, " "); \ 476ef1023bdSBarry Smith } while (0) 4779371c9d4SSatish Balay #define PetscCallVoid(...) \ 4789371c9d4SSatish Balay do { \ 4793ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_void_; \ 480e33ced7fSLisandro Dalcin PetscStackUpdateLine; \ 4813ba16761SJacob Faibussowitsch ierr_petsc_call_void_ = __VA_ARGS__; \ 4823ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_void_ != PETSC_SUCCESS)) { \ 4833ba16761SJacob Faibussowitsch ierr_petsc_call_void_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_void_, PETSC_ERROR_REPEAT, " "); \ 4843ba16761SJacob Faibussowitsch (void)ierr_petsc_call_void_; \ 4859371c9d4SSatish Balay return; \ 4869371c9d4SSatish Balay } \ 4879566063dSJacob Faibussowitsch } while (0) 4889566063dSJacob Faibussowitsch #endif 4899566063dSJacob Faibussowitsch 4909566063dSJacob Faibussowitsch /*MC 4919566063dSJacob Faibussowitsch CHKERRQ - Checks error code returned from PETSc function 49230de9b25SBarry Smith 49330de9b25SBarry Smith Synopsis: 494aaa7dc30SBarry Smith #include <petscsys.h> 4959566063dSJacob Faibussowitsch void CHKERRQ(PetscErrorCode ierr) 4969566063dSJacob Faibussowitsch 4979566063dSJacob Faibussowitsch Not Collective 4989566063dSJacob Faibussowitsch 4992fe279fdSBarry Smith Input Parameter: 5009566063dSJacob Faibussowitsch . ierr - nonzero error code 5019566063dSJacob Faibussowitsch 5029566063dSJacob Faibussowitsch Level: deprecated 5039566063dSJacob Faibussowitsch 504667f096bSBarry Smith Note: 505667f096bSBarry Smith Deprecated in favor of `PetscCall()`. This routine behaves identically to it. 506667f096bSBarry Smith 507db781477SPatrick Sanan .seealso: `PetscCall()` 5089566063dSJacob Faibussowitsch M*/ 5099566063dSJacob Faibussowitsch #define CHKERRQ(...) PetscCall(__VA_ARGS__) 5109566063dSJacob Faibussowitsch #define CHKERRV(...) PetscCallVoid(__VA_ARGS__) 5119566063dSJacob Faibussowitsch 512db9cea48SBarry Smith PETSC_EXTERN void PetscMPIErrorString(PetscMPIInt, char *); 513db9cea48SBarry Smith 5149566063dSJacob Faibussowitsch /*MC 5159566063dSJacob Faibussowitsch PetscCallMPI - Checks error code returned from MPI calls, if non-zero it calls the error 5169566063dSJacob Faibussowitsch handler and then returns 5179566063dSJacob Faibussowitsch 5189566063dSJacob Faibussowitsch Synopsis: 5199566063dSJacob Faibussowitsch #include <petscerror.h> 52049c86fc7SBarry Smith void PetscCallMPI(MPI_Function(args)) 52130de9b25SBarry Smith 522eca87e8dSBarry Smith Not Collective 52330de9b25SBarry Smith 5242fe279fdSBarry Smith Input Parameter: 52549c86fc7SBarry Smith . MPI_Function - an MPI function that returns an MPI error code 52630de9b25SBarry Smith 527667f096bSBarry Smith Level: beginner 528667f096bSBarry Smith 5299566063dSJacob Faibussowitsch Notes: 53087497f52SBarry Smith Always returns the error code `PETSC_ERR_MPI`; the MPI error code and string are embedded in 5319566063dSJacob Faibussowitsch the string error message. Do not use this to call any other routines (for example PETSc 5323ba16761SJacob Faibussowitsch routines), it should only be used for direct MPI calls. The user may configure PETSc with the 5333ba16761SJacob Faibussowitsch `--with-strict-petscerrorcode` option to check this at compile-time, otherwise they must 5349566063dSJacob Faibussowitsch check this themselves. 5359566063dSJacob Faibussowitsch 536aaa8cc7dSPierre Jolivet This routine can only be used in functions returning `PetscErrorCode` themselves. If the 5373ba16761SJacob Faibussowitsch calling function returns a different type, use `PetscCallMPIAbort()` instead. 5383ba16761SJacob Faibussowitsch 5399566063dSJacob Faibussowitsch Example Usage: 5409566063dSJacob Faibussowitsch .vb 5419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(...)); // OK, calling MPI function 5429566063dSJacob Faibussowitsch 5439566063dSJacob Faibussowitsch PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead! 5449566063dSJacob Faibussowitsch .ve 5459566063dSJacob Faibussowitsch 54649c86fc7SBarry Smith Fortran Notes: 54787497f52SBarry Smith The Fortran function from which this is used must declare a variable `PetscErrorCode` ierr and ierr must be 54849c86fc7SBarry Smith the final argument to the MPI function being called. 54949c86fc7SBarry Smith 55049c86fc7SBarry Smith In the main program and in Fortran subroutines that do not have ierr as the final return parameter one 55187497f52SBarry Smith should use `PetscCallMPIA()` 55249c86fc7SBarry Smith 55349c86fc7SBarry Smith Fortran Usage: 55449c86fc7SBarry Smith .vb 55549c86fc7SBarry Smith PetscErrorCode ierr or integer ierr 55649c86fc7SBarry Smith ... 55749c86fc7SBarry Smith PetscCallMPI(MPI_Comm_size(...,ierr)) 55849c86fc7SBarry Smith PetscCallMPIA(MPI_Comm_size(...,ierr)) ! Will abort after calling error handler 55949c86fc7SBarry Smith 56049c86fc7SBarry Smith PetscCallMPI(MPI_Comm_size(...,eflag)) ! ERROR, final argument must be ierr 56149c86fc7SBarry Smith .ve 56249c86fc7SBarry Smith 563db781477SPatrick Sanan .seealso: `SETERRMPI()`, `PetscCall()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`, 5643ba16761SJacob Faibussowitsch `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 5653ba16761SJacob Faibussowitsch `PetscError()`, `CHKMEMQ` 5663ba16761SJacob Faibussowitsch M*/ 5673ba16761SJacob Faibussowitsch 5683ba16761SJacob Faibussowitsch /*MC 5693ba16761SJacob Faibussowitsch PetscCallMPIAbort - Like `PetscCallMPI()` but calls `MPI_Abort()` on error 5703ba16761SJacob Faibussowitsch 5713ba16761SJacob Faibussowitsch Synopsis: 5723ba16761SJacob Faibussowitsch #include <petscerror.h> 5733ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm comm, MPI_Function(args)) 5743ba16761SJacob Faibussowitsch 5753ba16761SJacob Faibussowitsch Not Collective 5763ba16761SJacob Faibussowitsch 5773ba16761SJacob Faibussowitsch Input Parameters: 5783ba16761SJacob Faibussowitsch + comm - the MPI communicator to abort on 5793ba16761SJacob Faibussowitsch - MPI_Function - an MPI function that returns an MPI error code 5803ba16761SJacob Faibussowitsch 581667f096bSBarry Smith Level: beginner 582667f096bSBarry Smith 5833ba16761SJacob Faibussowitsch Notes: 5843ba16761SJacob Faibussowitsch Usage is identical to `PetscCallMPI()`. See `PetscCallMPI()` for detailed discussion. 5853ba16761SJacob Faibussowitsch 5863ba16761SJacob Faibussowitsch This routine may be used in functions returning `void` or other non-`PetscErrorCode` types. 5873ba16761SJacob Faibussowitsch 588989c49a2SBarry Smith Fortran Note: 589989c49a2SBarry Smith In Fortran this is called `PetscCallMPIA()` and is intended to be used in the main program while `PetscCallMPI()` is 590989c49a2SBarry Smith used in Fortran subroutines. 591989c49a2SBarry Smith 592989c49a2SBarry Smith Developer Note: 593989c49a2SBarry Smith This should have the same name in Fortran. 594989c49a2SBarry Smith 5953ba16761SJacob Faibussowitsch .seealso: `PetscCallMPI()`, `PetscCallAbort()`, `SETERRABORT()` 59630de9b25SBarry Smith M*/ 5973fcd9f07SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 59863a3b9bcSJacob Faibussowitsch void PetscCallMPI(PetscMPIInt); 5993ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm, PetscMPIInt); 600064a246eSJacob Faibussowitsch #else 6013ba16761SJacob Faibussowitsch #define PetscCallMPI_Private(__PETSC_STACK_POP_FUNC__, __SETERR_FUNC__, __COMM__, ...) \ 6029371c9d4SSatish Balay do { \ 6033ba16761SJacob Faibussowitsch PetscMPIInt ierr_petsc_call_mpi_; \ 604ef1023bdSBarry Smith PetscStackUpdateLine; \ 605792fecdfSBarry Smith PetscStackPushExternal("MPI function"); \ 606d71ae5a4SJacob Faibussowitsch { \ 6073ba16761SJacob Faibussowitsch ierr_petsc_call_mpi_ = __VA_ARGS__; \ 608d71ae5a4SJacob Faibussowitsch } \ 6093ba16761SJacob Faibussowitsch __PETSC_STACK_POP_FUNC__; \ 6103ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_mpi_ != MPI_SUCCESS)) { \ 6113ba16761SJacob Faibussowitsch char petsc_mpi_7_errorstring[2 * MPI_MAX_ERROR_STRING]; \ 6123ba16761SJacob Faibussowitsch PetscMPIErrorString(ierr_petsc_call_mpi_, (char *)petsc_mpi_7_errorstring); \ 6133ba16761SJacob Faibussowitsch __SETERR_FUNC__(__COMM__, PETSC_ERR_MPI, "MPI error %d %s", (int)ierr_petsc_call_mpi_, petsc_mpi_7_errorstring); \ 6143fcd9f07SJacob Faibussowitsch } \ 6153fcd9f07SJacob Faibussowitsch } while (0) 6163ba16761SJacob Faibussowitsch 6173ba16761SJacob Faibussowitsch #define PetscCallMPI(...) PetscCallMPI_Private(PetscStackPop, SETERRQ, PETSC_COMM_SELF, __VA_ARGS__) 6183ba16761SJacob Faibussowitsch #define PetscCallMPIAbort(comm, ...) PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRABORT, comm, __VA_ARGS__) 6199566063dSJacob Faibussowitsch #endif 620064a246eSJacob Faibussowitsch 6217037b0edSPatrick Sanan /*MC 6229566063dSJacob Faibussowitsch CHKERRMPI - Checks error code returned from MPI calls, if non-zero it calls the error 6239566063dSJacob Faibussowitsch handler and then returns 6249566063dSJacob Faibussowitsch 6259566063dSJacob Faibussowitsch Synopsis: 6269566063dSJacob Faibussowitsch #include <petscerror.h> 6279566063dSJacob Faibussowitsch void CHKERRMPI(PetscErrorCode ierr) 6289566063dSJacob Faibussowitsch 6299566063dSJacob Faibussowitsch Not Collective 6309566063dSJacob Faibussowitsch 6319566063dSJacob Faibussowitsch Input Parameter: 6329566063dSJacob Faibussowitsch . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 6339566063dSJacob Faibussowitsch 6349566063dSJacob Faibussowitsch Level: deprecated 6359566063dSJacob Faibussowitsch 636667f096bSBarry Smith Note: 637667f096bSBarry Smith Deprecated in favor of `PetscCallMPI()`. This routine behaves identically to it. 638667f096bSBarry Smith 639db781477SPatrick Sanan .seealso: `PetscCallMPI()` 6409566063dSJacob Faibussowitsch M*/ 6419566063dSJacob Faibussowitsch #define CHKERRMPI(...) PetscCallMPI(__VA_ARGS__) 6429566063dSJacob Faibussowitsch 6439566063dSJacob Faibussowitsch /*MC 644989c49a2SBarry Smith PetscCallAbort - Checks error code returned from PETSc function, if non-zero it aborts immediately by calling `MPI_Abort()` 6459566063dSJacob Faibussowitsch 6469566063dSJacob Faibussowitsch Synopsis: 6479566063dSJacob Faibussowitsch #include <petscerror.h> 6489566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm comm, PetscErrorCode ierr) 6499566063dSJacob Faibussowitsch 650c3339decSBarry Smith Collective 6519566063dSJacob Faibussowitsch 6529566063dSJacob Faibussowitsch Input Parameters: 6539566063dSJacob Faibussowitsch + comm - the MPI communicator on which to abort 6549566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 6559566063dSJacob Faibussowitsch 656667f096bSBarry Smith Level: intermediate 657667f096bSBarry Smith 6589566063dSJacob Faibussowitsch Notes: 65987497f52SBarry Smith This macro has identical type and usage semantics to `PetscCall()` with the important caveat 6609566063dSJacob Faibussowitsch that this macro does not return. Instead, if ierr is nonzero it calls the PETSc error handler 66187497f52SBarry Smith and then immediately calls `MPI_Abort()`. It can therefore be used anywhere. 6629566063dSJacob Faibussowitsch 66387497f52SBarry Smith As per `MPI_Abort()` semantics the communicator passed must be valid, although there is currently 66487497f52SBarry Smith no attempt made at handling any potential errors from `MPI_Abort()`. Note that while 66587497f52SBarry Smith `MPI_Abort()` is required to terminate only those processes which reside on comm, it is often 66687497f52SBarry Smith the case that `MPI_Abort()` terminates *all* processes. 6679566063dSJacob Faibussowitsch 6689566063dSJacob Faibussowitsch Example Usage: 6699566063dSJacob Faibussowitsch .vb 6709566063dSJacob Faibussowitsch PetscErrorCode boom(void) { return PETSC_ERR_MEM; } 6719566063dSJacob Faibussowitsch 6729566063dSJacob Faibussowitsch void foo(void) 6739566063dSJacob Faibussowitsch { 6749566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type 6759566063dSJacob Faibussowitsch } 6769566063dSJacob Faibussowitsch 6779566063dSJacob Faibussowitsch double bar(void) 6789566063dSJacob Faibussowitsch { 6799566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type 6809566063dSJacob Faibussowitsch } 6819566063dSJacob Faibussowitsch 6829566063dSJacob Faibussowitsch PetscCallAbort(MPI_COMM_NULL,boom()); // ERROR, communicator should be valid 6839566063dSJacob Faibussowitsch 6849566063dSJacob Faibussowitsch struct baz 6859566063dSJacob Faibussowitsch { 6869566063dSJacob Faibussowitsch baz() 6879566063dSJacob Faibussowitsch { 6889566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,boom()); // OK 6899566063dSJacob Faibussowitsch } 6909566063dSJacob Faibussowitsch 6919566063dSJacob Faibussowitsch ~baz() 6929566063dSJacob Faibussowitsch { 6939566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,boom()); // OK (in fact the only way to handle PETSc errors) 6949566063dSJacob Faibussowitsch } 6959566063dSJacob Faibussowitsch }; 6969566063dSJacob Faibussowitsch .ve 6979566063dSJacob Faibussowitsch 698989c49a2SBarry Smith Fortran Note: 699989c49a2SBarry Smith Use `PetscCallA()`. 700989c49a2SBarry Smith 701989c49a2SBarry Smith Developer Note: 702989c49a2SBarry Smith This should have the same name in Fortran as in C. 703989c49a2SBarry Smith 704db781477SPatrick Sanan .seealso: `SETERRABORT()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, 7055687f061SJacob Faibussowitsch `SETERRQ()`, `CHKMEMQ`, `PetscCallMPI()`, `PetscCallCXXAbort()` 7069566063dSJacob Faibussowitsch M*/ 7079566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 7089566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm, PetscErrorCode); 7099566063dSJacob Faibussowitsch void PetscCallContinue(PetscErrorCode); 7109566063dSJacob Faibussowitsch #else 7119371c9d4SSatish Balay #define PetscCallAbort(comm, ...) \ 7129371c9d4SSatish Balay do { \ 7137da6b3ccSLisandro Dalcin PetscErrorCode ierr_petsc_call_abort_; \ 7143ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 7157da6b3ccSLisandro Dalcin ierr_petsc_call_abort_ = __VA_ARGS__; \ 7163ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_abort_ != PETSC_SUCCESS)) { \ 7173ba16761SJacob Faibussowitsch ierr_petsc_call_abort_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_abort_, PETSC_ERROR_REPEAT, " "); \ 7183ba16761SJacob Faibussowitsch (void)MPI_Abort(comm, (PetscMPIInt)ierr_petsc_call_abort_); \ 7199566063dSJacob Faibussowitsch } \ 7209566063dSJacob Faibussowitsch } while (0) 7219371c9d4SSatish Balay #define PetscCallContinue(...) \ 7229371c9d4SSatish Balay do { \ 7237da6b3ccSLisandro Dalcin PetscErrorCode ierr_petsc_call_continue_; \ 7243ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 7257da6b3ccSLisandro Dalcin ierr_petsc_call_continue_ = __VA_ARGS__; \ 7263ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_continue_ != PETSC_SUCCESS)) { \ 7273ba16761SJacob Faibussowitsch ierr_petsc_call_continue_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_continue_, PETSC_ERROR_REPEAT, " "); \ 7283ba16761SJacob Faibussowitsch (void)ierr_petsc_call_continue_; \ 7293ba16761SJacob Faibussowitsch } \ 7309566063dSJacob Faibussowitsch } while (0) 7319566063dSJacob Faibussowitsch #endif 7329566063dSJacob Faibussowitsch 7339566063dSJacob Faibussowitsch /*MC 7349566063dSJacob Faibussowitsch CHKERRABORT - Checks error code returned from PETSc function. If non-zero it aborts immediately. 7359566063dSJacob Faibussowitsch 7369566063dSJacob Faibussowitsch Synopsis: 7379566063dSJacob Faibussowitsch #include <petscerror.h> 7389566063dSJacob Faibussowitsch void CHKERRABORT(MPI_Comm comm, PetscErrorCode ierr) 7399566063dSJacob Faibussowitsch 7409566063dSJacob Faibussowitsch Not Collective 7419566063dSJacob Faibussowitsch 7429566063dSJacob Faibussowitsch Input Parameters: 7439566063dSJacob Faibussowitsch + comm - the MPI communicator 7449566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 7459566063dSJacob Faibussowitsch 7469566063dSJacob Faibussowitsch Level: deprecated 7479566063dSJacob Faibussowitsch 748667f096bSBarry Smith Note: 749667f096bSBarry Smith Deprecated in favor of `PetscCallAbort()`. This routine behaves identically to it. 750667f096bSBarry Smith 751*af27ebaaSBarry Smith .seealso: `PetscCallAbort()`, `PetscErrorCode` 7529566063dSJacob Faibussowitsch M*/ 7539566063dSJacob Faibussowitsch #define CHKERRABORT(comm, ...) PetscCallAbort(comm, __VA_ARGS__) 7549566063dSJacob Faibussowitsch #define CHKERRCONTINUE(...) PetscCallContinue(__VA_ARGS__) 7559566063dSJacob Faibussowitsch 7569566063dSJacob Faibussowitsch /*MC 75787497f52SBarry Smith CHKERRA - Fortran-only replacement for use of `CHKERRQ()` in the main program, which aborts immediately 758f388eb8bSPatrick Sanan 759f388eb8bSPatrick Sanan Synopsis: 760f388eb8bSPatrick Sanan #include <petscsys.h> 761f388eb8bSPatrick Sanan PetscErrorCode CHKERRA(PetscErrorCode ierr) 762f388eb8bSPatrick Sanan 763f388eb8bSPatrick Sanan Not Collective 764f388eb8bSPatrick Sanan 7652fe279fdSBarry Smith Input Parameter: 766f388eb8bSPatrick Sanan . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 767f388eb8bSPatrick Sanan 76887497f52SBarry Smith Level: deprecated 769f388eb8bSPatrick Sanan 77087497f52SBarry Smith Note: 77187497f52SBarry Smith This macro is rarely needed, normal usage is `PetscCallA()` in the main Fortran program. 772f388eb8bSPatrick Sanan 773989c49a2SBarry Smith Developer Note: 774989c49a2SBarry Smith Why isn't this named `CHKERRABORT()` in Fortran? 775989c49a2SBarry Smith 77687497f52SBarry Smith .seealso: `PetscCall()`, `PetscCallA()`, `PetscCallAbort()`, `CHKERRQ()`, `SETERRA()`, `SETERRQ()`, `SETERRABORT()` 777f388eb8bSPatrick Sanan M*/ 778f388eb8bSPatrick Sanan 7791e4f893aSSatish Balay PETSC_EXTERN PetscBool petscwaitonerrorflg; 7801e4f893aSSatish Balay PETSC_EXTERN PetscBool petscindebugger; 78135f00c14SToby Isaac PETSC_EXTERN PetscBool petscabortmpifinalize; 7827c66cc67SJunchao Zhang 7837c66cc67SJunchao Zhang /*MC 784667f096bSBarry Smith PETSCABORT - Call `MPI_Abort()` with an informative error code 7857c66cc67SJunchao Zhang 7867c66cc67SJunchao Zhang Synopsis: 7877c66cc67SJunchao Zhang #include <petscsys.h> 7887c66cc67SJunchao Zhang PETSCABORT(MPI_Comm comm, PetscErrorCode ierr) 7897c66cc67SJunchao Zhang 7907cdbe19fSJose E. Roman Collective; No Fortran Support 7917c66cc67SJunchao Zhang 7927c66cc67SJunchao Zhang Input Parameters: 79395bd0b28SBarry Smith + comm - An MPI communicator, so that the error can be collective 7947c66cc67SJunchao Zhang - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 7957c66cc67SJunchao Zhang 796bf4d2887SBarry Smith Level: advanced 7977c66cc67SJunchao Zhang 798bf4d2887SBarry Smith Notes: 799989c49a2SBarry Smith If the option `-start_in_debugger` was used then this calls `abort()` to stop the program in the debugger. 800bf4d2887SBarry Smith 801989c49a2SBarry Smith if `PetscCIEnabledPortableErrorOutput` is set, which means the code is running in the PETSc test harness (make test), 802989c49a2SBarry Smith and `comm` is `MPI_COMM_WORLD` it strives to exit cleanly without calling `MPI_Abort()` and instead calling `MPI_Finalize()`. 803660278c0SBarry Smith 804989c49a2SBarry Smith This is currently only used when an error propagates up to the C `main()` program and is detected by a `PetscCall()`, `PetscCallMPI()`, 805989c49a2SBarry Smith or is set in `main()` with `SETERRQ()`. Abort calls such as `SETERRABORT()`, 806989c49a2SBarry Smith `PetscCheckAbort()`, `PetscCallMPIAbort()`, and `PetscCallAbort()` always call `MPI_Abort()` and do not have any special 807989c49a2SBarry Smith handling for the test harness. 808989c49a2SBarry Smith 809989c49a2SBarry Smith Developer Note: 810989c49a2SBarry Smith Should the other abort calls also pass through this call instead of calling `MPI_Abort()` directly? 811989c49a2SBarry Smith 812989c49a2SBarry Smith .seealso: `PetscError()`, `PetscCall()`, `SETERRABORT()`, `PetscCheckAbort()`, `PetscCallMPIAbort()`, `PetscCall()`, `PetscCallMPI()`, 813*af27ebaaSBarry Smith `PetscCallAbort()`, `MPI_Abort()`, `PetscErrorCode` 8147c66cc67SJunchao Zhang M*/ 81529f88feaSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 81629f88feaSJacob Faibussowitsch void PETSCABORT(MPI_Comm, PetscErrorCode); 81729f88feaSJacob Faibussowitsch #else 8189371c9d4SSatish Balay #define PETSCABORT(comm, ...) \ 8199371c9d4SSatish Balay do { \ 8203ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_abort_; \ 8213ba16761SJacob Faibussowitsch if (petscwaitonerrorflg) { ierr_petsc_abort_ = PetscSleep(1000); } \ 8223ba16761SJacob Faibussowitsch if (petscindebugger) { \ 8233ba16761SJacob Faibussowitsch abort(); \ 8243ba16761SJacob Faibussowitsch } else { \ 8257da6b3ccSLisandro Dalcin PetscMPIInt size_; \ 8263ba16761SJacob Faibussowitsch ierr_petsc_abort_ = __VA_ARGS__; \ 8277da6b3ccSLisandro Dalcin MPI_Comm_size(comm, &size_); \ 82835f00c14SToby Isaac if (PetscCIEnabledPortableErrorOutput && (size_ == PetscGlobalSize || petscabortmpifinalize) && ierr_petsc_abort_ != PETSC_ERR_SIG) { \ 8299371c9d4SSatish Balay MPI_Finalize(); \ 8309371c9d4SSatish Balay exit(0); \ 831660278c0SBarry Smith } else if (PetscCIEnabledPortableErrorOutput && PetscGlobalSize == 1) { \ 832660278c0SBarry Smith exit(0); \ 833660278c0SBarry Smith } else { \ 834660278c0SBarry Smith MPI_Abort(comm, (PetscMPIInt)ierr_petsc_abort_); \ 835660278c0SBarry Smith } \ 8363fcd9f07SJacob Faibussowitsch } \ 8377c66cc67SJunchao Zhang } while (0) 83829f88feaSJacob Faibussowitsch #endif 839986eef2eSBarry Smith 8409566063dSJacob Faibussowitsch #ifdef PETSC_CLANGUAGE_CXX 841986eef2eSBarry Smith /*MC 8429566063dSJacob Faibussowitsch PetscCallThrow - Checks error code, if non-zero it calls the C++ error handler which throws 8439566063dSJacob Faibussowitsch an exception 844986eef2eSBarry Smith 845986eef2eSBarry Smith Synopsis: 8469566063dSJacob Faibussowitsch #include <petscerror.h> 8479566063dSJacob Faibussowitsch void PetscCallThrow(PetscErrorCode ierr) 848986eef2eSBarry Smith 849986eef2eSBarry Smith Not Collective 850986eef2eSBarry Smith 8519566063dSJacob Faibussowitsch Input Parameter: 852986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 853986eef2eSBarry Smith 854667f096bSBarry Smith Level: beginner 855667f096bSBarry Smith 856986eef2eSBarry Smith Notes: 857*af27ebaaSBarry Smith Requires PETSc to be configured with clanguage of c++. Throws a std::runtime_error() on error. 858986eef2eSBarry Smith 85987497f52SBarry Smith Once the error handler throws the exception you can use `PetscCallVoid()` which returns without 86087497f52SBarry Smith an error code (bad idea since the error is ignored) or `PetscCallAbort()` to have `MPI_Abort()` 8619566063dSJacob Faibussowitsch called immediately. 8629566063dSJacob Faibussowitsch 863db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscCall()`, `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, 864db781477SPatrick Sanan `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ` 865986eef2eSBarry Smith M*/ 8669371c9d4SSatish Balay #define PetscCallThrow(...) \ 8679371c9d4SSatish Balay do { \ 8683ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 8693ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_throw_ = __VA_ARGS__; \ 8703ba16761SJacob 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); \ 871ffc4695bSBarry Smith } while (0) 87285614651SBarry Smith 873cc26af49SMatthew Knepley /*MC 874cc26af49SMatthew Knepley CHKERRXX - Checks error code, if non-zero it calls the C++ error handler which throws an exception 875cc26af49SMatthew Knepley 876cc26af49SMatthew Knepley Synopsis: 8779566063dSJacob Faibussowitsch #include <petscerror.h> 8783af045c5SBarry Smith void CHKERRXX(PetscErrorCode ierr) 879cc26af49SMatthew Knepley 880eca87e8dSBarry Smith Not Collective 881cc26af49SMatthew Knepley 8829566063dSJacob Faibussowitsch Input Parameter: 8833af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 884cc26af49SMatthew Knepley 8859566063dSJacob Faibussowitsch Level: deprecated 886cc26af49SMatthew Knepley 887667f096bSBarry Smith Note: 888667f096bSBarry Smith Deprecated in favor of `PetscCallThrow()`. This routine behaves identically to it. 889667f096bSBarry Smith 890db781477SPatrick Sanan .seealso: `PetscCallThrow()` 891cc26af49SMatthew Knepley M*/ 8929566063dSJacob Faibussowitsch #define CHKERRXX(...) PetscCallThrow(__VA_ARGS__) 893fd705b32SMatthew Knepley #endif 894fd705b32SMatthew Knepley 8953ba16761SJacob Faibussowitsch #define PetscCallCXX_Private(__SETERR_FUNC__, __COMM__, ...) \ 8963ba16761SJacob Faibussowitsch do { \ 8973ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 8983ba16761SJacob Faibussowitsch try { \ 8993ba16761SJacob Faibussowitsch __VA_ARGS__; \ 9003ba16761SJacob Faibussowitsch } catch (const std::exception &e) { \ 9013ba16761SJacob Faibussowitsch __SETERR_FUNC__(__COMM__, PETSC_ERR_LIB, "%s", e.what()); \ 9023ba16761SJacob Faibussowitsch } \ 9033ba16761SJacob Faibussowitsch } while (0) 9043ba16761SJacob Faibussowitsch 9053f520e80SJunchao Zhang /*MC 9069566063dSJacob Faibussowitsch PetscCallCXX - Checks C++ function calls and if they throw an exception, catch it and then 9079566063dSJacob Faibussowitsch return a PETSc error code 9083f520e80SJunchao Zhang 9093f520e80SJunchao Zhang Synopsis: 9109566063dSJacob Faibussowitsch #include <petscerror.h> 9115687f061SJacob Faibussowitsch void PetscCallCXX(...) noexcept; 9123f520e80SJunchao Zhang 9133f520e80SJunchao Zhang Not Collective 9143f520e80SJunchao Zhang 9159566063dSJacob Faibussowitsch Input Parameter: 9165687f061SJacob Faibussowitsch . __VA_ARGS__ - An arbitrary expression 9175687f061SJacob Faibussowitsch 9185687f061SJacob Faibussowitsch Level: beginner 9193f520e80SJunchao Zhang 9203f520e80SJunchao Zhang Notes: 9215687f061SJacob Faibussowitsch `PetscCallCXX(...)` is a macro replacement for 9229566063dSJacob Faibussowitsch .vb 9239566063dSJacob Faibussowitsch try { 9245687f061SJacob Faibussowitsch __VA_ARGS__; 9259566063dSJacob Faibussowitsch } catch (const std::exception& e) { 9269566063dSJacob Faibussowitsch return ConvertToPetscErrorCode(e); 9279566063dSJacob Faibussowitsch } 9289566063dSJacob Faibussowitsch .ve 9299566063dSJacob Faibussowitsch Due to the fact that it catches any (reasonable) exception, it is essentially noexcept. 9303f520e80SJunchao Zhang 9315687f061SJacob Faibussowitsch If you cannot return a `PetscErrorCode` use `PetscCallCXXAbort()` instead. 9325687f061SJacob Faibussowitsch 9339566063dSJacob Faibussowitsch Example Usage: 9349566063dSJacob Faibussowitsch .vb 9359566063dSJacob Faibussowitsch void foo(void) { throw std::runtime_error("error"); } 9363f520e80SJunchao Zhang 9379566063dSJacob Faibussowitsch void bar() 9389566063dSJacob Faibussowitsch { 939e8952933SJacob Faibussowitsch PetscCallCXX(foo()); // ERROR bar() does not return PetscErrorCode 9409566063dSJacob Faibussowitsch } 9419566063dSJacob Faibussowitsch 9429566063dSJacob Faibussowitsch PetscErrorCode baz() 9439566063dSJacob Faibussowitsch { 9449566063dSJacob Faibussowitsch PetscCallCXX(foo()); // OK 9459566063dSJacob Faibussowitsch 9469566063dSJacob Faibussowitsch PetscCallCXX( 9479566063dSJacob Faibussowitsch bar(); 948d5b43468SJose E. Roman foo(); // OK multiple statements allowed 9499566063dSJacob Faibussowitsch ); 9509566063dSJacob Faibussowitsch } 9519566063dSJacob Faibussowitsch 9529566063dSJacob Faibussowitsch struct bop 9539566063dSJacob Faibussowitsch { 9549566063dSJacob Faibussowitsch bop() 9559566063dSJacob Faibussowitsch { 956e8952933SJacob Faibussowitsch PetscCallCXX(foo()); // ERROR returns PetscErrorCode, cannot be used in constructors 9579566063dSJacob Faibussowitsch } 9589566063dSJacob Faibussowitsch }; 9599566063dSJacob Faibussowitsch 960e8952933SJacob Faibussowitsch // ERROR contains do-while, cannot be used as function-try block 9619566063dSJacob Faibussowitsch PetscErrorCode qux() PetscCallCXX( 9629566063dSJacob Faibussowitsch bar(); 9639566063dSJacob Faibussowitsch baz(); 9649566063dSJacob Faibussowitsch foo(); 9659566063dSJacob Faibussowitsch return 0; 9669566063dSJacob Faibussowitsch ) 9679566063dSJacob Faibussowitsch .ve 9689566063dSJacob Faibussowitsch 9695687f061SJacob Faibussowitsch .seealso: `PetscCallCXXAbort()`, `PetscCallThrow()`, `SETERRQ()`, `PetscCall()`, 9705687f061SJacob Faibussowitsch `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 9715687f061SJacob Faibussowitsch `PetscError()`, `CHKMEMQ` 9723f520e80SJunchao Zhang M*/ 9733ba16761SJacob Faibussowitsch #define PetscCallCXX(...) PetscCallCXX_Private(SETERRQ, PETSC_COMM_SELF, __VA_ARGS__) 9745687f061SJacob Faibussowitsch 9755687f061SJacob Faibussowitsch /*MC 9765687f061SJacob Faibussowitsch PetscCallCXXAbort - Like `PetscCallCXX()` but calls `MPI_Abort()` instead of returning an 9775687f061SJacob Faibussowitsch error-code 9785687f061SJacob Faibussowitsch 9795687f061SJacob Faibussowitsch Synopsis: 9805687f061SJacob Faibussowitsch #include <petscerror.h> 9815687f061SJacob Faibussowitsch void PetscCallCXXAbort(MPI_Comm comm, ...) noexcept; 9825687f061SJacob Faibussowitsch 98320f4b53cSBarry Smith Collective; No Fortran Support 9845687f061SJacob Faibussowitsch 9855687f061SJacob Faibussowitsch Input Parameters: 9865687f061SJacob Faibussowitsch + comm - The MPI communicator to abort on 9875687f061SJacob Faibussowitsch - __VA_ARGS__ - An arbitrary expression 9885687f061SJacob Faibussowitsch 9895687f061SJacob Faibussowitsch Level: beginner 9905687f061SJacob Faibussowitsch 9915687f061SJacob Faibussowitsch Notes: 9925687f061SJacob Faibussowitsch This macro may be used to check C++ expressions for exceptions in cases where you cannot 9935687f061SJacob Faibussowitsch return an error code. This includes constructors, destructors, copy/move assignment functions 9945687f061SJacob Faibussowitsch or constructors among others. 9955687f061SJacob Faibussowitsch 9965687f061SJacob Faibussowitsch If an exception is caught, the macro calls `SETERRABORT()` on `comm`. The exception must 9975687f061SJacob Faibussowitsch derive from `std::exception` in order to be caught. 9985687f061SJacob Faibussowitsch 9995687f061SJacob Faibussowitsch If the routine _can_ return an error-code it is highly advised to use `PetscCallCXX()` 10005687f061SJacob Faibussowitsch instead. 10015687f061SJacob Faibussowitsch 10025687f061SJacob Faibussowitsch See `PetscCallCXX()` for additional discussion. 10035687f061SJacob Faibussowitsch 10045687f061SJacob Faibussowitsch Example Usage: 10055687f061SJacob Faibussowitsch .vb 10065687f061SJacob Faibussowitsch class Foo 10075687f061SJacob Faibussowitsch { 10085687f061SJacob Faibussowitsch std::vector<int> data_; 10095687f061SJacob Faibussowitsch 10105687f061SJacob Faibussowitsch public: 10115687f061SJacob Faibussowitsch // normally std::vector::reserve() may raise an exception, but since we handle it with 10125687f061SJacob Faibussowitsch // PetscCallCXXAbort() we may mark this routine as noexcept! 10135687f061SJacob Faibussowitsch Foo() noexcept 10145687f061SJacob Faibussowitsch { 10155687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, data_.reserve(10)); 10165687f061SJacob Faibussowitsch } 10175687f061SJacob Faibussowitsch }; 10185687f061SJacob Faibussowitsch 10195687f061SJacob Faibussowitsch std::vector<int> bar() 10205687f061SJacob Faibussowitsch { 10215687f061SJacob Faibussowitsch std::vector<int> v; 10225687f061SJacob Faibussowitsch 10235687f061SJacob Faibussowitsch PetscFunctionBegin; 10245687f061SJacob Faibussowitsch // OK! 10255687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1)); 10265687f061SJacob Faibussowitsch PetscFunctionReturn(v); 10275687f061SJacob Faibussowitsch } 10285687f061SJacob Faibussowitsch 10295687f061SJacob Faibussowitsch PetscErrorCode baz() 10305687f061SJacob Faibussowitsch { 10315687f061SJacob Faibussowitsch std::vector<int> v; 10325687f061SJacob Faibussowitsch 10335687f061SJacob Faibussowitsch PetscFunctionBegin; 10345687f061SJacob Faibussowitsch // WRONG! baz() returns a PetscErrorCode, prefer PetscCallCXX() instead 10355687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1)); 10363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10375687f061SJacob Faibussowitsch } 10385687f061SJacob Faibussowitsch .ve 10395687f061SJacob Faibussowitsch 10405687f061SJacob Faibussowitsch .seealso: `PetscCallCXX()`, `SETERRABORT()`, `PetscCallAbort()` 10415687f061SJacob Faibussowitsch M*/ 10423ba16761SJacob Faibussowitsch #define PetscCallCXXAbort(comm, ...) PetscCallCXX_Private(SETERRABORT, comm, __VA_ARGS__) 10433f520e80SJunchao Zhang 104430de9b25SBarry Smith /*MC 10459566063dSJacob Faibussowitsch CHKERRCXX - Checks C++ function calls and if they throw an exception, catch it and then 10469566063dSJacob Faibussowitsch return a PETSc error code 10479566063dSJacob Faibussowitsch 10489566063dSJacob Faibussowitsch Synopsis: 10499566063dSJacob Faibussowitsch #include <petscerror.h> 10509566063dSJacob Faibussowitsch void CHKERRCXX(func) noexcept; 10519566063dSJacob Faibussowitsch 10529566063dSJacob Faibussowitsch Not Collective 10539566063dSJacob Faibussowitsch 10549566063dSJacob Faibussowitsch Input Parameter: 10559566063dSJacob Faibussowitsch . func - C++ function calls 10569566063dSJacob Faibussowitsch 10579566063dSJacob Faibussowitsch Level: deprecated 10589566063dSJacob Faibussowitsch 1059667f096bSBarry Smith Note: 1060667f096bSBarry Smith Deprecated in favor of `PetscCallCXX()`. This routine behaves identically to it. 1061667f096bSBarry Smith 1062db781477SPatrick Sanan .seealso: `PetscCallCXX()` 10639566063dSJacob Faibussowitsch M*/ 10649566063dSJacob Faibussowitsch #define CHKERRCXX(...) PetscCallCXX(__VA_ARGS__) 10659566063dSJacob Faibussowitsch 10669566063dSJacob Faibussowitsch /*MC 106730de9b25SBarry Smith CHKMEMQ - Checks the memory for corruption, calls error handler if any is detected 106830de9b25SBarry Smith 106930de9b25SBarry Smith Synopsis: 1070aaa7dc30SBarry Smith #include <petscsys.h> 107191d3bdf4SKris Buschelman CHKMEMQ; 107230de9b25SBarry Smith 1073eca87e8dSBarry Smith Not Collective 1074eca87e8dSBarry Smith 107530de9b25SBarry Smith Level: beginner 107630de9b25SBarry Smith 107730de9b25SBarry Smith Notes: 1078*af27ebaaSBarry Smith We recommend using Valgrind <https://petsc.org/release/faq/#valgrind> or for NVIDIA CUDA systems 1079*af27ebaaSBarry Smith <https://docs.nvidia.com/cuda/cuda-memcheck/index.html> for finding memory problems. The ``CHKMEMQ`` macro is useful on systems that 10805ed36255SBarry Smith do not have valgrind, but is not as good as valgrind or cuda-memcheck. 10811957e957SBarry Smith 1082667f096bSBarry Smith Must run with the option `-malloc_debug` (`-malloc_test` in debug mode; or if `PetscMallocSetDebug()` called) to enable this option 108330de9b25SBarry Smith 108430de9b25SBarry Smith Once the error handler is called the calling function is then returned from with the given error code. 108530de9b25SBarry Smith 108630de9b25SBarry Smith By defaults prints location where memory that is corrupted was allocated. 108730de9b25SBarry Smith 108887497f52SBarry Smith Use `CHKMEMA` for functions that return void 1089f621e05eSBarry Smith 1090db781477SPatrick Sanan .seealso: `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `SETERRQ()`, `PetscMallocValidate()` 109130de9b25SBarry Smith M*/ 10926d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 1093064a246eSJacob Faibussowitsch #define CHKMEMQ 1094064a246eSJacob Faibussowitsch #define CHKMEMA 10956d210af2SJacob Faibussowitsch #else 10969371c9d4SSatish Balay #define CHKMEMQ \ 10979371c9d4SSatish Balay do { \ 10983ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_memq_ = PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__); \ 10993ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_memq_ != PETSC_SUCCESS)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_memq_, PETSC_ERROR_REPEAT, " "); \ 110086d09637SLisandro Dalcin } while (0) 11016d210af2SJacob Faibussowitsch #define CHKMEMA PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__) 1102064a246eSJacob Faibussowitsch #endif 11039566063dSJacob Faibussowitsch 1104668f157eSBarry Smith /*E 1105668f157eSBarry Smith PetscErrorType - passed to the PETSc error handling routines indicating if this is the first or a later call to the error handlers 1106668f157eSBarry Smith 1107668f157eSBarry Smith Level: advanced 1108668f157eSBarry Smith 1109667f096bSBarry Smith Note: 111087497f52SBarry Smith `PETSC_ERROR_IN_CXX` indicates the error was detected in C++ and an exception should be generated 1111d736bfebSBarry Smith 1112667f096bSBarry Smith Developer Note: 1113667f096bSBarry Smith This is currently used to decide when to print the detailed information about the run in `PetscTraceBackErrorHandler()` 1114668f157eSBarry Smith 111587497f52SBarry Smith .seealso: `PetscError()`, `SETERRQ()` 1116668f157eSBarry Smith E*/ 11179371c9d4SSatish Balay typedef enum { 11189371c9d4SSatish Balay PETSC_ERROR_INITIAL = 0, 11199371c9d4SSatish Balay PETSC_ERROR_REPEAT = 1, 11209371c9d4SSatish Balay PETSC_ERROR_IN_CXX = 2 11219371c9d4SSatish Balay } PetscErrorType; 11224b209cf6SBarry Smith 1123eb9e708aSLisandro Dalcin #if defined(__clang_analyzer__) 1124eb9e708aSLisandro Dalcin __attribute__((analyzer_noreturn)) 1125eb9e708aSLisandro Dalcin #endif 11269371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode 11279371c9d4SSatish Balay PetscError(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, ...) PETSC_ATTRIBUTE_COLD PETSC_ATTRIBUTE_FORMAT(7, 8); 1128eb9e708aSLisandro Dalcin 1129014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscErrorPrintfInitialize(void); 11303ba16761SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscErrorMessage(PetscErrorCode, const char *[], char **); 1131d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscTraceBackErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1132d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscIgnoreErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1133d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscEmacsClientErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1134d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscMPIAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1135d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1136d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAttachDebuggerErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1137d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscReturnErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1138efca3c55SSatish Balay PETSC_EXTERN PetscErrorCode PetscPushErrorHandler(PetscErrorCode (*handler)(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *), void *); 1139014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopErrorHandler(void); 11408d359177SBarry Smith PETSC_EXTERN PetscErrorCode PetscSignalHandlerDefault(int, void *); 1141014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPushSignalHandler(PetscErrorCode (*)(int, void *), void *); 1142014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopSignalHandler(void); 114328559dc8SJed Brown PETSC_EXTERN PetscErrorCode PetscCheckPointerSetIntensity(PetscInt); 1144c2a741eeSJunchao Zhang PETSC_EXTERN void PetscSignalSegvCheckPointerOrMpi(void); 1145edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 13, 0, "PetscSignalSegvCheckPointerOrMpi()", ) static inline void PetscSignalSegvCheckPointer(void) 1146d71ae5a4SJacob Faibussowitsch { 11479371c9d4SSatish Balay PetscSignalSegvCheckPointerOrMpi(); 11489371c9d4SSatish Balay } 1149329f5518SBarry Smith 1150639ff905SBarry Smith /*MC 1151639ff905SBarry Smith PetscErrorPrintf - Prints error messages. 1152639ff905SBarry Smith 1153639ff905SBarry Smith Synopsis: 1154aaa7dc30SBarry Smith #include <petscsys.h> 1155639ff905SBarry Smith PetscErrorCode (*PetscErrorPrintf)(const char format[],...); 1156639ff905SBarry Smith 11577cdbe19fSJose E. Roman Not Collective; No Fortran Support 11587cdbe19fSJose E. Roman 1159f899ff85SJose E. Roman Input Parameter: 1160cd05f99aSJacob Faibussowitsch . format - the usual `printf()` format string 1161639ff905SBarry Smith 1162639ff905SBarry Smith Options Database Keys: 11631957e957SBarry Smith + -error_output_stdout - cause error messages to be printed to stdout instead of the (default) stderr 1164e1bc860dSBarry Smith - -error_output_none - to turn off all printing of error messages (does not change the way the error is handled.) 1165639ff905SBarry Smith 1166cf53795eSBarry Smith Level: developer 1167cf53795eSBarry Smith 116895452b02SPatrick Sanan Notes: 116995452b02SPatrick Sanan Use 1170667f096bSBarry Smith .vb 117195bd0b28SBarry Smith PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the error is handled) and 1172667f096bSBarry Smith PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on or you can use your own function 1173667f096bSBarry Smith .ve 1174639ff905SBarry Smith Use 1175667f096bSBarry Smith .vb 117687497f52SBarry Smith `PETSC_STDERR` = FILE* obtained from a file open etc. to have stderr printed to the file. 117787497f52SBarry Smith `PETSC_STDOUT` = FILE* obtained from a file open etc. to have stdout printed to the file. 1178667f096bSBarry Smith .ve 1179639ff905SBarry Smith Use 1180*af27ebaaSBarry Smith .vb 118187497f52SBarry Smith `PetscPushErrorHandler()` to provide your own error handler that determines what kind of messages to print 1182*af27ebaaSBarry Smith .ve 1183639ff905SBarry Smith 1184db781477SPatrick Sanan .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscHelpPrintf()`, `PetscPrintf()`, `PetscPushErrorHandler()`, `PetscVFPrintf()`, `PetscHelpPrintf()` 1185639ff905SBarry Smith M*/ 11863ca90d2dSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode (*PetscErrorPrintf)(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2); 1187639ff905SBarry Smith 1188cf0818bdSBarry Smith /*E 1189cf0818bdSBarry Smith PetscFPTrap - types of floating point exceptions that may be trapped 1190cf0818bdSBarry Smith 119187497f52SBarry Smith Currently only `PETSC_FP_TRAP_OFF` and `PETSC_FP_TRAP_ON` are handled. All others are treated as `PETSC_FP_TRAP_ON`. 1192cf0818bdSBarry Smith 1193cf0818bdSBarry Smith Level: intermediate 1194cf0818bdSBarry Smith 11957de69702SBarry Smith .seealso: `PetscSetFPTrap()`, `PetscFPTrapPush()` 1196cf0818bdSBarry Smith E*/ 11979371c9d4SSatish Balay typedef enum { 11989371c9d4SSatish Balay PETSC_FP_TRAP_OFF = 0, 11999371c9d4SSatish Balay PETSC_FP_TRAP_INDIV = 1, 12009371c9d4SSatish Balay PETSC_FP_TRAP_FLTOPERR = 2, 12019371c9d4SSatish Balay PETSC_FP_TRAP_FLTOVF = 4, 12029371c9d4SSatish Balay PETSC_FP_TRAP_FLTUND = 8, 12039371c9d4SSatish Balay PETSC_FP_TRAP_FLTDIV = 16, 12049371c9d4SSatish Balay PETSC_FP_TRAP_FLTINEX = 32 12059371c9d4SSatish Balay } PetscFPTrap; 1206bd2b07b1SBarry 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) 1207014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSetFPTrap(PetscFPTrap); 1208014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPush(PetscFPTrap); 1209014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPop(void); 1210aba4c478SBarry Smith PETSC_EXTERN PetscErrorCode PetscDetermineInitialFPTrap(void); 121154a8ef01SBarry Smith 12123a40ed3dSBarry Smith /* 12133a40ed3dSBarry Smith Allows the code to build a stack frame as it runs 12143a40ed3dSBarry Smith */ 12153a40ed3dSBarry Smith 121699cd645aSJed Brown #define PETSCSTACKSIZE 64 12173a40ed3dSBarry Smith typedef struct { 12180e33f6ddSBarry Smith const char *function[PETSCSTACKSIZE]; 12190e33f6ddSBarry Smith const char *file[PETSCSTACKSIZE]; 1220184914b5SBarry Smith int line[PETSCSTACKSIZE]; 1221362febeeSStefano Zampini int petscroutine[PETSCSTACKSIZE]; /* 0 external called from petsc, 1 petsc functions, 2 petsc user functions */ 1222184914b5SBarry Smith int currentsize; 1223a2f94806SJed Brown int hotdepth; 12244be741a6SBarry Smith PetscBool check; /* option to check for correct Push/Pop semantics, true for default petscstack but not other stacks */ 12253a40ed3dSBarry Smith } PetscStack; 1226dfb7d7afSStefano Zampini #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY) 122727104ee2SJacob Faibussowitsch PETSC_EXTERN PetscStack petscstack; 122827104ee2SJacob Faibussowitsch #endif 12293a40ed3dSBarry Smith 12305d12eec7SSatish Balay #if defined(PETSC_SERIALIZE_FUNCTIONS) 12315d12eec7SSatish Balay #include <petsc/private/petscfptimpl.h> 12325d12eec7SSatish Balay /* 12335d12eec7SSatish Balay Registers the current function into the global function pointer to function name table 12345d12eec7SSatish Balay 12355d12eec7SSatish Balay Have to fix this to handle errors but cannot return error since used in PETSC_VIEWER_DRAW_() etc 12365d12eec7SSatish Balay */ 12379371c9d4SSatish Balay #define PetscRegister__FUNCT__() \ 12389371c9d4SSatish Balay do { \ 12395d12eec7SSatish Balay static PetscBool __chked = PETSC_FALSE; \ 12405d12eec7SSatish Balay if (!__chked) { \ 12419371c9d4SSatish Balay void *ptr; \ 12423ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscDLSym(NULL, PETSC_FUNCTION_NAME, &ptr)); \ 12435d12eec7SSatish Balay __chked = PETSC_TRUE; \ 12449371c9d4SSatish Balay } \ 12459371c9d4SSatish Balay } while (0) 12465d12eec7SSatish Balay #else 12475d12eec7SSatish Balay #define PetscRegister__FUNCT__() 12485d12eec7SSatish Balay #endif 12495d12eec7SSatish Balay 1250eae3dc7dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) || defined(__clang_analyzer__) 12516d210af2SJacob Faibussowitsch #define PetscStackPushNoCheck(funct, petsc_routine, hot) 125237154d10SBarry Smith #define PetscStackUpdateLine 1253792fecdfSBarry Smith #define PetscStackPushExternal(funct) 12546d210af2SJacob Faibussowitsch #define PetscStackPopNoCheck 12556d210af2SJacob Faibussowitsch #define PetscStackClearTop 12566d210af2SJacob Faibussowitsch #define PetscFunctionBegin 12576d210af2SJacob Faibussowitsch #define PetscFunctionBeginUser 12586d210af2SJacob Faibussowitsch #define PetscFunctionBeginHot 12590a57284eSJacob Faibussowitsch #define PetscFunctionReturn(...) return __VA_ARGS__ 12606d210af2SJacob Faibussowitsch #define PetscFunctionReturnVoid() return 12616d210af2SJacob Faibussowitsch #define PetscStackPop 12626d210af2SJacob Faibussowitsch #define PetscStackPush(f) 1263dfb7d7afSStefano Zampini #elif defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY) 1264660278c0SBarry Smith 12659371c9d4SSatish Balay #define PetscStackPush_Private(stack__, file__, func__, line__, petsc_routine__, hot__) \ 12669371c9d4SSatish Balay do { \ 12675a96b57dSJacob Faibussowitsch if (stack__.currentsize < PETSCSTACKSIZE) { \ 12685a96b57dSJacob Faibussowitsch stack__.function[stack__.currentsize] = func__; \ 1269ef1023bdSBarry Smith if (petsc_routine__) { \ 1270ef1023bdSBarry Smith stack__.file[stack__.currentsize] = file__; \ 12715a96b57dSJacob Faibussowitsch stack__.line[stack__.currentsize] = line__; \ 1272ef1023bdSBarry Smith } else { \ 1273648bc8c4SBarry Smith stack__.file[stack__.currentsize] = PETSC_NULLPTR; \ 1274ef1023bdSBarry Smith stack__.line[stack__.currentsize] = 0; \ 1275ef1023bdSBarry Smith } \ 12765a96b57dSJacob Faibussowitsch stack__.petscroutine[stack__.currentsize] = petsc_routine__; \ 12775a96b57dSJacob Faibussowitsch } \ 12785a96b57dSJacob Faibussowitsch ++stack__.currentsize; \ 12795a96b57dSJacob Faibussowitsch stack__.hotdepth += (hot__ || stack__.hotdepth); \ 12805a96b57dSJacob Faibussowitsch } while (0) 12815a96b57dSJacob Faibussowitsch 12824be741a6SBarry Smith /* uses PetscCheckAbort() because may be used in a function that does not return an error code */ 12839371c9d4SSatish Balay #define PetscStackPop_Private(stack__, func__) \ 12849371c9d4SSatish Balay do { \ 12854be741a6SBarry 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__); \ 12865a96b57dSJacob Faibussowitsch if (--stack__.currentsize < PETSCSTACKSIZE) { \ 12879371c9d4SSatish 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", \ 12889371c9d4SSatish Balay stack__.function[stack__.currentsize], stack__.file[stack__.currentsize], stack__.line[stack__.currentsize], func__, __FILE__, __LINE__); \ 12895a96b57dSJacob Faibussowitsch stack__.function[stack__.currentsize] = PETSC_NULLPTR; \ 12905a96b57dSJacob Faibussowitsch stack__.file[stack__.currentsize] = PETSC_NULLPTR; \ 12915a96b57dSJacob Faibussowitsch stack__.line[stack__.currentsize] = 0; \ 12925a96b57dSJacob Faibussowitsch stack__.petscroutine[stack__.currentsize] = 0; \ 12935a96b57dSJacob Faibussowitsch } \ 12945a96b57dSJacob Faibussowitsch stack__.hotdepth = PetscMax(stack__.hotdepth - 1, 0); \ 12955a96b57dSJacob Faibussowitsch } while (0) 12965a96b57dSJacob Faibussowitsch 1297586f9135SBarry Smith /*MC 1298586f9135SBarry Smith PetscStackPushNoCheck - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is 1299586f9135SBarry Smith currently in the source code. 1300586f9135SBarry Smith 1301586f9135SBarry Smith Synopsis: 1302586f9135SBarry Smith #include <petscsys.h> 1303586f9135SBarry Smith void PetscStackPushNoCheck(char *funct,int petsc_routine,PetscBool hot); 1304586f9135SBarry Smith 13057cdbe19fSJose E. Roman Not Collective 13067cdbe19fSJose E. Roman 1307586f9135SBarry Smith Input Parameters: 1308586f9135SBarry Smith + funct - the function name 1309586f9135SBarry Smith . petsc_routine - 2 user function, 1 PETSc function, 0 some other function 1310586f9135SBarry Smith - hot - indicates that the function may be called often so expensive error checking should be turned off inside the function 1311586f9135SBarry Smith 1312586f9135SBarry Smith Level: developer 1313586f9135SBarry Smith 1314586f9135SBarry Smith Notes: 1315586f9135SBarry 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 131687497f52SBarry 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 1317586f9135SBarry Smith help debug the problem. 1318586f9135SBarry Smith 1319ef1023bdSBarry Smith This version does not check the memory corruption (an expensive operation), use `PetscStackPush()` to check the memory. 1320ef1023bdSBarry Smith 1321792fecdfSBarry Smith Use `PetscStackPushExternal()` for a function call that is about to be made to a non-PETSc or user function (such as BLAS etc). 1322ef1023bdSBarry Smith 1323586f9135SBarry Smith The default stack is a global variable called `petscstack`. 1324586f9135SBarry Smith 1325586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1326ef1023bdSBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPush()`, `PetscStackPop`, 1327792fecdfSBarry Smith `PetscStackPushExternal()` 1328586f9135SBarry Smith M*/ 13299371c9d4SSatish Balay #define PetscStackPushNoCheck(funct, petsc_routine, hot) \ 13309371c9d4SSatish Balay do { \ 1331e04113cfSBarry Smith PetscStackSAWsTakeAccess(); \ 13325a96b57dSJacob Faibussowitsch PetscStackPush_Private(petscstack, __FILE__, funct, __LINE__, petsc_routine, hot); \ 1333e04113cfSBarry Smith PetscStackSAWsGrantAccess(); \ 1334441dd030SJed Brown } while (0) 1335441dd030SJed Brown 1336586f9135SBarry Smith /*MC 133787497f52SBarry Smith PetscStackUpdateLine - in a function that has a `PetscFunctionBegin` or `PetscFunctionBeginUser` updates the stack line number to the 133837154d10SBarry Smith current line number. 133937154d10SBarry Smith 134037154d10SBarry Smith Synopsis: 134137154d10SBarry Smith #include <petscsys.h> 134237154d10SBarry Smith void PetscStackUpdateLine 134337154d10SBarry Smith 13447cdbe19fSJose E. Roman Not Collective 13457cdbe19fSJose E. Roman 134637154d10SBarry Smith Level: developer 134737154d10SBarry Smith 134837154d10SBarry Smith Notes: 134987497f52SBarry Smith Using `PetscCall()` and friends automatically handles this process 135087497f52SBarry Smith 135137154d10SBarry 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 135237154d10SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 135337154d10SBarry Smith help debug the problem. 135437154d10SBarry Smith 135595bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 135637154d10SBarry Smith 135737154d10SBarry Smith This is used by `PetscCall()` and is otherwise not like to be needed 135837154d10SBarry Smith 135937154d10SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`, `PetscCall()` 136037154d10SBarry Smith M*/ 136137154d10SBarry Smith #define PetscStackUpdateLine \ 13623ba16761SJacob Faibussowitsch do { \ 13633ba16761SJacob Faibussowitsch if (petscstack.currentsize > 0 && petscstack.function[petscstack.currentsize - 1] == PETSC_FUNCTION_NAME) { petscstack.line[petscstack.currentsize - 1] = __LINE__; } \ 13643ba16761SJacob Faibussowitsch } while (0) 136537154d10SBarry Smith 136637154d10SBarry Smith /*MC 1367792fecdfSBarry Smith PetscStackPushExternal - Pushes a new function name onto the PETSc default stack that tracks where the running program is 1368ef1023bdSBarry Smith currently in the source code. Does not include the filename or line number since this is called by the calling routine 1369ef1023bdSBarry Smith for non-PETSc or user functions. 1370ef1023bdSBarry Smith 1371ef1023bdSBarry Smith Synopsis: 1372ef1023bdSBarry Smith #include <petscsys.h> 1373792fecdfSBarry Smith void PetscStackPushExternal(char *funct); 1374ef1023bdSBarry Smith 13757cdbe19fSJose E. Roman Not Collective 13767cdbe19fSJose E. Roman 13772fe279fdSBarry Smith Input Parameter: 1378ef1023bdSBarry Smith . funct - the function name 1379ef1023bdSBarry Smith 1380ef1023bdSBarry Smith Level: developer 1381ef1023bdSBarry Smith 1382ef1023bdSBarry Smith Notes: 138387497f52SBarry Smith Using `PetscCallExternal()` and friends automatically handles this process 138487497f52SBarry Smith 1385ef1023bdSBarry 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 1386ef1023bdSBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1387ef1023bdSBarry Smith help debug the problem. 1388ef1023bdSBarry Smith 1389ef1023bdSBarry Smith The default stack is a global variable called `petscstack`. 1390ef1023bdSBarry Smith 1391ef1023bdSBarry Smith This is to be used when calling an external package function such as a BLAS function. 1392ef1023bdSBarry Smith 1393ef1023bdSBarry Smith This also updates the stack line number for the current stack function. 1394ef1023bdSBarry Smith 1395ef1023bdSBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1396ef1023bdSBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1397ef1023bdSBarry Smith M*/ 13989371c9d4SSatish Balay #define PetscStackPushExternal(funct) \ 13999371c9d4SSatish Balay do { \ 14009371c9d4SSatish Balay PetscStackUpdateLine; \ 14019371c9d4SSatish Balay PetscStackPushNoCheck(funct, 0, PETSC_TRUE); \ 1402a8f51744SPierre Jolivet } while (0) 1403ef1023bdSBarry Smith 1404ef1023bdSBarry Smith /*MC 1405586f9135SBarry Smith PetscStackPopNoCheck - Pops a function name from the PETSc default stack that tracks where the running program is 1406586f9135SBarry Smith currently in the source code. 1407586f9135SBarry Smith 1408586f9135SBarry Smith Synopsis: 1409586f9135SBarry Smith #include <petscsys.h> 1410586f9135SBarry Smith void PetscStackPopNoCheck(char *funct); 1411586f9135SBarry Smith 14127cdbe19fSJose E. Roman Not Collective 14137cdbe19fSJose E. Roman 1414586f9135SBarry Smith Input Parameter: 1415586f9135SBarry Smith . funct - the function name 1416586f9135SBarry Smith 1417586f9135SBarry Smith Level: developer 1418586f9135SBarry Smith 1419586f9135SBarry Smith Notes: 142087497f52SBarry Smith Using `PetscCall()`, `PetscCallExternal()`, `PetscCallBack()` and friends negates the need to call this 142187497f52SBarry Smith 1422586f9135SBarry 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 1423586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1424586f9135SBarry Smith help debug the problem. 1425586f9135SBarry Smith 142695bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 1427586f9135SBarry Smith 1428586f9135SBarry Smith Developer Note: 1429586f9135SBarry Smith `PetscStackPopNoCheck()` takes a function argument while `PetscStackPop` does not, this difference is likely just historical. 1430586f9135SBarry Smith 1431586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1432586f9135SBarry Smith M*/ 14339371c9d4SSatish Balay #define PetscStackPopNoCheck(funct) \ 14349371c9d4SSatish Balay do { \ 1435362febeeSStefano Zampini PetscStackSAWsTakeAccess(); \ 14365a96b57dSJacob Faibussowitsch PetscStackPop_Private(petscstack, funct); \ 1437362febeeSStefano Zampini PetscStackSAWsGrantAccess(); \ 1438362febeeSStefano Zampini } while (0) 1439362febeeSStefano Zampini 14409371c9d4SSatish Balay #define PetscStackClearTop \ 14419371c9d4SSatish Balay do { \ 1442e04113cfSBarry Smith PetscStackSAWsTakeAccess(); \ 14439371c9d4SSatish Balay if (petscstack.currentsize > 0 && --petscstack.currentsize < PETSCSTACKSIZE) { \ 144427104ee2SJacob Faibussowitsch petscstack.function[petscstack.currentsize] = PETSC_NULLPTR; \ 144527104ee2SJacob Faibussowitsch petscstack.file[petscstack.currentsize] = PETSC_NULLPTR; \ 144627104ee2SJacob Faibussowitsch petscstack.line[petscstack.currentsize] = 0; \ 144727104ee2SJacob Faibussowitsch petscstack.petscroutine[petscstack.currentsize] = 0; \ 1448441dd030SJed Brown } \ 144927104ee2SJacob Faibussowitsch petscstack.hotdepth = PetscMax(petscstack.hotdepth - 1, 0); \ 1450e04113cfSBarry Smith PetscStackSAWsGrantAccess(); \ 1451441dd030SJed Brown } while (0) 1452441dd030SJed Brown 145330de9b25SBarry Smith /*MC 14541957e957SBarry Smith PetscFunctionBegin - First executable line of each PETSc function, used for error handling. Final 145587497f52SBarry Smith line of PETSc functions should be `PetscFunctionReturn`(0); 145630de9b25SBarry Smith 145730de9b25SBarry Smith Synopsis: 1458aaa7dc30SBarry Smith #include <petscsys.h> 145930de9b25SBarry Smith void PetscFunctionBegin; 146030de9b25SBarry Smith 146120f4b53cSBarry Smith Not Collective; No Fortran Support 1462eca87e8dSBarry Smith 146330de9b25SBarry Smith Usage: 146430de9b25SBarry Smith .vb 146530de9b25SBarry Smith int something; 146630de9b25SBarry Smith 146730de9b25SBarry Smith PetscFunctionBegin; 146830de9b25SBarry Smith .ve 146930de9b25SBarry Smith 147030de9b25SBarry Smith Level: developer 147130de9b25SBarry Smith 147220f4b53cSBarry Smith Note: 147320f4b53cSBarry Smith Use `PetscFunctionBeginUser` for application codes. 147420f4b53cSBarry Smith 1475586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()` 147630de9b25SBarry Smith 147730de9b25SBarry Smith M*/ 14789371c9d4SSatish Balay #define PetscFunctionBegin \ 14799371c9d4SSatish Balay do { \ 1480362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_FALSE); \ 1481a2f94806SJed Brown PetscRegister__FUNCT__(); \ 1482a2f94806SJed Brown } while (0) 1483a2f94806SJed Brown 1484a2f94806SJed Brown /*MC 148587497f52SBarry Smith PetscFunctionBeginHot - Substitute for `PetscFunctionBegin` to be used in functions that are called in 1486a2f94806SJed Brown performance-critical circumstances. Use of this function allows for lighter profiling by default. 1487a2f94806SJed Brown 1488a2f94806SJed Brown Synopsis: 1489aaa7dc30SBarry Smith #include <petscsys.h> 1490a2f94806SJed Brown void PetscFunctionBeginHot; 1491a2f94806SJed Brown 149220f4b53cSBarry Smith Not Collective; No Fortran Support 1493a2f94806SJed Brown 1494a2f94806SJed Brown Usage: 1495a2f94806SJed Brown .vb 1496a2f94806SJed Brown int something; 1497a2f94806SJed Brown 1498a2f94806SJed Brown PetscFunctionBeginHot; 1499a2f94806SJed Brown .ve 1500a2f94806SJed Brown 1501a2f94806SJed Brown Level: developer 1502a2f94806SJed Brown 1503586f9135SBarry Smith .seealso: `PetscFunctionBegin`, `PetscFunctionReturn()`, `PetscStackPushNoCheck()` 1504a2f94806SJed Brown 1505a2f94806SJed Brown M*/ 15069371c9d4SSatish Balay #define PetscFunctionBeginHot \ 15079371c9d4SSatish Balay do { \ 1508362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_TRUE); \ 15092d53ad75SBarry Smith PetscRegister__FUNCT__(); \ 151053c77d0aSJed Brown } while (0) 151153c77d0aSJed Brown 1512a8d2bbe5SBarry Smith /*MC 1513530d85c1SBarry Smith PetscFunctionBeginUser - First executable line of user provided routines 1514a8d2bbe5SBarry Smith 1515a8d2bbe5SBarry Smith Synopsis: 1516aaa7dc30SBarry Smith #include <petscsys.h> 1517a8d2bbe5SBarry Smith void PetscFunctionBeginUser; 1518a8d2bbe5SBarry Smith 151920f4b53cSBarry Smith Not Collective; No Fortran Support 1520a8d2bbe5SBarry Smith 1521a8d2bbe5SBarry Smith Usage: 1522a8d2bbe5SBarry Smith .vb 1523a8d2bbe5SBarry Smith int something; 1524a8d2bbe5SBarry Smith 1525ac285190SBarry Smith PetscFunctionBeginUser; 1526a8d2bbe5SBarry Smith .ve 1527a8d2bbe5SBarry Smith 152820f4b53cSBarry Smith Level: intermediate 152920f4b53cSBarry Smith 1530a8d2bbe5SBarry Smith Notes: 1531530d85c1SBarry Smith Functions that incorporate this must call `PetscFunctionReturn()` instead of return except for main(). 1532530d85c1SBarry Smith 1533530d85c1SBarry Smith May be used before `PetscInitialize()` 15341957e957SBarry Smith 1535530d85c1SBarry Smith This is identical to `PetscFunctionBegin` except it labels the routine as a user 1536ac285190SBarry Smith routine instead of as a PETSc library routine. 1537ac285190SBarry Smith 1538586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, `PetscFunctionBeginHot`, `PetscStackPushNoCheck()` 1539a8d2bbe5SBarry Smith M*/ 15409371c9d4SSatish Balay #define PetscFunctionBeginUser \ 15419371c9d4SSatish Balay do { \ 1542362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 2, PETSC_FALSE); \ 1543a8d2bbe5SBarry Smith PetscRegister__FUNCT__(); \ 1544a8d2bbe5SBarry Smith } while (0) 1545a8d2bbe5SBarry Smith 1546586f9135SBarry Smith /*MC 1547586f9135SBarry Smith PetscStackPush - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is 1548586f9135SBarry Smith currently in the source code and verifies the memory is not corrupted. 1549586f9135SBarry Smith 1550586f9135SBarry Smith Synopsis: 1551586f9135SBarry Smith #include <petscsys.h> 1552586f9135SBarry Smith void PetscStackPush(char *funct) 1553586f9135SBarry Smith 15547cdbe19fSJose E. Roman Not Collective 15557cdbe19fSJose E. Roman 1556586f9135SBarry Smith Input Parameter: 1557586f9135SBarry Smith . funct - the function name 1558586f9135SBarry Smith 1559586f9135SBarry Smith Level: developer 1560586f9135SBarry Smith 1561586f9135SBarry Smith Notes: 1562586f9135SBarry 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 1563586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1564586f9135SBarry Smith help debug the problem. 1565586f9135SBarry Smith 156695bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 1567586f9135SBarry Smith 1568586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1569586f9135SBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1570586f9135SBarry Smith M*/ 15719371c9d4SSatish Balay #define PetscStackPush(n) \ 15729371c9d4SSatish Balay do { \ 1573362febeeSStefano Zampini PetscStackPushNoCheck(n, 0, PETSC_FALSE); \ 157415681b3cSBarry Smith CHKMEMQ; \ 157515681b3cSBarry Smith } while (0) 15763a40ed3dSBarry Smith 1577586f9135SBarry Smith /*MC 1578586f9135SBarry Smith PetscStackPop - Pops a function name from the PETSc default stack that tracks where the running program is 1579586f9135SBarry Smith currently in the source code and verifies the memory is not corrupted. 1580586f9135SBarry Smith 1581586f9135SBarry Smith Synopsis: 1582586f9135SBarry Smith #include <petscsys.h> 1583586f9135SBarry Smith void PetscStackPop 1584586f9135SBarry Smith 15857cdbe19fSJose E. Roman Not Collective 15867cdbe19fSJose E. Roman 1587586f9135SBarry Smith Level: developer 1588586f9135SBarry Smith 1589586f9135SBarry Smith Notes: 1590586f9135SBarry 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 1591586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1592586f9135SBarry Smith help debug the problem. 1593586f9135SBarry Smith 159495bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 1595586f9135SBarry Smith 1596586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPopNoCheck()`, `PetscStackPush()` 1597586f9135SBarry Smith M*/ 15989371c9d4SSatish Balay #define PetscStackPop \ 15999371c9d4SSatish Balay do { \ 1600441dd030SJed Brown CHKMEMQ; \ 1601362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 160215681b3cSBarry Smith } while (0) 1603d64ed03dSBarry Smith 160430de9b25SBarry Smith /*MC 16050a57284eSJacob Faibussowitsch PetscFunctionReturn - Last executable line of each PETSc function used for error 16060a57284eSJacob Faibussowitsch handling. Replaces `return()`. 160730de9b25SBarry Smith 160830de9b25SBarry Smith Synopsis: 16090a57284eSJacob Faibussowitsch #include <petscerror.h> 16100a57284eSJacob Faibussowitsch void PetscFunctionReturn(...) 161130de9b25SBarry Smith 1612667f096bSBarry Smith Not Collective; No Fortran Support 1613eca87e8dSBarry Smith 16145687f061SJacob Faibussowitsch Level: beginner 16155687f061SJacob Faibussowitsch 16160a57284eSJacob Faibussowitsch Notes: 16170a57284eSJacob Faibussowitsch This routine is a macro, so while it does not "return" anything itself, it does return from 16180a57284eSJacob Faibussowitsch the function in the literal sense. 16190a57284eSJacob Faibussowitsch 16200a57284eSJacob Faibussowitsch Usually the return value is the integer literal `0` (for example in any function returning 16210a57284eSJacob Faibussowitsch `PetscErrorCode`), however it is possible to return any arbitrary type. The arguments of 16220a57284eSJacob Faibussowitsch this macro are placed before the `return` statement as-is. 16230a57284eSJacob Faibussowitsch 16240a57284eSJacob Faibussowitsch Any routine which returns via `PetscFunctionReturn()` must begin with a corresponding 16250a57284eSJacob Faibussowitsch `PetscFunctionBegin`. 16260a57284eSJacob Faibussowitsch 16270a57284eSJacob Faibussowitsch For routines which return `void` use `PetscFunctionReturnVoid()` instead. 16280a57284eSJacob Faibussowitsch 16290a57284eSJacob Faibussowitsch Example Usage: 163030de9b25SBarry Smith .vb 16310a57284eSJacob Faibussowitsch PetscErrorCode foo(int *x) 16320a57284eSJacob Faibussowitsch { 16330a57284eSJacob Faibussowitsch PetscFunctionBegin; // don't forget the begin! 16340a57284eSJacob Faibussowitsch *x = 10; 16353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 163630de9b25SBarry Smith } 163730de9b25SBarry Smith .ve 163830de9b25SBarry Smith 16390a57284eSJacob Faibussowitsch May return any arbitrary type\: 16400a57284eSJacob Faibussowitsch .vb 16410a57284eSJacob Faibussowitsch struct Foo 16420a57284eSJacob Faibussowitsch { 16430a57284eSJacob Faibussowitsch int x; 16440a57284eSJacob Faibussowitsch }; 16450a57284eSJacob Faibussowitsch 16460a57284eSJacob Faibussowitsch struct Foo make_foo(int value) 16470a57284eSJacob Faibussowitsch { 16480a57284eSJacob Faibussowitsch struct Foo f; 16490a57284eSJacob Faibussowitsch 16500a57284eSJacob Faibussowitsch PetscFunctionBegin; 16510a57284eSJacob Faibussowitsch f.x = value; 16520a57284eSJacob Faibussowitsch PetscFunctionReturn(f); 16530a57284eSJacob Faibussowitsch } 16540a57284eSJacob Faibussowitsch .ve 16550a57284eSJacob Faibussowitsch 16560a57284eSJacob Faibussowitsch .seealso: `PetscFunctionBegin`, `PetscFunctionBeginUser`, `PetscFunctionReturnVoid()`, 16570a57284eSJacob Faibussowitsch `PetscStackPopNoCheck()` 165830de9b25SBarry Smith M*/ 16595687f061SJacob Faibussowitsch #define PetscFunctionReturn(...) \ 16609371c9d4SSatish Balay do { \ 1661362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 16625687f061SJacob Faibussowitsch return __VA_ARGS__; \ 166327104ee2SJacob Faibussowitsch } while (0) 1664d64ed03dSBarry Smith 16650a57284eSJacob Faibussowitsch /*MC 16660a57284eSJacob Faibussowitsch PetscFunctionReturnVoid - Like `PetscFunctionReturn()` but returns `void` 16670a57284eSJacob Faibussowitsch 16680a57284eSJacob Faibussowitsch Synopsis: 16690a57284eSJacob Faibussowitsch #include <petscerror.h> 16700a57284eSJacob Faibussowitsch void PetscFunctionReturnVoid() 16710a57284eSJacob Faibussowitsch 16720a57284eSJacob Faibussowitsch Not Collective 16730a57284eSJacob Faibussowitsch 16745687f061SJacob Faibussowitsch Level: beginner 16755687f061SJacob Faibussowitsch 16765687f061SJacob Faibussowitsch Note: 16770a57284eSJacob Faibussowitsch Behaves identically to `PetscFunctionReturn()` except that it returns `void`. That is, this 16785687f061SJacob Faibussowitsch macro culminates with `return`. 16790a57284eSJacob Faibussowitsch 16800a57284eSJacob Faibussowitsch Example Usage: 16810a57284eSJacob Faibussowitsch .vb 16820a57284eSJacob Faibussowitsch void foo() 16830a57284eSJacob Faibussowitsch { 16840a57284eSJacob Faibussowitsch PetscFunctionBegin; // must start with PetscFunctionBegin! 16850a57284eSJacob Faibussowitsch bar(); 16860a57284eSJacob Faibussowitsch baz(); 16870a57284eSJacob Faibussowitsch PetscFunctionReturnVoid(); 16880a57284eSJacob Faibussowitsch } 16890a57284eSJacob Faibussowitsch .ve 16900a57284eSJacob Faibussowitsch 16910a57284eSJacob Faibussowitsch .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, PetscFunctionBeginUser` 16920a57284eSJacob Faibussowitsch M*/ 16939371c9d4SSatish Balay #define PetscFunctionReturnVoid() \ 16949371c9d4SSatish Balay do { \ 1695362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 169627104ee2SJacob Faibussowitsch return; \ 169727104ee2SJacob Faibussowitsch } while (0) 169827104ee2SJacob Faibussowitsch #else /* PETSC_USE_DEBUG */ 169927104ee2SJacob Faibussowitsch #define PetscStackPushNoCheck(funct, petsc_routine, hot) 170037154d10SBarry Smith #define PetscStackUpdateLine 1701792fecdfSBarry Smith #define PetscStackPushExternal(funct) 17023ba16761SJacob Faibussowitsch #define PetscStackPopNoCheck(...) 170327104ee2SJacob Faibussowitsch #define PetscStackClearTop 17043a40ed3dSBarry Smith #define PetscFunctionBegin 17050bdf7c52SPeter Brune #define PetscFunctionBeginUser 1706a2f94806SJed Brown #define PetscFunctionBeginHot 17070a57284eSJacob Faibussowitsch #define PetscFunctionReturn(...) return __VA_ARGS__ 17085665465eSBarry Smith #define PetscFunctionReturnVoid() return 1709812af9f3SBarry Smith #define PetscStackPop CHKMEMQ 1710812af9f3SBarry Smith #define PetscStackPush(f) CHKMEMQ 171127104ee2SJacob Faibussowitsch #endif /* PETSC_USE_DEBUG */ 17123a40ed3dSBarry Smith 17136d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 17143ba16761SJacob Faibussowitsch #define PetscStackCallExternalVoid(...) 17153ba16761SJacob Faibussowitsch template <typename F, typename... Args> 17163ba16761SJacob Faibussowitsch void PetscCallExternal(F, Args...); 17176d210af2SJacob Faibussowitsch #else 1718586f9135SBarry Smith /*MC 1719e77caa6dSBarry Smith PetscStackCallExternalVoid - Calls an external library routine or user function after pushing the name of the routine on the stack. 1720eb6b5d47SBarry Smith 1721eb6b5d47SBarry Smith Input Parameters: 1722eb6b5d47SBarry Smith + name - string that gives the name of the function being called 1723586f9135SBarry Smith - routine - actual call to the routine, for example, functionname(a,b) 1724fd3f9acdSBarry Smith 1725586f9135SBarry Smith Level: developer 1726eb6b5d47SBarry Smith 172795bd0b28SBarry Smith Notes: 1728792fecdfSBarry Smith Often one should use `PetscCallExternal()` instead. This routine is intended for external library routines that DO NOT return error codes 1729eb6b5d47SBarry Smith 1730586f9135SBarry Smith In debug mode this also checks the memory for corruption at the end of the function call. 1731586f9135SBarry Smith 173295bd0b28SBarry Smith Certain external packages, such as BLAS/LAPACK may have their own macros, `PetscCallBLAS()` for managing the call, error checking, etc. 1733586f9135SBarry Smith 1734586f9135SBarry Smith Developer Note: 1735586f9135SBarry 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. 1736586f9135SBarry Smith 1737792fecdfSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscCallExternal()`, `PetscCallBLAS()` 1738586f9135SBarry Smith @*/ 17393ba16761SJacob Faibussowitsch #define PetscStackCallExternalVoid(name, ...) \ 17409371c9d4SSatish Balay do { \ 174128770333SStefano Zampini PetscStackPushExternal(name); \ 17423ba16761SJacob Faibussowitsch __VA_ARGS__; \ 17439371c9d4SSatish Balay PetscStackPop; \ 17449371c9d4SSatish Balay } while (0) 1745eb6b5d47SBarry Smith 1746586f9135SBarry Smith /*MC 1747792fecdfSBarry Smith PetscCallExternal - Calls an external library routine that returns an error code after pushing the name of the routine on the stack. 1748fd3f9acdSBarry Smith 1749fd3f9acdSBarry Smith Input Parameters: 1750fd3f9acdSBarry Smith + func - name of the routine 1751586f9135SBarry Smith - args - arguments to the routine 1752586f9135SBarry Smith 1753586f9135SBarry Smith Level: developer 1754fd3f9acdSBarry Smith 175595452b02SPatrick Sanan Notes: 1756e77caa6dSBarry Smith This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not. 1757dbf62e16SBarry Smith 1758586f9135SBarry Smith In debug mode this also checks the memory for corruption at the end of the function call. 1759fd3f9acdSBarry Smith 176087497f52SBarry Smith Assumes the error return code of the function is an integer and that a value of 0 indicates success 176187497f52SBarry Smith 1762586f9135SBarry Smith Developer Note: 1763d5b43468SJose E. Roman This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc. 1764586f9135SBarry Smith 1765e77caa6dSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()` 1766586f9135SBarry Smith M*/ 17679371c9d4SSatish Balay #define PetscCallExternal(func, ...) \ 17689371c9d4SSatish Balay do { \ 1769a74df02fSJacob Faibussowitsch PetscStackPush(PetscStringize(func)); \ 17703ba16761SJacob Faibussowitsch int ierr_petsc_call_external_ = func(__VA_ARGS__); \ 17711d4906efSStefano Zampini PetscStackPop; \ 17723ba16761SJacob 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_); \ 1773fd3f9acdSBarry Smith } while (0) 17746d210af2SJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */ 1775