154a8ef01SBarry Smith /* 2f621e05eSBarry Smith Contains all error handling interfaces for PETSc. 354a8ef01SBarry Smith */ 4a4963045SJacob Faibussowitsch #pragma once 56c7e564aSBarry Smith 6e1bf4ed2SJacob Faibussowitsch #include <petscmacros.h> 7e1bf4ed2SJacob Faibussowitsch #include <petscsystypes.h> 8e1bf4ed2SJacob Faibussowitsch 9af75f258SJacob Faibussowitsch #if defined(__cplusplus) 10af75f258SJacob Faibussowitsch #include <exception> // std::exception 11af75f258SJacob Faibussowitsch #endif 12af75f258SJacob Faibussowitsch 13ac09b921SBarry Smith /* SUBMANSEC = Sys */ 14ac09b921SBarry Smith 15edd03b47SJacob Faibussowitsch #define SETERRQ1(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 16edd03b47SJacob Faibussowitsch #define SETERRQ2(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 17edd03b47SJacob Faibussowitsch #define SETERRQ3(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 18edd03b47SJacob Faibussowitsch #define SETERRQ4(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 19edd03b47SJacob Faibussowitsch #define SETERRQ5(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 20edd03b47SJacob Faibussowitsch #define SETERRQ6(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 21edd03b47SJacob Faibussowitsch #define SETERRQ7(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 22edd03b47SJacob Faibussowitsch #define SETERRQ8(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 23edd03b47SJacob Faibussowitsch #define SETERRQ9(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__) 2498921bdaSJacob Faibussowitsch 2530de9b25SBarry Smith /*MC 261957e957SBarry Smith SETERRQ - Macro to be called when an error has been detected, 2730de9b25SBarry Smith 2830de9b25SBarry Smith Synopsis: 29aaa7dc30SBarry Smith #include <petscsys.h> 3098921bdaSJacob Faibussowitsch PetscErrorCode SETERRQ(MPI_Comm comm, PetscErrorCode ierr, char *message, ...) 3130de9b25SBarry Smith 32d083f849SBarry Smith Collective 3330de9b25SBarry Smith 3430de9b25SBarry Smith Input Parameters: 3595bd0b28SBarry Smith + comm - An MPI communicator, use `PETSC_COMM_SELF` unless you know all ranks of another communicator will detect the error 363af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 3730de9b25SBarry Smith - message - error message 3830de9b25SBarry Smith 3930de9b25SBarry Smith Level: beginner 4030de9b25SBarry Smith 4130de9b25SBarry Smith Notes: 4287497f52SBarry Smith This is rarely needed, one should use `PetscCheck()` and `PetscCall()` and friends to automatically handle error conditions. 4330de9b25SBarry Smith Once the error handler is called the calling function is then returned from with the given error code. 4430de9b25SBarry Smith 4587497f52SBarry Smith Experienced users can set the error handler with `PetscPushErrorHandler()`. 4630de9b25SBarry Smith 4716a05f60SBarry Smith Fortran Note: 48989c49a2SBarry Smith `SETERRQ()` may be called from Fortran subroutines but `SETERRA()` must be called from the 493b1008b8SBarry Smith Fortran main program. 503b1008b8SBarry Smith 51db781477SPatrick Sanan .seealso: `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 52af27ebaaSBarry Smith `PetscError()`, `PetscCall()`, `CHKMEMQ`, `CHKERRA()`, `PetscCallMPI()`, `PetscErrorCode` 5330de9b25SBarry Smith M*/ 54e809461dSJacob Faibussowitsch #define SETERRQ(comm, ierr, ...) \ 55e809461dSJacob Faibussowitsch do { \ 56e809461dSJacob Faibussowitsch PetscErrorCode ierr_seterrq_petsc_ = PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \ 57e809461dSJacob Faibussowitsch return ierr_seterrq_petsc_ ? ierr_seterrq_petsc_ : PETSC_ERR_RETURN; \ 58e809461dSJacob Faibussowitsch } while (0) 59986eef2eSBarry Smith 60648c30bcSBarry Smith #define SETERRQNULL(comm, ierr, ...) \ 61648c30bcSBarry Smith do { \ 62648c30bcSBarry Smith (void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \ 63648c30bcSBarry Smith return NULL; \ 64648c30bcSBarry Smith } while (0) 65648c30bcSBarry Smith 66ffc4695bSBarry Smith /* 67ffc4695bSBarry Smith Returned from PETSc functions that are called from MPI, such as related to attributes 68ffc4695bSBarry Smith Do not confuse PETSC_MPI_ERROR_CODE and PETSC_ERR_MPI, the first is registered with MPI and returned to MPI as 69888a1fb3SBarry Smith an error code, the latter is a regular PETSc error code passed within PETSc code indicating an error was detected in an MPI call. 70ffc4695bSBarry Smith */ 71ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CLASS; 72ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CODE; 73ffc4695bSBarry Smith 74986eef2eSBarry Smith /*MC 75986eef2eSBarry Smith SETERRMPI - Macro to be called when an error has been detected within an MPI callback function 76986eef2eSBarry Smith 77989c49a2SBarry Smith No Fortran Support 78989c49a2SBarry Smith 79986eef2eSBarry Smith Synopsis: 80986eef2eSBarry Smith #include <petscsys.h> 8198921bdaSJacob Faibussowitsch PetscErrorCode SETERRMPI(MPI_Comm comm, PetscErrorCode ierr, char *message, ...) 82986eef2eSBarry Smith 83d083f849SBarry Smith Collective 84986eef2eSBarry Smith 85986eef2eSBarry Smith Input Parameters: 8695bd0b28SBarry Smith + comm - An MPI communicator, use `PETSC_COMM_SELF` unless you know all ranks of another communicator will detect the error 87986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 88986eef2eSBarry Smith - message - error message 89986eef2eSBarry Smith 90986eef2eSBarry Smith Level: developer 91986eef2eSBarry Smith 9295bd0b28SBarry Smith Note: 9387497f52SBarry Smith This macro is FOR USE IN MPI CALLBACK FUNCTIONS ONLY, such as those passed to `MPI_Comm_create_keyval()`. It always returns the error code `PETSC_MPI_ERROR_CODE` 9487497f52SBarry Smith which is registered with `MPI_Add_error_code()` when PETSc is initialized. 95986eef2eSBarry Smith 96af27ebaaSBarry Smith .seealso: `SETERRQ()`, `PetscCall()`, `PetscCallMPI()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `PetscErrorCode` 97986eef2eSBarry Smith M*/ 983ba16761SJacob Faibussowitsch #define SETERRMPI(comm, ierr, ...) return ((void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__), PETSC_MPI_ERROR_CODE) 99ee8199e6SMatthew G. Knepley 100ee8199e6SMatthew G. Knepley /*MC 101f388eb8bSPatrick Sanan SETERRA - Fortran-only macro that can be called when an error has been detected from the main program 102f388eb8bSPatrick Sanan 103f388eb8bSPatrick Sanan Synopsis: 104f388eb8bSPatrick Sanan #include <petscsys.h> 105f388eb8bSPatrick Sanan PetscErrorCode SETERRA(MPI_Comm comm, PetscErrorCode ierr, char *message) 106f388eb8bSPatrick Sanan 107f388eb8bSPatrick Sanan Collective 108f388eb8bSPatrick Sanan 109f388eb8bSPatrick Sanan Input Parameters: 11095bd0b28SBarry Smith + comm - An MPI communicator, so that the error can be collective 111f388eb8bSPatrick Sanan . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 112e55b9c91SBarry Smith - message - error message in the `printf()` format 113f388eb8bSPatrick Sanan 114f388eb8bSPatrick Sanan Level: beginner 115f388eb8bSPatrick Sanan 116f388eb8bSPatrick Sanan Notes: 11787497f52SBarry Smith This should only be used with Fortran. With C/C++, use `SETERRQ()`. 118f388eb8bSPatrick Sanan 11987497f52SBarry Smith `SETERRQ()` may be called from Fortran subroutines but `SETERRA()` must be called from the 120f388eb8bSPatrick Sanan Fortran main program. 121f388eb8bSPatrick Sanan 122af27ebaaSBarry Smith .seealso: `SETERRQ()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()`, `PetscErrorCode` 123f388eb8bSPatrick Sanan M*/ 124f388eb8bSPatrick Sanan 125f388eb8bSPatrick Sanan /*MC 126fa190f98SMatthew G. Knepley SETERRABORT - Macro that can be called when an error has been detected, 127fa190f98SMatthew G. Knepley 128fa190f98SMatthew G. Knepley Synopsis: 129fa190f98SMatthew G. Knepley #include <petscsys.h> 13098921bdaSJacob Faibussowitsch PetscErrorCode SETERRABORT(MPI_Comm comm, PetscErrorCode ierr, char *message, ...) 131fa190f98SMatthew G. Knepley 132d083f849SBarry Smith Collective 133fa190f98SMatthew G. Knepley 134fa190f98SMatthew G. Knepley Input Parameters: 13595bd0b28SBarry Smith + comm - An MPI communicator, so that the error can be collective 1363af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 137e55b9c91SBarry Smith - message - error message in the `printf()` format 138fa190f98SMatthew G. Knepley 139fa190f98SMatthew G. Knepley Level: beginner 140fa190f98SMatthew G. Knepley 141fa190f98SMatthew G. Knepley Notes: 14287497f52SBarry Smith This function just calls `MPI_Abort()`. 14387497f52SBarry Smith 14487497f52SBarry Smith This should only be called in routines that cannot return an error code, such as in C++ constructors. 145fa190f98SMatthew G. Knepley 146989c49a2SBarry Smith Fortran Note: 147989c49a2SBarry Smith Use `SETERRA()` in Fortran main program and `SETERRQ()` in Fortran subroutines 148989c49a2SBarry Smith 149989c49a2SBarry Smith Developer Note: 150989c49a2SBarry Smith In Fortran `SETERRA()` could be called `SETERRABORT()` since they serve the same purpose 151989c49a2SBarry Smith 152af27ebaaSBarry Smith .seealso: `SETERRQ()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `PetscCall()`, `CHKMEMQ`, `PetscErrorCode` 153fa190f98SMatthew G. Knepley M*/ 1549371c9d4SSatish Balay #define SETERRABORT(comm, ierr, ...) \ 1559371c9d4SSatish Balay do { \ 1563ba16761SJacob Faibussowitsch (void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \ 157dd460d27SBarry Smith (void)MPI_Abort(comm, ierr); \ 15898921bdaSJacob Faibussowitsch } while (0) 1599a00fa46SSatish Balay 16030de9b25SBarry Smith /*MC 161139c8099SSatish Balay PetscCheck - Checks that a particular condition is true; if not true, then returns the provided error code 1622c71b3e2SJacob Faibussowitsch 1632c71b3e2SJacob Faibussowitsch Synopsis: 1642c71b3e2SJacob Faibussowitsch #include <petscerror.h> 1652c71b3e2SJacob Faibussowitsch void PetscCheck(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 1662c71b3e2SJacob Faibussowitsch 1677cdbe19fSJose E. Roman Collective; No Fortran Support 1682c71b3e2SJacob Faibussowitsch 1692c71b3e2SJacob Faibussowitsch Input Parameters: 1702c71b3e2SJacob Faibussowitsch + cond - The boolean condition 1712c71b3e2SJacob Faibussowitsch . comm - The communicator on which the check can be collective on 1722c71b3e2SJacob Faibussowitsch . ierr - A nonzero error code, see include/petscerror.h for the complete list 173e55b9c91SBarry Smith - message - Error message in the `printf()` format 1742c71b3e2SJacob Faibussowitsch 175667f096bSBarry Smith Level: beginner 176667f096bSBarry Smith 1772c71b3e2SJacob Faibussowitsch Notes: 1782c71b3e2SJacob Faibussowitsch Enabled in both optimized and debug builds. 1792c71b3e2SJacob Faibussowitsch 180a2454668SBarry Smith As a general rule, `PetscCheck()` is used to check "usage error" (for example, passing an incorrect value as a function argument), 181a2454668SBarry Smith `PetscAssert()` is used to "check for bugs in PETSc" (for example, is a value in a PETSc data structure nonsensical). 182a2454668SBarry Smith However, for functions that are called in a "hot spot", for example, thousands of times in a loop, `PetscAssert()` should be used instead 183a2454668SBarry Smith of `PetscCheck()` since the former is compiled out in PETSc's optimization code. 184a2454668SBarry Smith 18587497f52SBarry Smith Calls `SETERRQ()` if the assertion fails, so can only be called from functions returning a 18687497f52SBarry Smith `PetscErrorCode` (or equivalent type after conversion). 1872c71b3e2SJacob Faibussowitsch 1887c5b2466SBarry Smith .seealso: `PetscAssert()`, `PetscCheckReturnMPI()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheckAbort()`, `PetscErrorCode` 1899566063dSJacob Faibussowitsch M*/ 1909371c9d4SSatish Balay #define PetscCheck(cond, comm, ierr, ...) \ 1919371c9d4SSatish Balay do { \ 1929371c9d4SSatish Balay if (PetscUnlikely(!(cond))) SETERRQ(comm, ierr, __VA_ARGS__); \ 1939371c9d4SSatish Balay } while (0) 1942c71b3e2SJacob Faibussowitsch 1952c71b3e2SJacob Faibussowitsch /*MC 1967c5b2466SBarry Smith PetscCheckReturnMPI - Checks that a particular condition is true; if not true, then returns an MPI error code. 1977c5b2466SBarry Smith To check for errors in PETSc-provided MPI callbacks. 1987c5b2466SBarry Smith 1997c5b2466SBarry Smith Synopsis: 2007c5b2466SBarry Smith #include <petscerror.h> 2017c5b2466SBarry Smith void PetscCheckReturnMPI(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 2027c5b2466SBarry Smith 2037c5b2466SBarry Smith Collective; No Fortran Support 2047c5b2466SBarry Smith 2057c5b2466SBarry Smith Input Parameters: 2067c5b2466SBarry Smith + cond - The boolean condition 2077c5b2466SBarry Smith . comm - The communicator on which the check can be collective on 2087c5b2466SBarry Smith . ierr - A nonzero error code, see include/petscerror.h for the complete list 2097c5b2466SBarry Smith - message - Error message in the `printf()` format 2107c5b2466SBarry Smith 2117c5b2466SBarry Smith Level: beginner 2127c5b2466SBarry Smith 2137c5b2466SBarry Smith Note: 2147c5b2466SBarry Smith Enabled in both optimized and debug builds. 2157c5b2466SBarry Smith 2167c5b2466SBarry Smith .seealso: `PetscCheck()`, `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheckAbort()`, `PetscErrorCode` 2177c5b2466SBarry Smith M*/ 2187c5b2466SBarry Smith #define PetscCheckReturnMPI(cond, comm, ierr, ...) \ 2197c5b2466SBarry Smith do { \ 2207c5b2466SBarry Smith if (PetscUnlikely(!(cond))) SETERRMPI(comm, ierr, __VA_ARGS__); \ 2217c5b2466SBarry Smith } while (0) 2227c5b2466SBarry Smith 2237c5b2466SBarry Smith /*MC 2244be741a6SBarry Smith PetscCheckAbort - Check that a particular condition is true, otherwise prints error and aborts 2254be741a6SBarry Smith 2264be741a6SBarry Smith Synopsis: 2274be741a6SBarry Smith #include <petscerror.h> 2284be741a6SBarry Smith void PetscCheckAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 2294be741a6SBarry Smith 2307cdbe19fSJose E. Roman Collective; No Fortran Support 2314be741a6SBarry Smith 2324be741a6SBarry Smith Input Parameters: 2334be741a6SBarry Smith + cond - The boolean condition 2344be741a6SBarry Smith . comm - The communicator on which the check can be collective on 2354be741a6SBarry Smith . ierr - A nonzero error code, see include/petscerror.h for the complete list 236e55b9c91SBarry Smith - message - Error message in the `printf()` format 2374be741a6SBarry Smith 238667f096bSBarry Smith Level: developer 239667f096bSBarry Smith 2404be741a6SBarry Smith Notes: 2414be741a6SBarry Smith Enabled in both optimized and debug builds. 2424be741a6SBarry Smith 24387497f52SBarry Smith Calls `SETERRABORT()` if the assertion fails, can be called from a function that does not return an 24487497f52SBarry Smith error code, such as a C++ constructor. usually `PetscCheck()` should be used. 2454be741a6SBarry Smith 246af27ebaaSBarry Smith .seealso: `PetscAssertAbort()`, `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheck()`, `SETERRABORT()`, `PetscErrorCode` 2474be741a6SBarry Smith M*/ 2489371c9d4SSatish Balay #define PetscCheckAbort(cond, comm, ierr, ...) \ 2490e6b6b59SJacob Faibussowitsch do { \ 2500e6b6b59SJacob Faibussowitsch if (PetscUnlikely(!(cond))) SETERRABORT(comm, ierr, __VA_ARGS__); \ 251c7481402SJacob Faibussowitsch } while (0) 2524be741a6SBarry Smith 2534be741a6SBarry Smith /*MC 2549ace16cdSJacob Faibussowitsch PetscAssert - Assert that a particular condition is true 2559ace16cdSJacob Faibussowitsch 2569ace16cdSJacob Faibussowitsch Synopsis: 2579ace16cdSJacob Faibussowitsch #include <petscerror.h> 2589ace16cdSJacob Faibussowitsch void PetscAssert(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 2599ace16cdSJacob Faibussowitsch 2607cdbe19fSJose E. Roman Collective; No Fortran Support 2619ace16cdSJacob Faibussowitsch 2629ace16cdSJacob Faibussowitsch Input Parameters: 2639ace16cdSJacob Faibussowitsch + cond - The boolean condition 2649ace16cdSJacob Faibussowitsch . comm - The communicator on which the check can be collective on 2659ace16cdSJacob Faibussowitsch . ierr - A nonzero error code, see include/petscerror.h for the complete list 266af27ebaaSBarry Smith - message - Error message in `printf()` format 2679ace16cdSJacob Faibussowitsch 268667f096bSBarry Smith Level: beginner 269667f096bSBarry Smith 2709ace16cdSJacob Faibussowitsch Notes: 271c7481402SJacob Faibussowitsch Equivalent to `PetscCheck()` if debugging is enabled, and `PetscAssume(cond)` otherwise. 2722c71b3e2SJacob Faibussowitsch 27387497f52SBarry Smith See `PetscCheck()` for usage and behaviour. 27487497f52SBarry Smith 27587497f52SBarry Smith This is needed instead of simply using `assert()` because this correctly handles the collective nature of errors under MPI 2769ace16cdSJacob Faibussowitsch 277af27ebaaSBarry Smith .seealso: `PetscCheck()`, `SETERRQ()`, `PetscError()`, `PetscAssertAbort()`, `PetscErrorCode` 2789566063dSJacob Faibussowitsch M*/ 279c7481402SJacob Faibussowitsch #if PetscDefined(USE_DEBUG) 280c7481402SJacob Faibussowitsch #define PetscAssert(cond, comm, ierr, ...) PetscCheck(cond, comm, ierr, __VA_ARGS__) 281c7481402SJacob Faibussowitsch #else 282c7481402SJacob Faibussowitsch #define PetscAssert(cond, ...) PetscAssume(cond) 283c7481402SJacob Faibussowitsch #endif 2849ace16cdSJacob Faibussowitsch 2859ace16cdSJacob Faibussowitsch /*MC 2860e6b6b59SJacob Faibussowitsch PetscAssertAbort - Assert that a particular condition is true, otherwise prints error and aborts 2870e6b6b59SJacob Faibussowitsch 2880e6b6b59SJacob Faibussowitsch Synopsis: 2890e6b6b59SJacob Faibussowitsch #include <petscerror.h> 2900e6b6b59SJacob Faibussowitsch void PetscAssertAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...) 2910e6b6b59SJacob Faibussowitsch 2927cdbe19fSJose E. Roman Collective; No Fortran Support 2930e6b6b59SJacob Faibussowitsch 2940e6b6b59SJacob Faibussowitsch Input Parameters: 2950e6b6b59SJacob Faibussowitsch + cond - The boolean condition 2960e6b6b59SJacob Faibussowitsch . comm - The communicator on which the check can be collective on 2970e6b6b59SJacob Faibussowitsch . ierr - A nonzero error code, see include/petscerror.h for the complete list 298e55b9c91SBarry Smith - message - Error message in the `printf()` format 2990e6b6b59SJacob Faibussowitsch 300667f096bSBarry Smith Level: beginner 301667f096bSBarry Smith 30295bd0b28SBarry Smith Note: 3030e6b6b59SJacob Faibussowitsch Enabled only in debug builds. See `PetscCheckAbort()` for usage. 3040e6b6b59SJacob Faibussowitsch 3050e6b6b59SJacob Faibussowitsch .seealso: `PetscCheckAbort()`, `PetscAssert()`, `PetscCheck()`, `SETERRABORT()`, `PetscError()` 3060e6b6b59SJacob Faibussowitsch M*/ 3073ba16761SJacob Faibussowitsch #if PetscDefined(USE_DEBUG) 3083ba16761SJacob Faibussowitsch #define PetscAssertAbort(cond, comm, ierr, ...) PetscCheckAbort(cond, comm, ierr, __VA_ARGS__) 3093ba16761SJacob Faibussowitsch #else 3103ba16761SJacob Faibussowitsch #define PetscAssertAbort(cond, comm, ierr, ...) PetscAssume(cond) 3113ba16761SJacob Faibussowitsch #endif 3120e6b6b59SJacob Faibussowitsch 3130e6b6b59SJacob Faibussowitsch /*MC 3143ba16761SJacob Faibussowitsch PetscCall - Calls a PETSc function and then checks the resulting error code, if it is 3153ba16761SJacob Faibussowitsch non-zero it calls the error handler and returns from the current function with the error 3163ba16761SJacob Faibussowitsch code. 3179566063dSJacob Faibussowitsch 3189566063dSJacob Faibussowitsch Synopsis: 3199566063dSJacob Faibussowitsch #include <petscerror.h> 32049c86fc7SBarry Smith void PetscCall(PetscFunction(args)) 3219566063dSJacob Faibussowitsch 3229566063dSJacob Faibussowitsch Not Collective 3239566063dSJacob Faibussowitsch 3249566063dSJacob Faibussowitsch Input Parameter: 32549c86fc7SBarry Smith . PetscFunction - any PETSc function that returns an error code 3269566063dSJacob Faibussowitsch 327667f096bSBarry Smith Level: beginner 328667f096bSBarry Smith 3299566063dSJacob Faibussowitsch Notes: 3309566063dSJacob Faibussowitsch Once the error handler is called the calling function is then returned from with the given 33187497f52SBarry Smith error code. Experienced users can set the error handler with `PetscPushErrorHandler()`. 3329566063dSJacob Faibussowitsch 33387497f52SBarry Smith `PetscCall()` cannot be used in functions returning a datatype not convertible to 3348af9f190SBarry Smith `PetscErrorCode`. For example, `PetscCall()` may not be used in functions returning `void`, use 3353ba16761SJacob Faibussowitsch `PetscCallAbort()` or `PetscCallVoid()` in this case. 3369566063dSJacob Faibussowitsch 3379566063dSJacob Faibussowitsch Example Usage: 3389566063dSJacob Faibussowitsch .vb 3399566063dSJacob Faibussowitsch PetscCall(PetscInitiailize(...)); // OK to call even when PETSc is not yet initialized! 3409566063dSJacob Faibussowitsch 3419566063dSJacob Faibussowitsch struct my_struct 3429566063dSJacob Faibussowitsch { 3439566063dSJacob Faibussowitsch void *data; 3449566063dSJacob Faibussowitsch } my_complex_type; 3459566063dSJacob Faibussowitsch 3469566063dSJacob Faibussowitsch struct my_struct bar(void) 3479566063dSJacob Faibussowitsch { 3486aad120cSJose E. Roman PetscCall(foo(15)); // ERROR PetscErrorCode not convertible to struct my_struct! 3499566063dSJacob Faibussowitsch } 3509566063dSJacob Faibussowitsch 3516aad120cSJose E. Roman PetscCall(bar()) // ERROR input not convertible to PetscErrorCode 3529566063dSJacob Faibussowitsch .ve 3539566063dSJacob Faibussowitsch 35487497f52SBarry Smith It is also possible to call this directly on a `PetscErrorCode` variable 35549c86fc7SBarry Smith .vb 35649c86fc7SBarry Smith PetscCall(ierr); // check if ierr is nonzero 35749c86fc7SBarry Smith .ve 35849c86fc7SBarry Smith 359792fecdfSBarry Smith Should not be used to call callback functions provided by users, `PetscCallBack()` should be used in that situation. 360ef1023bdSBarry Smith 3616a8be23eSBarry Smith `PetscUseTypeMethod()` or `PetscTryTypeMethod()` should be used when calling functions pointers contained in a PETSc object's `ops` array 3626a8be23eSBarry Smith 36349c86fc7SBarry Smith Fortran Notes: 36498d50953SPierre Jolivet The Fortran function in which this is used must declare a `PetscErrorCode` variable necessarily named `ierr`, and `ierr` must be 36587497f52SBarry Smith the final argument to the PETSc function being called. 36649c86fc7SBarry Smith 36798d50953SPierre Jolivet In the main program and in Fortran subroutines that do not have `ierr` as the final return parameter, one 36887497f52SBarry Smith should use `PetscCallA()` 36949c86fc7SBarry Smith 37049c86fc7SBarry Smith Example Fortran Usage: 37149c86fc7SBarry Smith .vb 37249c86fc7SBarry Smith PetscErrorCode ierr 37349c86fc7SBarry Smith Vec v 37449c86fc7SBarry Smith 37549c86fc7SBarry Smith ... 37649c86fc7SBarry Smith PetscCall(VecShift(v, 1.0, ierr)) 37749c86fc7SBarry Smith PetscCallA(VecShift(v, 1.0, ierr)) 37849c86fc7SBarry Smith .ve 37949c86fc7SBarry Smith 380667f096bSBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`, 381667f096bSBarry Smith `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, 382648c30bcSBarry Smith `CHKERRMPI()`, `PetscCallBack()`, `PetscCallAbort()`, `PetscCallVoid()`, `PetscCallNull()` 383648c30bcSBarry Smith M*/ 384648c30bcSBarry Smith 385648c30bcSBarry Smith /*MC 386648c30bcSBarry Smith PetscCallNull - Calls a PETSc function and then checks the resulting error code, if it is 387648c30bcSBarry Smith non-zero it calls the error handler and returns a `NULL` 388648c30bcSBarry Smith 389648c30bcSBarry Smith Synopsis: 390648c30bcSBarry Smith #include <petscerror.h> 391648c30bcSBarry Smith void PetscCallNull(PetscFunction(args)) 392648c30bcSBarry Smith 393648c30bcSBarry Smith Not Collective; No Fortran Support 394648c30bcSBarry Smith 395648c30bcSBarry Smith Input Parameter: 396648c30bcSBarry Smith . PetscFunction - any PETSc function that returns something that can be returned as a `NULL` 397648c30bcSBarry Smith 398648c30bcSBarry Smith Level: developer 399648c30bcSBarry Smith 400648c30bcSBarry Smith .seealso: `PetscCall()`, `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`, 401648c30bcSBarry Smith `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, 402648c30bcSBarry Smith `CHKERRMPI()`, `PetscCallBack()`, `PetscCallAbort()`, `PetscCallVoid()`, `PetscCall()` 4039566063dSJacob Faibussowitsch M*/ 404ef1023bdSBarry Smith 405ef1023bdSBarry Smith /*MC 40698d50953SPierre Jolivet PetscCallA - Fortran-only macro that should be used in the main program and subroutines that do not have `ierr` as the final return parameter, to call PETSc functions instead of using 40798d50953SPierre Jolivet `PetscCall()` which should be used in other Fortran subroutines 408989c49a2SBarry Smith 409989c49a2SBarry Smith Synopsis: 410989c49a2SBarry Smith #include <petscsys.h> 411989c49a2SBarry Smith PetscErrorCode PetscCallA(PetscFunction(arguments, ierr)) 412989c49a2SBarry Smith 413989c49a2SBarry Smith Collective 414989c49a2SBarry Smith 415989c49a2SBarry Smith Input Parameter: 416989c49a2SBarry Smith . PetscFunction(arguments,ierr) - the call to the function 417989c49a2SBarry Smith 418989c49a2SBarry Smith Level: beginner 419989c49a2SBarry Smith 420989c49a2SBarry Smith Notes: 421989c49a2SBarry Smith This should only be used with Fortran. With C/C++, use `PetscCall()` always. 422989c49a2SBarry Smith 42398d50953SPierre Jolivet The Fortran function in which this is used must declare a `PetscErrorCode` variable necessarily named `ierr` 424989c49a2SBarry Smith Use `SETERRA()` to set an error in a Fortran main program and `SETERRQ()` in Fortran subroutines 425989c49a2SBarry Smith 426989c49a2SBarry Smith .seealso: `SETERRQ()`, `SETERRA()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()` 427989c49a2SBarry Smith M*/ 428989c49a2SBarry Smith 429989c49a2SBarry Smith /*MC 430792fecdfSBarry 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 431ef1023bdSBarry Smith handler and returns from the current function with the error code. 432ef1023bdSBarry Smith 433ef1023bdSBarry Smith Synopsis: 434ef1023bdSBarry Smith #include <petscerror.h> 435792fecdfSBarry Smith void PetscCallBack(const char *functionname, PetscFunction(args)) 436ef1023bdSBarry Smith 4377cdbe19fSJose E. Roman Not Collective; No Fortran Support 438ef1023bdSBarry Smith 439ef1023bdSBarry Smith Input Parameters: 440ef1023bdSBarry Smith + functionname - the name of the function being called, this can be a string with spaces that describes the meaning of the callback 44187497f52SBarry Smith - PetscFunction - user provided callback function that returns an error code 442ef1023bdSBarry Smith 443ef1023bdSBarry Smith Example Usage: 444ef1023bdSBarry Smith .vb 445792fecdfSBarry Smith PetscCallBack("XXX callback to do something", a->callback(...)); 446ef1023bdSBarry Smith .ve 447ef1023bdSBarry Smith 448ef1023bdSBarry Smith Level: developer 449ef1023bdSBarry Smith 450667f096bSBarry Smith Notes: 4517e6f8dd6SBarry Smith `PetscUseTypeMethod()` and ` PetscTryTypeMethod()` are the preferred API for this functionality. But when the callback functions are associated with a 4527e6f8dd6SBarry Smith `DMSNES` or `DMTS` this API must be used. 4537e6f8dd6SBarry Smith 454667f096bSBarry Smith Once the error handler is called the calling function is then returned from with the given 455667f096bSBarry Smith error code. Experienced users can set the error handler with `PetscPushErrorHandler()`. 456667f096bSBarry Smith 457667f096bSBarry Smith `PetscCallBack()` should only be called in PETSc when a call is being made to a user provided call-back routine. 458667f096bSBarry Smith 4597e6f8dd6SBarry Smith Developer Note: 4607e6f8dd6SBarry Smith It would be good to provide a new API for when the callbacks are associated with `DMSNES` or `DMTS` so this routine could be used less 4617e6f8dd6SBarry Smith 46287497f52SBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscCall()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()` 4637e6f8dd6SBarry Smith `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, `CHKERRMPI()`, `PetscCall()`, `PetscUseTypeMethod()`, `PetscTryTypeMethod()` 464ef1023bdSBarry Smith M*/ 465ef1023bdSBarry Smith 4663ba16761SJacob Faibussowitsch /*MC 4678af9f190SBarry Smith PetscCallVoid - Like `PetscCall()` but for use in functions that return `void` 4683ba16761SJacob Faibussowitsch 4693ba16761SJacob Faibussowitsch Synopsis: 4703ba16761SJacob Faibussowitsch #include <petscerror.h> 4718af9f190SBarry Smith void PetscCallVoid(PetscFunction(args)) 4723ba16761SJacob Faibussowitsch 4737cdbe19fSJose E. Roman Not Collective; No Fortran Support 4743ba16761SJacob Faibussowitsch 4753ba16761SJacob Faibussowitsch Input Parameter: 4763ba16761SJacob Faibussowitsch . PetscFunction - any PETSc function that returns an error code 4773ba16761SJacob Faibussowitsch 4783ba16761SJacob Faibussowitsch Example Usage: 4793ba16761SJacob Faibussowitsch .vb 4803ba16761SJacob Faibussowitsch void foo() 4813ba16761SJacob Faibussowitsch { 4823ba16761SJacob Faibussowitsch KSP ksp; 4833ba16761SJacob Faibussowitsch 4843ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 4853ba16761SJacob Faibussowitsch // OK, properly handles PETSc error codes 4863ba16761SJacob Faibussowitsch PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp)); 4878af9f190SBarry Smith PetscFunctionReturnVoid(); 4883ba16761SJacob Faibussowitsch } 4893ba16761SJacob Faibussowitsch 4903ba16761SJacob Faibussowitsch PetscErrorCode bar() 4913ba16761SJacob Faibussowitsch { 4923ba16761SJacob Faibussowitsch KSP ksp; 4933ba16761SJacob Faibussowitsch 4943ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 4953ba16761SJacob Faibussowitsch // ERROR, Non-void function 'bar' should return a value 4963ba16761SJacob Faibussowitsch PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp)); 4973ba16761SJacob Faibussowitsch // OK, returning PetscErrorCode 4983ba16761SJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp)); 4993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5003ba16761SJacob Faibussowitsch } 501667f096bSBarry Smith .ve 5023ba16761SJacob Faibussowitsch 5033ba16761SJacob Faibussowitsch Level: beginner 5043ba16761SJacob Faibussowitsch 505667f096bSBarry Smith Notes: 506667f096bSBarry Smith Has identical usage to `PetscCall()`, except that it returns `void` on error instead of a 507667f096bSBarry Smith `PetscErrorCode`. See `PetscCall()` for more detailed discussion. 508667f096bSBarry Smith 509667f096bSBarry Smith Note that users should prefer `PetscCallAbort()` to this routine. While this routine does 510667f096bSBarry Smith "handle" errors by returning from the enclosing function, it effectively gobbles the 511667f096bSBarry Smith error. Since the enclosing function itself returns `void`, its callers have no way of knowing 512667f096bSBarry Smith that the routine returned early due to an error. `PetscCallAbort()` at least ensures that the 513667f096bSBarry Smith program crashes gracefully. 514667f096bSBarry Smith 515648c30bcSBarry Smith .seealso: `PetscCall()`, `PetscErrorCode`, `PetscCallAbort()`, `PetscCallNull()` 5163ba16761SJacob Faibussowitsch M*/ 5177c5b2466SBarry Smith 5187c5b2466SBarry Smith /*MC 5197c5b2466SBarry Smith PetscCallReturnMPI - Calls a PETSc function and then checks the resulting error code, if it is 5207c5b2466SBarry Smith non-zero it calls the error handler and returns from the current function with an MPI error code. 5217c5b2466SBarry Smith To check for errors in PETSc provided MPI callbacks. 5227c5b2466SBarry Smith 5237c5b2466SBarry Smith Synopsis: 5247c5b2466SBarry Smith #include <petscerror.h> 5257c5b2466SBarry Smith void PetscCallReturnMPI(PetscFunction(args)) 5267c5b2466SBarry Smith 5277c5b2466SBarry Smith Not Collective 5287c5b2466SBarry Smith 5297c5b2466SBarry Smith Input Parameter: 5307c5b2466SBarry Smith . PetscFunction - any PETSc function that returns an error code 5317c5b2466SBarry Smith 5327c5b2466SBarry Smith Level: advanced 5337c5b2466SBarry Smith 5347c5b2466SBarry Smith Notes: 5357c5b2466SBarry Smith Note to be confused with `PetscCallMPI()`. 5367c5b2466SBarry Smith 5377c5b2466SBarry Smith This is be used in a PETSc-provided MPI callback function, such as `MPI_Comm_delete_attr_function function()`. 5387c5b2466SBarry Smith 5397c5b2466SBarry Smith Currently, it always returns `MPI_ERR_OTHER` on failure 5407c5b2466SBarry Smith 5417c5b2466SBarry Smith .seealso: `PetscCall()`, `PetscCallMPI()`, `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`, 5427c5b2466SBarry Smith `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, 5437c5b2466SBarry Smith `CHKERRMPI()`, `PetscCallBack()`, `PetscCallAbort()`, `PetscCallVoid()`, `PetscCallNull()` 5447c5b2466SBarry Smith M*/ 5457c5b2466SBarry Smith 5469566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 5479566063dSJacob Faibussowitsch void PetscCall(PetscErrorCode); 548792fecdfSBarry Smith void PetscCallBack(const char *, PetscErrorCode); 5499566063dSJacob Faibussowitsch void PetscCallVoid(PetscErrorCode); 550648c30bcSBarry Smith void PetscCallNull(PetscErrorCode); 5517c5b2466SBarry Smith void PetscCallReturnMPI(PetscErrorCode); 5529566063dSJacob Faibussowitsch #else 5539371c9d4SSatish Balay #define PetscCall(...) \ 5549371c9d4SSatish Balay do { \ 5553ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_q_; \ 55637154d10SBarry Smith PetscStackUpdateLine; \ 5573ba16761SJacob Faibussowitsch ierr_petsc_call_q_ = __VA_ARGS__; \ 5583ba16761SJacob 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, " "); \ 5599566063dSJacob Faibussowitsch } while (0) 560648c30bcSBarry Smith #define PetscCallNull(...) \ 561648c30bcSBarry Smith do { \ 562648c30bcSBarry Smith PetscErrorCode ierr_petsc_call_q_; \ 563648c30bcSBarry Smith PetscStackUpdateLine; \ 564648c30bcSBarry Smith ierr_petsc_call_q_ = __VA_ARGS__; \ 565648c30bcSBarry Smith if (PetscUnlikely(ierr_petsc_call_q_ != PETSC_SUCCESS)) { \ 566648c30bcSBarry Smith (void)PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); \ 567648c30bcSBarry Smith PetscFunctionReturn(NULL); \ 568648c30bcSBarry Smith } \ 569648c30bcSBarry Smith } while (0) 5709371c9d4SSatish Balay #define PetscCallBack(function, ...) \ 5719371c9d4SSatish Balay do { \ 5723ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_q_; \ 573e33ced7fSLisandro Dalcin PetscStackUpdateLine; \ 574e33ced7fSLisandro Dalcin PetscStackPushExternal(function); \ 5753ba16761SJacob Faibussowitsch ierr_petsc_call_q_ = __VA_ARGS__; \ 576ef1023bdSBarry Smith PetscStackPop; \ 5773ba16761SJacob 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, " "); \ 578ef1023bdSBarry Smith } while (0) 5799371c9d4SSatish Balay #define PetscCallVoid(...) \ 5809371c9d4SSatish Balay do { \ 5813ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_void_; \ 582e33ced7fSLisandro Dalcin PetscStackUpdateLine; \ 5833ba16761SJacob Faibussowitsch ierr_petsc_call_void_ = __VA_ARGS__; \ 5843ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_void_ != PETSC_SUCCESS)) { \ 585dd460d27SBarry Smith (void)PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_void_, PETSC_ERROR_REPEAT, " "); \ 5869371c9d4SSatish Balay return; \ 5879371c9d4SSatish Balay } \ 5889566063dSJacob Faibussowitsch } while (0) 5897c5b2466SBarry Smith #define PetscCallReturnMPI(...) \ 5907c5b2466SBarry Smith do { \ 5917c5b2466SBarry Smith PetscErrorCode ierr_petsc_call_q_; \ 5927c5b2466SBarry Smith PetscStackUpdateLine; \ 5937c5b2466SBarry Smith ierr_petsc_call_q_ = __VA_ARGS__; \ 5947c5b2466SBarry Smith if (PetscUnlikely(ierr_petsc_call_q_ != PETSC_SUCCESS)) { \ 5957c5b2466SBarry Smith (void)PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_q_, PETSC_ERROR_REPEAT, " "); \ 5967c5b2466SBarry Smith return MPI_ERR_OTHER; \ 5977c5b2466SBarry Smith } \ 5987c5b2466SBarry Smith } while (0) 5999566063dSJacob Faibussowitsch #endif 6009566063dSJacob Faibussowitsch 6019566063dSJacob Faibussowitsch /*MC 6029566063dSJacob Faibussowitsch CHKERRQ - Checks error code returned from PETSc function 60330de9b25SBarry Smith 60430de9b25SBarry Smith Synopsis: 605aaa7dc30SBarry Smith #include <petscsys.h> 6069566063dSJacob Faibussowitsch void CHKERRQ(PetscErrorCode ierr) 6079566063dSJacob Faibussowitsch 6089566063dSJacob Faibussowitsch Not Collective 6099566063dSJacob Faibussowitsch 6102fe279fdSBarry Smith Input Parameter: 6119566063dSJacob Faibussowitsch . ierr - nonzero error code 6129566063dSJacob Faibussowitsch 6139566063dSJacob Faibussowitsch Level: deprecated 6149566063dSJacob Faibussowitsch 615667f096bSBarry Smith Note: 616667f096bSBarry Smith Deprecated in favor of `PetscCall()`. This routine behaves identically to it. 617667f096bSBarry Smith 618db781477SPatrick Sanan .seealso: `PetscCall()` 6199566063dSJacob Faibussowitsch M*/ 6209566063dSJacob Faibussowitsch #define CHKERRQ(...) PetscCall(__VA_ARGS__) 6219566063dSJacob Faibussowitsch #define CHKERRV(...) PetscCallVoid(__VA_ARGS__) 6229566063dSJacob Faibussowitsch 623a1fd7ae3SBarry Smith PETSC_EXTERN void PetscMPIErrorString(PetscMPIInt, size_t, char *); 624db9cea48SBarry Smith 6259566063dSJacob Faibussowitsch /*MC 6269566063dSJacob Faibussowitsch PetscCallMPI - Checks error code returned from MPI calls, if non-zero it calls the error 627648c30bcSBarry Smith handler and then returns a `PetscErrorCode` 6289566063dSJacob Faibussowitsch 6299566063dSJacob Faibussowitsch Synopsis: 6309566063dSJacob Faibussowitsch #include <petscerror.h> 63149c86fc7SBarry Smith void PetscCallMPI(MPI_Function(args)) 63230de9b25SBarry Smith 633eca87e8dSBarry Smith Not Collective 63430de9b25SBarry Smith 6352fe279fdSBarry Smith Input Parameter: 63649c86fc7SBarry Smith . MPI_Function - an MPI function that returns an MPI error code 63730de9b25SBarry Smith 638667f096bSBarry Smith Level: beginner 639667f096bSBarry Smith 6409566063dSJacob Faibussowitsch Notes: 64187497f52SBarry Smith Always returns the error code `PETSC_ERR_MPI`; the MPI error code and string are embedded in 6429566063dSJacob Faibussowitsch the string error message. Do not use this to call any other routines (for example PETSc 6433ba16761SJacob Faibussowitsch routines), it should only be used for direct MPI calls. The user may configure PETSc with the 6443ba16761SJacob Faibussowitsch `--with-strict-petscerrorcode` option to check this at compile-time, otherwise they must 6459566063dSJacob Faibussowitsch check this themselves. 6469566063dSJacob Faibussowitsch 647aaa8cc7dSPierre Jolivet This routine can only be used in functions returning `PetscErrorCode` themselves. If the 6483ba16761SJacob Faibussowitsch calling function returns a different type, use `PetscCallMPIAbort()` instead. 6493ba16761SJacob Faibussowitsch 6509566063dSJacob Faibussowitsch Example Usage: 6519566063dSJacob Faibussowitsch .vb 6529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(...)); // OK, calling MPI function 6539566063dSJacob Faibussowitsch 6549566063dSJacob Faibussowitsch PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead! 6559566063dSJacob Faibussowitsch .ve 6569566063dSJacob Faibussowitsch 65749c86fc7SBarry Smith Fortran Notes: 65887497f52SBarry Smith The Fortran function from which this is used must declare a variable `PetscErrorCode` ierr and ierr must be 65949c86fc7SBarry Smith the final argument to the MPI function being called. 66049c86fc7SBarry Smith 66149c86fc7SBarry Smith In the main program and in Fortran subroutines that do not have ierr as the final return parameter one 66287497f52SBarry Smith should use `PetscCallMPIA()` 66349c86fc7SBarry Smith 66449c86fc7SBarry Smith Fortran Usage: 66549c86fc7SBarry Smith .vb 66649c86fc7SBarry Smith PetscErrorCode ierr or integer ierr 66749c86fc7SBarry Smith ... 66849c86fc7SBarry Smith PetscCallMPI(MPI_Comm_size(...,ierr)) 66949c86fc7SBarry Smith PetscCallMPIA(MPI_Comm_size(...,ierr)) ! Will abort after calling error handler 67049c86fc7SBarry Smith 67149c86fc7SBarry Smith PetscCallMPI(MPI_Comm_size(...,eflag)) ! ERROR, final argument must be ierr 67249c86fc7SBarry Smith .ve 67349c86fc7SBarry Smith 6747c5b2466SBarry Smith .seealso: `SETERRMPI()`, `PetscCall()`, `PetscCallReturnMPI()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`, 6757c5b2466SBarry Smith `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 6767c5b2466SBarry Smith `PetscError()`, `CHKMEMQ`, `PetscCallMPINull()` 6777c5b2466SBarry Smith M*/ 6787c5b2466SBarry Smith 6797c5b2466SBarry Smith /*MC 6807c5b2466SBarry Smith PetscCallMPIReturnMPI - Checks error code returned from MPI calls, if non-zero it calls the error 6817c5b2466SBarry Smith handler and then returns an MPI error code. To check for errors in PETSc-provided MPI callbacks. 6827c5b2466SBarry Smith 6837c5b2466SBarry Smith Synopsis: 6847c5b2466SBarry Smith #include <petscerror.h> 6857c5b2466SBarry Smith void PetscCallMPIReturnMPI(MPI_Function(args)) 6867c5b2466SBarry Smith 6877c5b2466SBarry Smith Not Collective 6887c5b2466SBarry Smith 6897c5b2466SBarry Smith Input Parameter: 6907c5b2466SBarry Smith . MPI_Function - an MPI function that returns an MPI error code 6917c5b2466SBarry Smith 6927c5b2466SBarry Smith Level: advanced 6937c5b2466SBarry Smith 6947c5b2466SBarry Smith .seealso: `SETERRMPI()`, `PetscCall()`, `PetscCallMPI()`, `PetscCallReturnMPI()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`, 6953ba16761SJacob Faibussowitsch `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 696648c30bcSBarry Smith `PetscError()`, `CHKMEMQ`, `PetscCallMPINull()` 697648c30bcSBarry Smith M*/ 698648c30bcSBarry Smith 699648c30bcSBarry Smith /*MC 700648c30bcSBarry Smith PetscCallMPINull - Checks error code returned from MPI calls, if non-zero it calls the error 701648c30bcSBarry Smith handler and then returns a `NULL` 702648c30bcSBarry Smith 703648c30bcSBarry Smith Synopsis: 704648c30bcSBarry Smith #include <petscerror.h> 705648c30bcSBarry Smith void PetscCallMPINull(MPI_Function(args)) 706648c30bcSBarry Smith 707648c30bcSBarry Smith Not Collective; No Fortran Support 708648c30bcSBarry Smith 709648c30bcSBarry Smith Input Parameter: 710648c30bcSBarry Smith . MPI_Function - an MPI function that returns an MPI error code 711648c30bcSBarry Smith 712648c30bcSBarry Smith Level: beginner 713648c30bcSBarry Smith 714648c30bcSBarry Smith Notes: 715648c30bcSBarry Smith Always passes the error code `PETSC_ERR_MPI` to the error handler `PetscError()`; the MPI error code and string are embedded in 716648c30bcSBarry Smith the string error message. Do not use this to call any other routines (for example PETSc 717648c30bcSBarry Smith routines), it should only be used for direct MPI calls. 718648c30bcSBarry Smith 719648c30bcSBarry Smith This routine can only be used in functions returning anything that can be returned as a `NULL` themselves. If the 720648c30bcSBarry Smith calling function returns a different type, use `PetscCallMPIAbort()` instead. 721648c30bcSBarry Smith 722648c30bcSBarry Smith Example Usage: 723648c30bcSBarry Smith .vb 724648c30bcSBarry Smith PetscCallMPINull(MPI_Comm_size(...)); // OK, calling MPI function 725648c30bcSBarry Smith 726648c30bcSBarry Smith PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead! 727648c30bcSBarry Smith .ve 728648c30bcSBarry Smith 729648c30bcSBarry Smith .seealso: `SETERRMPI()`, `PetscCall()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`, 730648c30bcSBarry Smith `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 731648c30bcSBarry Smith `PetscError()`, `CHKMEMQ`, `PetscCallMPI()` 7323ba16761SJacob Faibussowitsch M*/ 7333ba16761SJacob Faibussowitsch 7343ba16761SJacob Faibussowitsch /*MC 7353ba16761SJacob Faibussowitsch PetscCallMPIAbort - Like `PetscCallMPI()` but calls `MPI_Abort()` on error 7363ba16761SJacob Faibussowitsch 7373ba16761SJacob Faibussowitsch Synopsis: 7383ba16761SJacob Faibussowitsch #include <petscerror.h> 7393ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm comm, MPI_Function(args)) 7403ba16761SJacob Faibussowitsch 7413ba16761SJacob Faibussowitsch Not Collective 7423ba16761SJacob Faibussowitsch 7433ba16761SJacob Faibussowitsch Input Parameters: 7443ba16761SJacob Faibussowitsch + comm - the MPI communicator to abort on 7453ba16761SJacob Faibussowitsch - MPI_Function - an MPI function that returns an MPI error code 7463ba16761SJacob Faibussowitsch 747667f096bSBarry Smith Level: beginner 748667f096bSBarry Smith 7493ba16761SJacob Faibussowitsch Notes: 7503ba16761SJacob Faibussowitsch Usage is identical to `PetscCallMPI()`. See `PetscCallMPI()` for detailed discussion. 7513ba16761SJacob Faibussowitsch 7523ba16761SJacob Faibussowitsch This routine may be used in functions returning `void` or other non-`PetscErrorCode` types. 7533ba16761SJacob Faibussowitsch 754989c49a2SBarry Smith Fortran Note: 755989c49a2SBarry Smith In Fortran this is called `PetscCallMPIA()` and is intended to be used in the main program while `PetscCallMPI()` is 756989c49a2SBarry Smith used in Fortran subroutines. 757989c49a2SBarry Smith 758989c49a2SBarry Smith Developer Note: 759989c49a2SBarry Smith This should have the same name in Fortran. 760989c49a2SBarry Smith 7613ba16761SJacob Faibussowitsch .seealso: `PetscCallMPI()`, `PetscCallAbort()`, `SETERRABORT()` 76230de9b25SBarry Smith M*/ 7633fcd9f07SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 76463a3b9bcSJacob Faibussowitsch void PetscCallMPI(PetscMPIInt); 7653ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm, PetscMPIInt); 766648c30bcSBarry Smith void PetscCallMPINull(PetscMPIInt); 767064a246eSJacob Faibussowitsch #else 7683ba16761SJacob Faibussowitsch #define PetscCallMPI_Private(__PETSC_STACK_POP_FUNC__, __SETERR_FUNC__, __COMM__, ...) \ 7699371c9d4SSatish Balay do { \ 7703ba16761SJacob Faibussowitsch PetscMPIInt ierr_petsc_call_mpi_; \ 771ef1023bdSBarry Smith PetscStackUpdateLine; \ 772792fecdfSBarry Smith PetscStackPushExternal("MPI function"); \ 773d71ae5a4SJacob Faibussowitsch { \ 7743ba16761SJacob Faibussowitsch ierr_petsc_call_mpi_ = __VA_ARGS__; \ 775d71ae5a4SJacob Faibussowitsch } \ 7763ba16761SJacob Faibussowitsch __PETSC_STACK_POP_FUNC__; \ 7773ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_mpi_ != MPI_SUCCESS)) { \ 7783ba16761SJacob Faibussowitsch char petsc_mpi_7_errorstring[2 * MPI_MAX_ERROR_STRING]; \ 779a1fd7ae3SBarry Smith PetscMPIErrorString(ierr_petsc_call_mpi_, 2 * MPI_MAX_ERROR_STRING, (char *)petsc_mpi_7_errorstring); \ 780a1fd7ae3SBarry Smith __SETERR_FUNC__(__COMM__, PETSC_ERR_MPI, "MPI error %d %s", ierr_petsc_call_mpi_, petsc_mpi_7_errorstring); \ 7813fcd9f07SJacob Faibussowitsch } \ 7823fcd9f07SJacob Faibussowitsch } while (0) 7833ba16761SJacob Faibussowitsch 7843ba16761SJacob Faibussowitsch #define PetscCallMPI(...) PetscCallMPI_Private(PetscStackPop, SETERRQ, PETSC_COMM_SELF, __VA_ARGS__) 7857c5b2466SBarry Smith #define PetscCallMPIReturnMPI(...) PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRMPI, PETSC_COMM_SELF, __VA_ARGS__) 7863ba16761SJacob Faibussowitsch #define PetscCallMPIAbort(comm, ...) PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRABORT, comm, __VA_ARGS__) 787648c30bcSBarry Smith #define PetscCallMPINull(...) PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRQNULL, PETSC_COMM_SELF, __VA_ARGS__) 7889566063dSJacob Faibussowitsch #endif 789064a246eSJacob Faibussowitsch 7907037b0edSPatrick Sanan /*MC 7919566063dSJacob Faibussowitsch CHKERRMPI - Checks error code returned from MPI calls, if non-zero it calls the error 7929566063dSJacob Faibussowitsch handler and then returns 7939566063dSJacob Faibussowitsch 7949566063dSJacob Faibussowitsch Synopsis: 7959566063dSJacob Faibussowitsch #include <petscerror.h> 7969566063dSJacob Faibussowitsch void CHKERRMPI(PetscErrorCode ierr) 7979566063dSJacob Faibussowitsch 7989566063dSJacob Faibussowitsch Not Collective 7999566063dSJacob Faibussowitsch 8009566063dSJacob Faibussowitsch Input Parameter: 8019566063dSJacob Faibussowitsch . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 8029566063dSJacob Faibussowitsch 8039566063dSJacob Faibussowitsch Level: deprecated 8049566063dSJacob Faibussowitsch 805667f096bSBarry Smith Note: 806667f096bSBarry Smith Deprecated in favor of `PetscCallMPI()`. This routine behaves identically to it. 807667f096bSBarry Smith 808db781477SPatrick Sanan .seealso: `PetscCallMPI()` 8099566063dSJacob Faibussowitsch M*/ 8109566063dSJacob Faibussowitsch #define CHKERRMPI(...) PetscCallMPI(__VA_ARGS__) 8119566063dSJacob Faibussowitsch 8129566063dSJacob Faibussowitsch /*MC 813989c49a2SBarry Smith PetscCallAbort - Checks error code returned from PETSc function, if non-zero it aborts immediately by calling `MPI_Abort()` 8149566063dSJacob Faibussowitsch 8159566063dSJacob Faibussowitsch Synopsis: 8169566063dSJacob Faibussowitsch #include <petscerror.h> 8179566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm comm, PetscErrorCode ierr) 8189566063dSJacob Faibussowitsch 819c3339decSBarry Smith Collective 8209566063dSJacob Faibussowitsch 8219566063dSJacob Faibussowitsch Input Parameters: 8229566063dSJacob Faibussowitsch + comm - the MPI communicator on which to abort 8239566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 8249566063dSJacob Faibussowitsch 825667f096bSBarry Smith Level: intermediate 826667f096bSBarry Smith 8279566063dSJacob Faibussowitsch Notes: 82887497f52SBarry Smith This macro has identical type and usage semantics to `PetscCall()` with the important caveat 8299566063dSJacob Faibussowitsch that this macro does not return. Instead, if ierr is nonzero it calls the PETSc error handler 83087497f52SBarry Smith and then immediately calls `MPI_Abort()`. It can therefore be used anywhere. 8319566063dSJacob Faibussowitsch 83287497f52SBarry Smith As per `MPI_Abort()` semantics the communicator passed must be valid, although there is currently 83387497f52SBarry Smith no attempt made at handling any potential errors from `MPI_Abort()`. Note that while 83487497f52SBarry Smith `MPI_Abort()` is required to terminate only those processes which reside on comm, it is often 83587497f52SBarry Smith the case that `MPI_Abort()` terminates *all* processes. 8369566063dSJacob Faibussowitsch 8379566063dSJacob Faibussowitsch Example Usage: 8389566063dSJacob Faibussowitsch .vb 8399566063dSJacob Faibussowitsch PetscErrorCode boom(void) { return PETSC_ERR_MEM; } 8409566063dSJacob Faibussowitsch 8419566063dSJacob Faibussowitsch void foo(void) 8429566063dSJacob Faibussowitsch { 8439566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type 8449566063dSJacob Faibussowitsch } 8459566063dSJacob Faibussowitsch 8469566063dSJacob Faibussowitsch double bar(void) 8479566063dSJacob Faibussowitsch { 8489566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type 8499566063dSJacob Faibussowitsch } 8509566063dSJacob Faibussowitsch 8519566063dSJacob Faibussowitsch PetscCallAbort(MPI_COMM_NULL,boom()); // ERROR, communicator should be valid 8529566063dSJacob Faibussowitsch 8539566063dSJacob Faibussowitsch struct baz 8549566063dSJacob Faibussowitsch { 8559566063dSJacob Faibussowitsch baz() 8569566063dSJacob Faibussowitsch { 8579566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,boom()); // OK 8589566063dSJacob Faibussowitsch } 8599566063dSJacob Faibussowitsch 8609566063dSJacob Faibussowitsch ~baz() 8619566063dSJacob Faibussowitsch { 8629566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,boom()); // OK (in fact the only way to handle PETSc errors) 8639566063dSJacob Faibussowitsch } 8649566063dSJacob Faibussowitsch }; 8659566063dSJacob Faibussowitsch .ve 8669566063dSJacob Faibussowitsch 867989c49a2SBarry Smith Fortran Note: 868989c49a2SBarry Smith Use `PetscCallA()`. 869989c49a2SBarry Smith 870989c49a2SBarry Smith Developer Note: 871989c49a2SBarry Smith This should have the same name in Fortran as in C. 872989c49a2SBarry Smith 873db781477SPatrick Sanan .seealso: `SETERRABORT()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, 8745687f061SJacob Faibussowitsch `SETERRQ()`, `CHKMEMQ`, `PetscCallMPI()`, `PetscCallCXXAbort()` 8759566063dSJacob Faibussowitsch M*/ 8769566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 8779566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm, PetscErrorCode); 8789566063dSJacob Faibussowitsch void PetscCallContinue(PetscErrorCode); 8799566063dSJacob Faibussowitsch #else 8809371c9d4SSatish Balay #define PetscCallAbort(comm, ...) \ 8819371c9d4SSatish Balay do { \ 8827da6b3ccSLisandro Dalcin PetscErrorCode ierr_petsc_call_abort_; \ 8833ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 8847da6b3ccSLisandro Dalcin ierr_petsc_call_abort_ = __VA_ARGS__; \ 8853ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_abort_ != PETSC_SUCCESS)) { \ 8863ba16761SJacob Faibussowitsch ierr_petsc_call_abort_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_abort_, PETSC_ERROR_REPEAT, " "); \ 8873ba16761SJacob Faibussowitsch (void)MPI_Abort(comm, (PetscMPIInt)ierr_petsc_call_abort_); \ 8889566063dSJacob Faibussowitsch } \ 8899566063dSJacob Faibussowitsch } while (0) 8909371c9d4SSatish Balay #define PetscCallContinue(...) \ 8919371c9d4SSatish Balay do { \ 8927da6b3ccSLisandro Dalcin PetscErrorCode ierr_petsc_call_continue_; \ 8933ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 8947da6b3ccSLisandro Dalcin ierr_petsc_call_continue_ = __VA_ARGS__; \ 8953ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_call_continue_ != PETSC_SUCCESS)) { \ 8963ba16761SJacob Faibussowitsch ierr_petsc_call_continue_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_continue_, PETSC_ERROR_REPEAT, " "); \ 8973ba16761SJacob Faibussowitsch (void)ierr_petsc_call_continue_; \ 8983ba16761SJacob Faibussowitsch } \ 8999566063dSJacob Faibussowitsch } while (0) 9009566063dSJacob Faibussowitsch #endif 9019566063dSJacob Faibussowitsch 9029566063dSJacob Faibussowitsch /*MC 9039566063dSJacob Faibussowitsch CHKERRABORT - Checks error code returned from PETSc function. If non-zero it aborts immediately. 9049566063dSJacob Faibussowitsch 9059566063dSJacob Faibussowitsch Synopsis: 9069566063dSJacob Faibussowitsch #include <petscerror.h> 9079566063dSJacob Faibussowitsch void CHKERRABORT(MPI_Comm comm, PetscErrorCode ierr) 9089566063dSJacob Faibussowitsch 9099566063dSJacob Faibussowitsch Not Collective 9109566063dSJacob Faibussowitsch 9119566063dSJacob Faibussowitsch Input Parameters: 9129566063dSJacob Faibussowitsch + comm - the MPI communicator 9139566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 9149566063dSJacob Faibussowitsch 9159566063dSJacob Faibussowitsch Level: deprecated 9169566063dSJacob Faibussowitsch 917667f096bSBarry Smith Note: 918667f096bSBarry Smith Deprecated in favor of `PetscCallAbort()`. This routine behaves identically to it. 919667f096bSBarry Smith 920af27ebaaSBarry Smith .seealso: `PetscCallAbort()`, `PetscErrorCode` 9219566063dSJacob Faibussowitsch M*/ 9229566063dSJacob Faibussowitsch #define CHKERRABORT(comm, ...) PetscCallAbort(comm, __VA_ARGS__) 9239566063dSJacob Faibussowitsch #define CHKERRCONTINUE(...) PetscCallContinue(__VA_ARGS__) 9249566063dSJacob Faibussowitsch 9259566063dSJacob Faibussowitsch /*MC 92687497f52SBarry Smith CHKERRA - Fortran-only replacement for use of `CHKERRQ()` in the main program, which aborts immediately 927f388eb8bSPatrick Sanan 928f388eb8bSPatrick Sanan Synopsis: 929f388eb8bSPatrick Sanan #include <petscsys.h> 930f388eb8bSPatrick Sanan PetscErrorCode CHKERRA(PetscErrorCode ierr) 931f388eb8bSPatrick Sanan 932f388eb8bSPatrick Sanan Not Collective 933f388eb8bSPatrick Sanan 9342fe279fdSBarry Smith Input Parameter: 935f388eb8bSPatrick Sanan . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 936f388eb8bSPatrick Sanan 93787497f52SBarry Smith Level: deprecated 938f388eb8bSPatrick Sanan 93987497f52SBarry Smith Note: 94087497f52SBarry Smith This macro is rarely needed, normal usage is `PetscCallA()` in the main Fortran program. 941f388eb8bSPatrick Sanan 942989c49a2SBarry Smith Developer Note: 943989c49a2SBarry Smith Why isn't this named `CHKERRABORT()` in Fortran? 944989c49a2SBarry Smith 94587497f52SBarry Smith .seealso: `PetscCall()`, `PetscCallA()`, `PetscCallAbort()`, `CHKERRQ()`, `SETERRA()`, `SETERRQ()`, `SETERRABORT()` 946f388eb8bSPatrick Sanan M*/ 947f388eb8bSPatrick Sanan 9481e4f893aSSatish Balay PETSC_EXTERN PetscBool petscwaitonerrorflg; 9491e4f893aSSatish Balay PETSC_EXTERN PetscBool petscindebugger; 95035f00c14SToby Isaac PETSC_EXTERN PetscBool petscabortmpifinalize; 9517c66cc67SJunchao Zhang 9521fc6d13cSAlexander #if defined(PETSC_CLANG_STATIC_ANALYZER) 9531fc6d13cSAlexander void PETSCABORTWITHERR_Private(MPI_Comm, PetscErrorCode); 9541fc6d13cSAlexander #else 9551fc6d13cSAlexander #define PETSCABORTWITHIERR_Private(comm, ierr) \ 9561fc6d13cSAlexander do { \ 9571fc6d13cSAlexander PetscMPIInt size_; \ 958dd460d27SBarry Smith (void)MPI_Comm_size(comm, &size_); \ 9591fc6d13cSAlexander if (PetscCIEnabledPortableErrorOutput && (size_ == PetscGlobalSize || petscabortmpifinalize) && ierr != PETSC_ERR_SIG) { \ 960dd460d27SBarry Smith (void)MPI_Finalize(); \ 9611fc6d13cSAlexander exit(0); \ 9621fc6d13cSAlexander } else if (PetscCIEnabledPortableErrorOutput && PetscGlobalSize == 1) { \ 9631fc6d13cSAlexander exit(0); \ 9641fc6d13cSAlexander } else { \ 965dd460d27SBarry Smith (void)MPI_Abort(comm, (PetscMPIInt)ierr); \ 9661fc6d13cSAlexander } \ 9671fc6d13cSAlexander } while (0) 9681fc6d13cSAlexander #endif 9691fc6d13cSAlexander 9707c66cc67SJunchao Zhang /*MC 971667f096bSBarry Smith PETSCABORT - Call `MPI_Abort()` with an informative error code 9727c66cc67SJunchao Zhang 9737c66cc67SJunchao Zhang Synopsis: 9747c66cc67SJunchao Zhang #include <petscsys.h> 9757c66cc67SJunchao Zhang PETSCABORT(MPI_Comm comm, PetscErrorCode ierr) 9767c66cc67SJunchao Zhang 9777cdbe19fSJose E. Roman Collective; No Fortran Support 9787c66cc67SJunchao Zhang 9797c66cc67SJunchao Zhang Input Parameters: 98095bd0b28SBarry Smith + comm - An MPI communicator, so that the error can be collective 9817c66cc67SJunchao Zhang - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 9827c66cc67SJunchao Zhang 983bf4d2887SBarry Smith Level: advanced 9847c66cc67SJunchao Zhang 985bf4d2887SBarry Smith Notes: 986989c49a2SBarry Smith If the option `-start_in_debugger` was used then this calls `abort()` to stop the program in the debugger. 987bf4d2887SBarry Smith 988989c49a2SBarry Smith if `PetscCIEnabledPortableErrorOutput` is set, which means the code is running in the PETSc test harness (make test), 989989c49a2SBarry Smith and `comm` is `MPI_COMM_WORLD` it strives to exit cleanly without calling `MPI_Abort()` and instead calling `MPI_Finalize()`. 990660278c0SBarry Smith 991989c49a2SBarry Smith This is currently only used when an error propagates up to the C `main()` program and is detected by a `PetscCall()`, `PetscCallMPI()`, 992989c49a2SBarry Smith or is set in `main()` with `SETERRQ()`. Abort calls such as `SETERRABORT()`, 993989c49a2SBarry Smith `PetscCheckAbort()`, `PetscCallMPIAbort()`, and `PetscCallAbort()` always call `MPI_Abort()` and do not have any special 994989c49a2SBarry Smith handling for the test harness. 995989c49a2SBarry Smith 996989c49a2SBarry Smith Developer Note: 997989c49a2SBarry Smith Should the other abort calls also pass through this call instead of calling `MPI_Abort()` directly? 998989c49a2SBarry Smith 999989c49a2SBarry Smith .seealso: `PetscError()`, `PetscCall()`, `SETERRABORT()`, `PetscCheckAbort()`, `PetscCallMPIAbort()`, `PetscCall()`, `PetscCallMPI()`, 1000af27ebaaSBarry Smith `PetscCallAbort()`, `MPI_Abort()`, `PetscErrorCode` 10017c66cc67SJunchao Zhang M*/ 100229f88feaSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 100329f88feaSJacob Faibussowitsch void PETSCABORT(MPI_Comm, PetscErrorCode); 100429f88feaSJacob Faibussowitsch #else 10059371c9d4SSatish Balay #define PETSCABORT(comm, ...) \ 10069371c9d4SSatish Balay do { \ 10073ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_abort_; \ 10083ba16761SJacob Faibussowitsch if (petscwaitonerrorflg) { ierr_petsc_abort_ = PetscSleep(1000); } \ 10093ba16761SJacob Faibussowitsch if (petscindebugger) { \ 10103ba16761SJacob Faibussowitsch abort(); \ 10113ba16761SJacob Faibussowitsch } else { \ 10123ba16761SJacob Faibussowitsch ierr_petsc_abort_ = __VA_ARGS__; \ 10131fc6d13cSAlexander PETSCABORTWITHIERR_Private(comm, ierr_petsc_abort_); \ 10143fcd9f07SJacob Faibussowitsch } \ 10157c66cc67SJunchao Zhang } while (0) 101629f88feaSJacob Faibussowitsch #endif 1017986eef2eSBarry Smith 10189566063dSJacob Faibussowitsch #ifdef PETSC_CLANGUAGE_CXX 1019986eef2eSBarry Smith /*MC 10209566063dSJacob Faibussowitsch PetscCallThrow - Checks error code, if non-zero it calls the C++ error handler which throws 10219566063dSJacob Faibussowitsch an exception 1022986eef2eSBarry Smith 1023986eef2eSBarry Smith Synopsis: 10249566063dSJacob Faibussowitsch #include <petscerror.h> 10259566063dSJacob Faibussowitsch void PetscCallThrow(PetscErrorCode ierr) 1026986eef2eSBarry Smith 1027986eef2eSBarry Smith Not Collective 1028986eef2eSBarry Smith 10299566063dSJacob Faibussowitsch Input Parameter: 1030986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 1031986eef2eSBarry Smith 1032667f096bSBarry Smith Level: beginner 1033667f096bSBarry Smith 1034986eef2eSBarry Smith Notes: 1035af27ebaaSBarry Smith Requires PETSc to be configured with clanguage of c++. Throws a std::runtime_error() on error. 1036986eef2eSBarry Smith 103787497f52SBarry Smith Once the error handler throws the exception you can use `PetscCallVoid()` which returns without 103887497f52SBarry Smith an error code (bad idea since the error is ignored) or `PetscCallAbort()` to have `MPI_Abort()` 10399566063dSJacob Faibussowitsch called immediately. 10409566063dSJacob Faibussowitsch 1041db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscCall()`, `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, 1042db781477SPatrick Sanan `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ` 1043986eef2eSBarry Smith M*/ 10449371c9d4SSatish Balay #define PetscCallThrow(...) \ 10459371c9d4SSatish Balay do { \ 10463ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 10473ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_call_throw_ = __VA_ARGS__; \ 10483ba16761SJacob 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); \ 1049ffc4695bSBarry Smith } while (0) 105085614651SBarry Smith 1051cc26af49SMatthew Knepley /*MC 1052cc26af49SMatthew Knepley CHKERRXX - Checks error code, if non-zero it calls the C++ error handler which throws an exception 1053cc26af49SMatthew Knepley 1054cc26af49SMatthew Knepley Synopsis: 10559566063dSJacob Faibussowitsch #include <petscerror.h> 10563af045c5SBarry Smith void CHKERRXX(PetscErrorCode ierr) 1057cc26af49SMatthew Knepley 1058eca87e8dSBarry Smith Not Collective 1059cc26af49SMatthew Knepley 10609566063dSJacob Faibussowitsch Input Parameter: 10613af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h 1062cc26af49SMatthew Knepley 10639566063dSJacob Faibussowitsch Level: deprecated 1064cc26af49SMatthew Knepley 1065667f096bSBarry Smith Note: 1066667f096bSBarry Smith Deprecated in favor of `PetscCallThrow()`. This routine behaves identically to it. 1067667f096bSBarry Smith 1068db781477SPatrick Sanan .seealso: `PetscCallThrow()` 1069cc26af49SMatthew Knepley M*/ 10709566063dSJacob Faibussowitsch #define CHKERRXX(...) PetscCallThrow(__VA_ARGS__) 1071fd705b32SMatthew Knepley #endif 1072fd705b32SMatthew Knepley 10733ba16761SJacob Faibussowitsch #define PetscCallCXX_Private(__SETERR_FUNC__, __COMM__, ...) \ 10743ba16761SJacob Faibussowitsch do { \ 10753ba16761SJacob Faibussowitsch PetscStackUpdateLine; \ 10763ba16761SJacob Faibussowitsch try { \ 10773ba16761SJacob Faibussowitsch __VA_ARGS__; \ 10783ba16761SJacob Faibussowitsch } catch (const std::exception &e) { \ 10793ba16761SJacob Faibussowitsch __SETERR_FUNC__(__COMM__, PETSC_ERR_LIB, "%s", e.what()); \ 10803ba16761SJacob Faibussowitsch } \ 10813ba16761SJacob Faibussowitsch } while (0) 10823ba16761SJacob Faibussowitsch 10833f520e80SJunchao Zhang /*MC 10849566063dSJacob Faibussowitsch PetscCallCXX - Checks C++ function calls and if they throw an exception, catch it and then 10859566063dSJacob Faibussowitsch return a PETSc error code 10863f520e80SJunchao Zhang 10873f520e80SJunchao Zhang Synopsis: 10889566063dSJacob Faibussowitsch #include <petscerror.h> 10895687f061SJacob Faibussowitsch void PetscCallCXX(...) noexcept; 10903f520e80SJunchao Zhang 10913f520e80SJunchao Zhang Not Collective 10923f520e80SJunchao Zhang 10939566063dSJacob Faibussowitsch Input Parameter: 10945687f061SJacob Faibussowitsch . __VA_ARGS__ - An arbitrary expression 10955687f061SJacob Faibussowitsch 10965687f061SJacob Faibussowitsch Level: beginner 10973f520e80SJunchao Zhang 10983f520e80SJunchao Zhang Notes: 10995687f061SJacob Faibussowitsch `PetscCallCXX(...)` is a macro replacement for 11009566063dSJacob Faibussowitsch .vb 11019566063dSJacob Faibussowitsch try { 11025687f061SJacob Faibussowitsch __VA_ARGS__; 11039566063dSJacob Faibussowitsch } catch (const std::exception& e) { 11049566063dSJacob Faibussowitsch return ConvertToPetscErrorCode(e); 11059566063dSJacob Faibussowitsch } 11069566063dSJacob Faibussowitsch .ve 11079566063dSJacob Faibussowitsch Due to the fact that it catches any (reasonable) exception, it is essentially noexcept. 11083f520e80SJunchao Zhang 11095687f061SJacob Faibussowitsch If you cannot return a `PetscErrorCode` use `PetscCallCXXAbort()` instead. 11105687f061SJacob Faibussowitsch 11119566063dSJacob Faibussowitsch Example Usage: 11129566063dSJacob Faibussowitsch .vb 11139566063dSJacob Faibussowitsch void foo(void) { throw std::runtime_error("error"); } 11143f520e80SJunchao Zhang 11159566063dSJacob Faibussowitsch void bar() 11169566063dSJacob Faibussowitsch { 1117e8952933SJacob Faibussowitsch PetscCallCXX(foo()); // ERROR bar() does not return PetscErrorCode 11189566063dSJacob Faibussowitsch } 11199566063dSJacob Faibussowitsch 11209566063dSJacob Faibussowitsch PetscErrorCode baz() 11219566063dSJacob Faibussowitsch { 11229566063dSJacob Faibussowitsch PetscCallCXX(foo()); // OK 11239566063dSJacob Faibussowitsch 11249566063dSJacob Faibussowitsch PetscCallCXX( 11259566063dSJacob Faibussowitsch bar(); 1126d5b43468SJose E. Roman foo(); // OK multiple statements allowed 11279566063dSJacob Faibussowitsch ); 11289566063dSJacob Faibussowitsch } 11299566063dSJacob Faibussowitsch 11309566063dSJacob Faibussowitsch struct bop 11319566063dSJacob Faibussowitsch { 11329566063dSJacob Faibussowitsch bop() 11339566063dSJacob Faibussowitsch { 1134e8952933SJacob Faibussowitsch PetscCallCXX(foo()); // ERROR returns PetscErrorCode, cannot be used in constructors 11359566063dSJacob Faibussowitsch } 11369566063dSJacob Faibussowitsch }; 11379566063dSJacob Faibussowitsch 1138e8952933SJacob Faibussowitsch // ERROR contains do-while, cannot be used as function-try block 11399566063dSJacob Faibussowitsch PetscErrorCode qux() PetscCallCXX( 11409566063dSJacob Faibussowitsch bar(); 11419566063dSJacob Faibussowitsch baz(); 11429566063dSJacob Faibussowitsch foo(); 11439566063dSJacob Faibussowitsch return 0; 11449566063dSJacob Faibussowitsch ) 11459566063dSJacob Faibussowitsch .ve 11469566063dSJacob Faibussowitsch 11475687f061SJacob Faibussowitsch .seealso: `PetscCallCXXAbort()`, `PetscCallThrow()`, `SETERRQ()`, `PetscCall()`, 11485687f061SJacob Faibussowitsch `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, 11495687f061SJacob Faibussowitsch `PetscError()`, `CHKMEMQ` 11503f520e80SJunchao Zhang M*/ 11513ba16761SJacob Faibussowitsch #define PetscCallCXX(...) PetscCallCXX_Private(SETERRQ, PETSC_COMM_SELF, __VA_ARGS__) 11525687f061SJacob Faibussowitsch 11535687f061SJacob Faibussowitsch /*MC 11545687f061SJacob Faibussowitsch PetscCallCXXAbort - Like `PetscCallCXX()` but calls `MPI_Abort()` instead of returning an 11555687f061SJacob Faibussowitsch error-code 11565687f061SJacob Faibussowitsch 11575687f061SJacob Faibussowitsch Synopsis: 11585687f061SJacob Faibussowitsch #include <petscerror.h> 11595687f061SJacob Faibussowitsch void PetscCallCXXAbort(MPI_Comm comm, ...) noexcept; 11605687f061SJacob Faibussowitsch 116120f4b53cSBarry Smith Collective; No Fortran Support 11625687f061SJacob Faibussowitsch 11635687f061SJacob Faibussowitsch Input Parameters: 11645687f061SJacob Faibussowitsch + comm - The MPI communicator to abort on 11655687f061SJacob Faibussowitsch - __VA_ARGS__ - An arbitrary expression 11665687f061SJacob Faibussowitsch 11675687f061SJacob Faibussowitsch Level: beginner 11685687f061SJacob Faibussowitsch 11695687f061SJacob Faibussowitsch Notes: 11705687f061SJacob Faibussowitsch This macro may be used to check C++ expressions for exceptions in cases where you cannot 11715687f061SJacob Faibussowitsch return an error code. This includes constructors, destructors, copy/move assignment functions 11725687f061SJacob Faibussowitsch or constructors among others. 11735687f061SJacob Faibussowitsch 11745687f061SJacob Faibussowitsch If an exception is caught, the macro calls `SETERRABORT()` on `comm`. The exception must 11755687f061SJacob Faibussowitsch derive from `std::exception` in order to be caught. 11765687f061SJacob Faibussowitsch 11775687f061SJacob Faibussowitsch If the routine _can_ return an error-code it is highly advised to use `PetscCallCXX()` 11785687f061SJacob Faibussowitsch instead. 11795687f061SJacob Faibussowitsch 11805687f061SJacob Faibussowitsch See `PetscCallCXX()` for additional discussion. 11815687f061SJacob Faibussowitsch 11825687f061SJacob Faibussowitsch Example Usage: 11835687f061SJacob Faibussowitsch .vb 11845687f061SJacob Faibussowitsch class Foo 11855687f061SJacob Faibussowitsch { 11865687f061SJacob Faibussowitsch std::vector<int> data_; 11875687f061SJacob Faibussowitsch 11885687f061SJacob Faibussowitsch public: 11895687f061SJacob Faibussowitsch // normally std::vector::reserve() may raise an exception, but since we handle it with 11905687f061SJacob Faibussowitsch // PetscCallCXXAbort() we may mark this routine as noexcept! 11915687f061SJacob Faibussowitsch Foo() noexcept 11925687f061SJacob Faibussowitsch { 11935687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, data_.reserve(10)); 11945687f061SJacob Faibussowitsch } 11955687f061SJacob Faibussowitsch }; 11965687f061SJacob Faibussowitsch 11975687f061SJacob Faibussowitsch std::vector<int> bar() 11985687f061SJacob Faibussowitsch { 11995687f061SJacob Faibussowitsch std::vector<int> v; 12005687f061SJacob Faibussowitsch 12015687f061SJacob Faibussowitsch PetscFunctionBegin; 12025687f061SJacob Faibussowitsch // OK! 12035687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1)); 12045687f061SJacob Faibussowitsch PetscFunctionReturn(v); 12055687f061SJacob Faibussowitsch } 12065687f061SJacob Faibussowitsch 12075687f061SJacob Faibussowitsch PetscErrorCode baz() 12085687f061SJacob Faibussowitsch { 12095687f061SJacob Faibussowitsch std::vector<int> v; 12105687f061SJacob Faibussowitsch 12115687f061SJacob Faibussowitsch PetscFunctionBegin; 12125687f061SJacob Faibussowitsch // WRONG! baz() returns a PetscErrorCode, prefer PetscCallCXX() instead 12135687f061SJacob Faibussowitsch PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1)); 12143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12155687f061SJacob Faibussowitsch } 12165687f061SJacob Faibussowitsch .ve 12175687f061SJacob Faibussowitsch 12185687f061SJacob Faibussowitsch .seealso: `PetscCallCXX()`, `SETERRABORT()`, `PetscCallAbort()` 12195687f061SJacob Faibussowitsch M*/ 12203ba16761SJacob Faibussowitsch #define PetscCallCXXAbort(comm, ...) PetscCallCXX_Private(SETERRABORT, comm, __VA_ARGS__) 12213f520e80SJunchao Zhang 122230de9b25SBarry Smith /*MC 12239566063dSJacob Faibussowitsch CHKERRCXX - Checks C++ function calls and if they throw an exception, catch it and then 12249566063dSJacob Faibussowitsch return a PETSc error code 12259566063dSJacob Faibussowitsch 12269566063dSJacob Faibussowitsch Synopsis: 12279566063dSJacob Faibussowitsch #include <petscerror.h> 12289566063dSJacob Faibussowitsch void CHKERRCXX(func) noexcept; 12299566063dSJacob Faibussowitsch 12309566063dSJacob Faibussowitsch Not Collective 12319566063dSJacob Faibussowitsch 12329566063dSJacob Faibussowitsch Input Parameter: 12339566063dSJacob Faibussowitsch . func - C++ function calls 12349566063dSJacob Faibussowitsch 12359566063dSJacob Faibussowitsch Level: deprecated 12369566063dSJacob Faibussowitsch 1237667f096bSBarry Smith Note: 1238667f096bSBarry Smith Deprecated in favor of `PetscCallCXX()`. This routine behaves identically to it. 1239667f096bSBarry Smith 1240db781477SPatrick Sanan .seealso: `PetscCallCXX()` 12419566063dSJacob Faibussowitsch M*/ 12429566063dSJacob Faibussowitsch #define CHKERRCXX(...) PetscCallCXX(__VA_ARGS__) 12439566063dSJacob Faibussowitsch 12449566063dSJacob Faibussowitsch /*MC 124530de9b25SBarry Smith CHKMEMQ - Checks the memory for corruption, calls error handler if any is detected 124630de9b25SBarry Smith 124730de9b25SBarry Smith Synopsis: 1248aaa7dc30SBarry Smith #include <petscsys.h> 124991d3bdf4SKris Buschelman CHKMEMQ; 125030de9b25SBarry Smith 1251eca87e8dSBarry Smith Not Collective 1252eca87e8dSBarry Smith 125330de9b25SBarry Smith Level: beginner 125430de9b25SBarry Smith 125530de9b25SBarry Smith Notes: 1256af27ebaaSBarry Smith We recommend using Valgrind <https://petsc.org/release/faq/#valgrind> or for NVIDIA CUDA systems 1257af27ebaaSBarry Smith <https://docs.nvidia.com/cuda/cuda-memcheck/index.html> for finding memory problems. The ``CHKMEMQ`` macro is useful on systems that 12585ed36255SBarry Smith do not have valgrind, but is not as good as valgrind or cuda-memcheck. 12591957e957SBarry Smith 1260667f096bSBarry Smith Must run with the option `-malloc_debug` (`-malloc_test` in debug mode; or if `PetscMallocSetDebug()` called) to enable this option 126130de9b25SBarry Smith 126230de9b25SBarry Smith Once the error handler is called the calling function is then returned from with the given error code. 126330de9b25SBarry Smith 126430de9b25SBarry Smith By defaults prints location where memory that is corrupted was allocated. 126530de9b25SBarry Smith 12668af9f190SBarry Smith Use `CHKMEMA` for functions that return `void` 1267f621e05eSBarry Smith 1268db781477SPatrick Sanan .seealso: `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `SETERRQ()`, `PetscMallocValidate()` 126930de9b25SBarry Smith M*/ 12706d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 1271064a246eSJacob Faibussowitsch #define CHKMEMQ 1272064a246eSJacob Faibussowitsch #define CHKMEMA 12736d210af2SJacob Faibussowitsch #else 12749371c9d4SSatish Balay #define CHKMEMQ \ 12759371c9d4SSatish Balay do { \ 12763ba16761SJacob Faibussowitsch PetscErrorCode ierr_petsc_memq_ = PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__); \ 12773ba16761SJacob Faibussowitsch if (PetscUnlikely(ierr_petsc_memq_ != PETSC_SUCCESS)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_memq_, PETSC_ERROR_REPEAT, " "); \ 127886d09637SLisandro Dalcin } while (0) 12796d210af2SJacob Faibussowitsch #define CHKMEMA PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__) 1280064a246eSJacob Faibussowitsch #endif 12819566063dSJacob Faibussowitsch 1282668f157eSBarry Smith /*E 1283668f157eSBarry Smith PetscErrorType - passed to the PETSc error handling routines indicating if this is the first or a later call to the error handlers 1284668f157eSBarry Smith 1285668f157eSBarry Smith Level: advanced 1286668f157eSBarry Smith 1287667f096bSBarry Smith Note: 128887497f52SBarry Smith `PETSC_ERROR_IN_CXX` indicates the error was detected in C++ and an exception should be generated 1289d736bfebSBarry Smith 1290667f096bSBarry Smith Developer Note: 1291667f096bSBarry Smith This is currently used to decide when to print the detailed information about the run in `PetscTraceBackErrorHandler()` 1292668f157eSBarry Smith 129387497f52SBarry Smith .seealso: `PetscError()`, `SETERRQ()` 1294668f157eSBarry Smith E*/ 12959371c9d4SSatish Balay typedef enum { 12969371c9d4SSatish Balay PETSC_ERROR_INITIAL = 0, 12979371c9d4SSatish Balay PETSC_ERROR_REPEAT = 1, 12989371c9d4SSatish Balay PETSC_ERROR_IN_CXX = 2 12999371c9d4SSatish Balay } PetscErrorType; 13004b209cf6SBarry Smith 1301eb9e708aSLisandro Dalcin #if defined(__clang_analyzer__) 1302eb9e708aSLisandro Dalcin __attribute__((analyzer_noreturn)) 1303eb9e708aSLisandro Dalcin #endif 13049371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode 13059371c9d4SSatish Balay PetscError(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, ...) PETSC_ATTRIBUTE_COLD PETSC_ATTRIBUTE_FORMAT(7, 8); 1306eb9e708aSLisandro Dalcin 1307014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscErrorPrintfInitialize(void); 13083ba16761SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscErrorMessage(PetscErrorCode, const char *[], char **); 1309d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscTraceBackErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1310d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscIgnoreErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1311d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscEmacsClientErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1312d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscMPIAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1313d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1314d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAttachDebuggerErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1315d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscReturnErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD; 1316efca3c55SSatish Balay PETSC_EXTERN PetscErrorCode PetscPushErrorHandler(PetscErrorCode (*handler)(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *), void *); 1317014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopErrorHandler(void); 13188d359177SBarry Smith PETSC_EXTERN PetscErrorCode PetscSignalHandlerDefault(int, void *); 1319014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPushSignalHandler(PetscErrorCode (*)(int, void *), void *); 1320014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopSignalHandler(void); 132128559dc8SJed Brown PETSC_EXTERN PetscErrorCode PetscCheckPointerSetIntensity(PetscInt); 1322c2a741eeSJunchao Zhang PETSC_EXTERN void PetscSignalSegvCheckPointerOrMpi(void); 1323edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 13, 0, "PetscSignalSegvCheckPointerOrMpi()", ) static inline void PetscSignalSegvCheckPointer(void) 1324d71ae5a4SJacob Faibussowitsch { 13259371c9d4SSatish Balay PetscSignalSegvCheckPointerOrMpi(); 13269371c9d4SSatish Balay } 1327329f5518SBarry Smith 1328639ff905SBarry Smith /*MC 1329639ff905SBarry Smith PetscErrorPrintf - Prints error messages. 1330639ff905SBarry Smith 1331639ff905SBarry Smith Synopsis: 1332aaa7dc30SBarry Smith #include <petscsys.h> 1333639ff905SBarry Smith PetscErrorCode (*PetscErrorPrintf)(const char format[], ...); 1334639ff905SBarry Smith 13357cdbe19fSJose E. Roman Not Collective; No Fortran Support 13367cdbe19fSJose E. Roman 1337f899ff85SJose E. Roman Input Parameter: 1338cd05f99aSJacob Faibussowitsch . format - the usual `printf()` format string 1339639ff905SBarry Smith 1340639ff905SBarry Smith Options Database Keys: 13411957e957SBarry Smith + -error_output_stdout - cause error messages to be printed to stdout instead of the (default) stderr 1342e1bc860dSBarry Smith - -error_output_none - to turn off all printing of error messages (does not change the way the error is handled.) 1343639ff905SBarry Smith 1344cf53795eSBarry Smith Level: developer 1345cf53795eSBarry Smith 134695452b02SPatrick Sanan Notes: 134795452b02SPatrick Sanan Use 1348667f096bSBarry Smith .vb 134995bd0b28SBarry Smith PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the error is handled) and 1350667f096bSBarry Smith PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on or you can use your own function 1351667f096bSBarry Smith .ve 1352639ff905SBarry Smith Use 1353667f096bSBarry Smith .vb 135487497f52SBarry Smith `PETSC_STDERR` = FILE* obtained from a file open etc. to have stderr printed to the file. 135587497f52SBarry Smith `PETSC_STDOUT` = FILE* obtained from a file open etc. to have stdout printed to the file. 1356667f096bSBarry Smith .ve 1357639ff905SBarry Smith Use 1358af27ebaaSBarry Smith .vb 135987497f52SBarry Smith `PetscPushErrorHandler()` to provide your own error handler that determines what kind of messages to print 1360af27ebaaSBarry Smith .ve 1361639ff905SBarry Smith 1362db781477SPatrick Sanan .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscHelpPrintf()`, `PetscPrintf()`, `PetscPushErrorHandler()`, `PetscVFPrintf()`, `PetscHelpPrintf()` 1363639ff905SBarry Smith M*/ 13643ca90d2dSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode (*PetscErrorPrintf)(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2); 1365639ff905SBarry Smith 1366cf0818bdSBarry Smith /*E 1367cf0818bdSBarry Smith PetscFPTrap - types of floating point exceptions that may be trapped 1368cf0818bdSBarry Smith 136987497f52SBarry Smith Currently only `PETSC_FP_TRAP_OFF` and `PETSC_FP_TRAP_ON` are handled. All others are treated as `PETSC_FP_TRAP_ON`. 1370cf0818bdSBarry Smith 1371cf0818bdSBarry Smith Level: intermediate 1372cf0818bdSBarry Smith 13737de69702SBarry Smith .seealso: `PetscSetFPTrap()`, `PetscFPTrapPush()` 1374cf0818bdSBarry Smith E*/ 13759371c9d4SSatish Balay typedef enum { 13769371c9d4SSatish Balay PETSC_FP_TRAP_OFF = 0, 13779371c9d4SSatish Balay PETSC_FP_TRAP_INDIV = 1, 13789371c9d4SSatish Balay PETSC_FP_TRAP_FLTOPERR = 2, 13799371c9d4SSatish Balay PETSC_FP_TRAP_FLTOVF = 4, 13809371c9d4SSatish Balay PETSC_FP_TRAP_FLTUND = 8, 13819371c9d4SSatish Balay PETSC_FP_TRAP_FLTDIV = 16, 1382dd460d27SBarry Smith PETSC_FP_TRAP_FLTINEX = 32, 1383dd460d27SBarry Smith PETSC_FP_TRAP_ON = (PETSC_FP_TRAP_INDIV | PETSC_FP_TRAP_FLTOPERR | PETSC_FP_TRAP_FLTOVF | PETSC_FP_TRAP_FLTDIV | PETSC_FP_TRAP_FLTINEX) 13849371c9d4SSatish Balay } PetscFPTrap; 1385dd460d27SBarry Smith 1386014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSetFPTrap(PetscFPTrap); 1387014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPush(PetscFPTrap); 1388014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPop(void); 1389aba4c478SBarry Smith PETSC_EXTERN PetscErrorCode PetscDetermineInitialFPTrap(void); 139054a8ef01SBarry Smith 13913a40ed3dSBarry Smith /* 13923a40ed3dSBarry Smith Allows the code to build a stack frame as it runs 13933a40ed3dSBarry Smith */ 13943a40ed3dSBarry Smith 139599cd645aSJed Brown #define PETSCSTACKSIZE 64 13963a40ed3dSBarry Smith typedef struct { 13970e33f6ddSBarry Smith const char *function[PETSCSTACKSIZE]; 13980e33f6ddSBarry Smith const char *file[PETSCSTACKSIZE]; 1399184914b5SBarry Smith int line[PETSCSTACKSIZE]; 1400362febeeSStefano Zampini int petscroutine[PETSCSTACKSIZE]; /* 0 external called from petsc, 1 petsc functions, 2 petsc user functions */ 1401184914b5SBarry Smith int currentsize; 1402a2f94806SJed Brown int hotdepth; 14034be741a6SBarry Smith PetscBool check; /* option to check for correct Push/Pop semantics, true for default petscstack but not other stacks */ 14043a40ed3dSBarry Smith } PetscStack; 1405dfb7d7afSStefano Zampini #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY) 140627104ee2SJacob Faibussowitsch PETSC_EXTERN PetscStack petscstack; 140727104ee2SJacob Faibussowitsch #endif 14083a40ed3dSBarry Smith 14095d12eec7SSatish Balay #if defined(PETSC_SERIALIZE_FUNCTIONS) 14105d12eec7SSatish Balay #include <petsc/private/petscfptimpl.h> 14115d12eec7SSatish Balay /* 14125d12eec7SSatish Balay Registers the current function into the global function pointer to function name table 14135d12eec7SSatish Balay 14145d12eec7SSatish Balay Have to fix this to handle errors but cannot return error since used in PETSC_VIEWER_DRAW_() etc 14155d12eec7SSatish Balay */ 14169371c9d4SSatish Balay #define PetscRegister__FUNCT__() \ 14179371c9d4SSatish Balay do { \ 14185d12eec7SSatish Balay static PetscBool __chked = PETSC_FALSE; \ 14195d12eec7SSatish Balay if (!__chked) { \ 14209371c9d4SSatish Balay void *ptr; \ 14213ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscDLSym(NULL, PETSC_FUNCTION_NAME, &ptr)); \ 14225d12eec7SSatish Balay __chked = PETSC_TRUE; \ 14239371c9d4SSatish Balay } \ 14249371c9d4SSatish Balay } while (0) 14255d12eec7SSatish Balay #else 14265d12eec7SSatish Balay #define PetscRegister__FUNCT__() 14275d12eec7SSatish Balay #endif 14285d12eec7SSatish Balay 1429eae3dc7dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) || defined(__clang_analyzer__) 14306d210af2SJacob Faibussowitsch #define PetscStackPushNoCheck(funct, petsc_routine, hot) 143137154d10SBarry Smith #define PetscStackUpdateLine 1432792fecdfSBarry Smith #define PetscStackPushExternal(funct) 1433dd460d27SBarry Smith #define PetscStackPopNoCheck(funct) 14346d210af2SJacob Faibussowitsch #define PetscStackClearTop 14356d210af2SJacob Faibussowitsch #define PetscFunctionBegin 14366d210af2SJacob Faibussowitsch #define PetscFunctionBeginUser 14376d210af2SJacob Faibussowitsch #define PetscFunctionBeginHot 14380a57284eSJacob Faibussowitsch #define PetscFunctionReturn(...) return __VA_ARGS__ 14396d210af2SJacob Faibussowitsch #define PetscFunctionReturnVoid() return 14406d210af2SJacob Faibussowitsch #define PetscStackPop 14416d210af2SJacob Faibussowitsch #define PetscStackPush(f) 1442dd460d27SBarry Smith #define PetscStackPush_Private(stack__, file__, func__, line__, petsc_routine__, hot__) 1443dd460d27SBarry Smith #define PetscStackPop_Private(stack__, func__) 1444dfb7d7afSStefano Zampini #elif defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY) 1445660278c0SBarry Smith 14469371c9d4SSatish Balay #define PetscStackPush_Private(stack__, file__, func__, line__, petsc_routine__, hot__) \ 14479371c9d4SSatish Balay do { \ 14485a96b57dSJacob Faibussowitsch if (stack__.currentsize < PETSCSTACKSIZE) { \ 14495a96b57dSJacob Faibussowitsch stack__.function[stack__.currentsize] = func__; \ 1450ef1023bdSBarry Smith if (petsc_routine__) { \ 1451ef1023bdSBarry Smith stack__.file[stack__.currentsize] = file__; \ 14525a96b57dSJacob Faibussowitsch stack__.line[stack__.currentsize] = line__; \ 1453ef1023bdSBarry Smith } else { \ 1454648bc8c4SBarry Smith stack__.file[stack__.currentsize] = PETSC_NULLPTR; \ 1455ef1023bdSBarry Smith stack__.line[stack__.currentsize] = 0; \ 1456ef1023bdSBarry Smith } \ 14575a96b57dSJacob Faibussowitsch stack__.petscroutine[stack__.currentsize] = petsc_routine__; \ 14585a96b57dSJacob Faibussowitsch } \ 14595a96b57dSJacob Faibussowitsch ++stack__.currentsize; \ 14605a96b57dSJacob Faibussowitsch stack__.hotdepth += (hot__ || stack__.hotdepth); \ 14615a96b57dSJacob Faibussowitsch } while (0) 14625a96b57dSJacob Faibussowitsch 14634be741a6SBarry Smith /* uses PetscCheckAbort() because may be used in a function that does not return an error code */ 14649371c9d4SSatish Balay #define PetscStackPop_Private(stack__, func__) \ 14659371c9d4SSatish Balay do { \ 14664be741a6SBarry 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__); \ 14675a96b57dSJacob Faibussowitsch if (--stack__.currentsize < PETSCSTACKSIZE) { \ 14689371c9d4SSatish 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", \ 14699371c9d4SSatish Balay stack__.function[stack__.currentsize], stack__.file[stack__.currentsize], stack__.line[stack__.currentsize], func__, __FILE__, __LINE__); \ 14705a96b57dSJacob Faibussowitsch stack__.function[stack__.currentsize] = PETSC_NULLPTR; \ 14715a96b57dSJacob Faibussowitsch stack__.file[stack__.currentsize] = PETSC_NULLPTR; \ 14725a96b57dSJacob Faibussowitsch stack__.line[stack__.currentsize] = 0; \ 14735a96b57dSJacob Faibussowitsch stack__.petscroutine[stack__.currentsize] = 0; \ 14745a96b57dSJacob Faibussowitsch } \ 14755a96b57dSJacob Faibussowitsch stack__.hotdepth = PetscMax(stack__.hotdepth - 1, 0); \ 14765a96b57dSJacob Faibussowitsch } while (0) 14775a96b57dSJacob Faibussowitsch 1478586f9135SBarry Smith /*MC 1479586f9135SBarry Smith PetscStackPushNoCheck - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is 1480586f9135SBarry Smith currently in the source code. 1481586f9135SBarry Smith 1482586f9135SBarry Smith Synopsis: 1483586f9135SBarry Smith #include <petscsys.h> 1484586f9135SBarry Smith void PetscStackPushNoCheck(char *funct,int petsc_routine,PetscBool hot); 1485586f9135SBarry Smith 14867cdbe19fSJose E. Roman Not Collective 14877cdbe19fSJose E. Roman 1488586f9135SBarry Smith Input Parameters: 1489586f9135SBarry Smith + funct - the function name 1490586f9135SBarry Smith . petsc_routine - 2 user function, 1 PETSc function, 0 some other function 1491586f9135SBarry Smith - hot - indicates that the function may be called often so expensive error checking should be turned off inside the function 1492586f9135SBarry Smith 1493586f9135SBarry Smith Level: developer 1494586f9135SBarry Smith 1495586f9135SBarry Smith Notes: 1496586f9135SBarry 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 149787497f52SBarry 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 1498586f9135SBarry Smith help debug the problem. 1499586f9135SBarry Smith 1500ef1023bdSBarry Smith This version does not check the memory corruption (an expensive operation), use `PetscStackPush()` to check the memory. 1501ef1023bdSBarry Smith 1502792fecdfSBarry Smith Use `PetscStackPushExternal()` for a function call that is about to be made to a non-PETSc or user function (such as BLAS etc). 1503ef1023bdSBarry Smith 1504586f9135SBarry Smith The default stack is a global variable called `petscstack`. 1505586f9135SBarry Smith 1506586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1507ef1023bdSBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPush()`, `PetscStackPop`, 1508792fecdfSBarry Smith `PetscStackPushExternal()` 1509586f9135SBarry Smith M*/ 15109371c9d4SSatish Balay #define PetscStackPushNoCheck(funct, petsc_routine, hot) \ 15119371c9d4SSatish Balay do { \ 1512e04113cfSBarry Smith PetscStackSAWsTakeAccess(); \ 15135a96b57dSJacob Faibussowitsch PetscStackPush_Private(petscstack, __FILE__, funct, __LINE__, petsc_routine, hot); \ 1514e04113cfSBarry Smith PetscStackSAWsGrantAccess(); \ 1515441dd030SJed Brown } while (0) 1516441dd030SJed Brown 1517586f9135SBarry Smith /*MC 151887497f52SBarry Smith PetscStackUpdateLine - in a function that has a `PetscFunctionBegin` or `PetscFunctionBeginUser` updates the stack line number to the 151937154d10SBarry Smith current line number. 152037154d10SBarry Smith 152137154d10SBarry Smith Synopsis: 152237154d10SBarry Smith #include <petscsys.h> 152337154d10SBarry Smith void PetscStackUpdateLine 152437154d10SBarry Smith 15257cdbe19fSJose E. Roman Not Collective 15267cdbe19fSJose E. Roman 152737154d10SBarry Smith Level: developer 152837154d10SBarry Smith 152937154d10SBarry Smith Notes: 153087497f52SBarry Smith Using `PetscCall()` and friends automatically handles this process 153187497f52SBarry Smith 153237154d10SBarry 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 153337154d10SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 153437154d10SBarry Smith help debug the problem. 153537154d10SBarry Smith 153695bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 153737154d10SBarry Smith 153837154d10SBarry Smith This is used by `PetscCall()` and is otherwise not like to be needed 153937154d10SBarry Smith 154037154d10SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`, `PetscCall()` 154137154d10SBarry Smith M*/ 154237154d10SBarry Smith #define PetscStackUpdateLine \ 15433ba16761SJacob Faibussowitsch do { \ 1544f1e99121SPierre Jolivet if (petscstack.currentsize > 0 && petscstack.currentsize < PETSCSTACKSIZE && petscstack.function[petscstack.currentsize - 1] == PETSC_FUNCTION_NAME) { petscstack.line[petscstack.currentsize - 1] = __LINE__; } \ 15453ba16761SJacob Faibussowitsch } while (0) 154637154d10SBarry Smith 154737154d10SBarry Smith /*MC 1548792fecdfSBarry Smith PetscStackPushExternal - Pushes a new function name onto the PETSc default stack that tracks where the running program is 1549ef1023bdSBarry Smith currently in the source code. Does not include the filename or line number since this is called by the calling routine 1550ef1023bdSBarry Smith for non-PETSc or user functions. 1551ef1023bdSBarry Smith 1552ef1023bdSBarry Smith Synopsis: 1553ef1023bdSBarry Smith #include <petscsys.h> 1554792fecdfSBarry Smith void PetscStackPushExternal(char *funct); 1555ef1023bdSBarry Smith 15567cdbe19fSJose E. Roman Not Collective 15577cdbe19fSJose E. Roman 15582fe279fdSBarry Smith Input Parameter: 1559ef1023bdSBarry Smith . funct - the function name 1560ef1023bdSBarry Smith 1561ef1023bdSBarry Smith Level: developer 1562ef1023bdSBarry Smith 1563ef1023bdSBarry Smith Notes: 156487497f52SBarry Smith Using `PetscCallExternal()` and friends automatically handles this process 156587497f52SBarry Smith 1566ef1023bdSBarry 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 1567ef1023bdSBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1568ef1023bdSBarry Smith help debug the problem. 1569ef1023bdSBarry Smith 1570ef1023bdSBarry Smith The default stack is a global variable called `petscstack`. 1571ef1023bdSBarry Smith 1572ef1023bdSBarry Smith This is to be used when calling an external package function such as a BLAS function. 1573ef1023bdSBarry Smith 1574ef1023bdSBarry Smith This also updates the stack line number for the current stack function. 1575ef1023bdSBarry Smith 1576ef1023bdSBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1577ef1023bdSBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1578ef1023bdSBarry Smith M*/ 15799371c9d4SSatish Balay #define PetscStackPushExternal(funct) \ 15809371c9d4SSatish Balay do { \ 15819371c9d4SSatish Balay PetscStackUpdateLine; \ 15829371c9d4SSatish Balay PetscStackPushNoCheck(funct, 0, PETSC_TRUE); \ 1583a8f51744SPierre Jolivet } while (0) 1584ef1023bdSBarry Smith 1585ef1023bdSBarry Smith /*MC 1586586f9135SBarry Smith PetscStackPopNoCheck - Pops a function name from the PETSc default stack that tracks where the running program is 1587586f9135SBarry Smith currently in the source code. 1588586f9135SBarry Smith 1589586f9135SBarry Smith Synopsis: 1590586f9135SBarry Smith #include <petscsys.h> 1591586f9135SBarry Smith void PetscStackPopNoCheck(char *funct); 1592586f9135SBarry Smith 15937cdbe19fSJose E. Roman Not Collective 15947cdbe19fSJose E. Roman 1595586f9135SBarry Smith Input Parameter: 1596586f9135SBarry Smith . funct - the function name 1597586f9135SBarry Smith 1598586f9135SBarry Smith Level: developer 1599586f9135SBarry Smith 1600586f9135SBarry Smith Notes: 160187497f52SBarry Smith Using `PetscCall()`, `PetscCallExternal()`, `PetscCallBack()` and friends negates the need to call this 160287497f52SBarry Smith 1603586f9135SBarry 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 1604586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1605586f9135SBarry Smith help debug the problem. 1606586f9135SBarry Smith 160795bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 1608586f9135SBarry Smith 1609586f9135SBarry Smith Developer Note: 1610586f9135SBarry Smith `PetscStackPopNoCheck()` takes a function argument while `PetscStackPop` does not, this difference is likely just historical. 1611586f9135SBarry Smith 1612586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1613586f9135SBarry Smith M*/ 16149371c9d4SSatish Balay #define PetscStackPopNoCheck(funct) \ 16159371c9d4SSatish Balay do { \ 1616362febeeSStefano Zampini PetscStackSAWsTakeAccess(); \ 16175a96b57dSJacob Faibussowitsch PetscStackPop_Private(petscstack, funct); \ 1618362febeeSStefano Zampini PetscStackSAWsGrantAccess(); \ 1619362febeeSStefano Zampini } while (0) 1620362febeeSStefano Zampini 16219371c9d4SSatish Balay #define PetscStackClearTop \ 16229371c9d4SSatish Balay do { \ 1623e04113cfSBarry Smith PetscStackSAWsTakeAccess(); \ 16249371c9d4SSatish Balay if (petscstack.currentsize > 0 && --petscstack.currentsize < PETSCSTACKSIZE) { \ 162527104ee2SJacob Faibussowitsch petscstack.function[petscstack.currentsize] = PETSC_NULLPTR; \ 162627104ee2SJacob Faibussowitsch petscstack.file[petscstack.currentsize] = PETSC_NULLPTR; \ 162727104ee2SJacob Faibussowitsch petscstack.line[petscstack.currentsize] = 0; \ 162827104ee2SJacob Faibussowitsch petscstack.petscroutine[petscstack.currentsize] = 0; \ 1629441dd030SJed Brown } \ 163027104ee2SJacob Faibussowitsch petscstack.hotdepth = PetscMax(petscstack.hotdepth - 1, 0); \ 1631e04113cfSBarry Smith PetscStackSAWsGrantAccess(); \ 1632441dd030SJed Brown } while (0) 1633441dd030SJed Brown 163430de9b25SBarry Smith /*MC 16351957e957SBarry Smith PetscFunctionBegin - First executable line of each PETSc function, used for error handling. Final 163687497f52SBarry Smith line of PETSc functions should be `PetscFunctionReturn`(0); 163730de9b25SBarry Smith 163830de9b25SBarry Smith Synopsis: 1639aaa7dc30SBarry Smith #include <petscsys.h> 164030de9b25SBarry Smith void PetscFunctionBegin; 164130de9b25SBarry Smith 164220f4b53cSBarry Smith Not Collective; No Fortran Support 1643eca87e8dSBarry Smith 164430de9b25SBarry Smith Usage: 164530de9b25SBarry Smith .vb 164630de9b25SBarry Smith int something; 164730de9b25SBarry Smith 164830de9b25SBarry Smith PetscFunctionBegin; 164930de9b25SBarry Smith .ve 165030de9b25SBarry Smith 165130de9b25SBarry Smith Level: developer 165230de9b25SBarry Smith 165320f4b53cSBarry Smith Note: 165420f4b53cSBarry Smith Use `PetscFunctionBeginUser` for application codes. 165520f4b53cSBarry Smith 1656586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()` 165730de9b25SBarry Smith 165830de9b25SBarry Smith M*/ 16599371c9d4SSatish Balay #define PetscFunctionBegin \ 16609371c9d4SSatish Balay do { \ 1661362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_FALSE); \ 1662a2f94806SJed Brown PetscRegister__FUNCT__(); \ 1663a2f94806SJed Brown } while (0) 1664a2f94806SJed Brown 1665a2f94806SJed Brown /*MC 166687497f52SBarry Smith PetscFunctionBeginHot - Substitute for `PetscFunctionBegin` to be used in functions that are called in 1667a2f94806SJed Brown performance-critical circumstances. Use of this function allows for lighter profiling by default. 1668a2f94806SJed Brown 1669a2f94806SJed Brown Synopsis: 1670aaa7dc30SBarry Smith #include <petscsys.h> 1671a2f94806SJed Brown void PetscFunctionBeginHot; 1672a2f94806SJed Brown 167320f4b53cSBarry Smith Not Collective; No Fortran Support 1674a2f94806SJed Brown 1675a2f94806SJed Brown Usage: 1676a2f94806SJed Brown .vb 1677a2f94806SJed Brown int something; 1678a2f94806SJed Brown 1679a2f94806SJed Brown PetscFunctionBeginHot; 1680a2f94806SJed Brown .ve 1681a2f94806SJed Brown 1682a2f94806SJed Brown Level: developer 1683a2f94806SJed Brown 1684586f9135SBarry Smith .seealso: `PetscFunctionBegin`, `PetscFunctionReturn()`, `PetscStackPushNoCheck()` 1685a2f94806SJed Brown 1686a2f94806SJed Brown M*/ 16879371c9d4SSatish Balay #define PetscFunctionBeginHot \ 16889371c9d4SSatish Balay do { \ 1689362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_TRUE); \ 16902d53ad75SBarry Smith PetscRegister__FUNCT__(); \ 169153c77d0aSJed Brown } while (0) 169253c77d0aSJed Brown 1693a8d2bbe5SBarry Smith /*MC 1694530d85c1SBarry Smith PetscFunctionBeginUser - First executable line of user provided routines 1695a8d2bbe5SBarry Smith 1696a8d2bbe5SBarry Smith Synopsis: 1697aaa7dc30SBarry Smith #include <petscsys.h> 1698a8d2bbe5SBarry Smith void PetscFunctionBeginUser; 1699a8d2bbe5SBarry Smith 170020f4b53cSBarry Smith Not Collective; No Fortran Support 1701a8d2bbe5SBarry Smith 1702a8d2bbe5SBarry Smith Usage: 1703a8d2bbe5SBarry Smith .vb 1704a8d2bbe5SBarry Smith int something; 1705a8d2bbe5SBarry Smith 1706ac285190SBarry Smith PetscFunctionBeginUser; 1707a8d2bbe5SBarry Smith .ve 1708a8d2bbe5SBarry Smith 170920f4b53cSBarry Smith Level: intermediate 171020f4b53cSBarry Smith 1711a8d2bbe5SBarry Smith Notes: 1712530d85c1SBarry Smith Functions that incorporate this must call `PetscFunctionReturn()` instead of return except for main(). 1713530d85c1SBarry Smith 1714530d85c1SBarry Smith May be used before `PetscInitialize()` 17151957e957SBarry Smith 1716530d85c1SBarry Smith This is identical to `PetscFunctionBegin` except it labels the routine as a user 1717ac285190SBarry Smith routine instead of as a PETSc library routine. 1718ac285190SBarry Smith 1719586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, `PetscFunctionBeginHot`, `PetscStackPushNoCheck()` 1720a8d2bbe5SBarry Smith M*/ 17219371c9d4SSatish Balay #define PetscFunctionBeginUser \ 17229371c9d4SSatish Balay do { \ 1723362febeeSStefano Zampini PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 2, PETSC_FALSE); \ 1724a8d2bbe5SBarry Smith PetscRegister__FUNCT__(); \ 1725a8d2bbe5SBarry Smith } while (0) 1726a8d2bbe5SBarry Smith 1727586f9135SBarry Smith /*MC 1728586f9135SBarry Smith PetscStackPush - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is 1729586f9135SBarry Smith currently in the source code and verifies the memory is not corrupted. 1730586f9135SBarry Smith 1731586f9135SBarry Smith Synopsis: 1732586f9135SBarry Smith #include <petscsys.h> 1733586f9135SBarry Smith void PetscStackPush(char *funct) 1734586f9135SBarry Smith 17357cdbe19fSJose E. Roman Not Collective 17367cdbe19fSJose E. Roman 1737586f9135SBarry Smith Input Parameter: 1738586f9135SBarry Smith . funct - the function name 1739586f9135SBarry Smith 1740586f9135SBarry Smith Level: developer 1741586f9135SBarry Smith 1742586f9135SBarry Smith Notes: 1743586f9135SBarry 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 1744586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1745586f9135SBarry Smith help debug the problem. 1746586f9135SBarry Smith 174795bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 1748586f9135SBarry Smith 1749586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`, 1750586f9135SBarry Smith `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop` 1751586f9135SBarry Smith M*/ 17529371c9d4SSatish Balay #define PetscStackPush(n) \ 17539371c9d4SSatish Balay do { \ 1754362febeeSStefano Zampini PetscStackPushNoCheck(n, 0, PETSC_FALSE); \ 175515681b3cSBarry Smith CHKMEMQ; \ 175615681b3cSBarry Smith } while (0) 17573a40ed3dSBarry Smith 1758586f9135SBarry Smith /*MC 1759586f9135SBarry Smith PetscStackPop - Pops a function name from the PETSc default stack that tracks where the running program is 1760586f9135SBarry Smith currently in the source code and verifies the memory is not corrupted. 1761586f9135SBarry Smith 1762586f9135SBarry Smith Synopsis: 1763586f9135SBarry Smith #include <petscsys.h> 1764586f9135SBarry Smith void PetscStackPop 1765586f9135SBarry Smith 17667cdbe19fSJose E. Roman Not Collective 17677cdbe19fSJose E. Roman 1768586f9135SBarry Smith Level: developer 1769586f9135SBarry Smith 1770586f9135SBarry Smith Notes: 1771586f9135SBarry 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 1772586f9135SBarry Smith occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to 1773586f9135SBarry Smith help debug the problem. 1774586f9135SBarry Smith 177595bd0b28SBarry Smith The default stack is a global variable called `petscstack`. 1776586f9135SBarry Smith 1777586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPopNoCheck()`, `PetscStackPush()` 1778586f9135SBarry Smith M*/ 17799371c9d4SSatish Balay #define PetscStackPop \ 17809371c9d4SSatish Balay do { \ 1781441dd030SJed Brown CHKMEMQ; \ 1782362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 178315681b3cSBarry Smith } while (0) 1784d64ed03dSBarry Smith 178530de9b25SBarry Smith /*MC 17860a57284eSJacob Faibussowitsch PetscFunctionReturn - Last executable line of each PETSc function used for error 17870a57284eSJacob Faibussowitsch handling. Replaces `return()`. 178830de9b25SBarry Smith 178930de9b25SBarry Smith Synopsis: 17900a57284eSJacob Faibussowitsch #include <petscerror.h> 17910a57284eSJacob Faibussowitsch void PetscFunctionReturn(...) 179230de9b25SBarry Smith 1793667f096bSBarry Smith Not Collective; No Fortran Support 1794eca87e8dSBarry Smith 17955687f061SJacob Faibussowitsch Level: beginner 17965687f061SJacob Faibussowitsch 17970a57284eSJacob Faibussowitsch Notes: 17980a57284eSJacob Faibussowitsch This routine is a macro, so while it does not "return" anything itself, it does return from 17990a57284eSJacob Faibussowitsch the function in the literal sense. 18000a57284eSJacob Faibussowitsch 18010a57284eSJacob Faibussowitsch Usually the return value is the integer literal `0` (for example in any function returning 18020a57284eSJacob Faibussowitsch `PetscErrorCode`), however it is possible to return any arbitrary type. The arguments of 18030a57284eSJacob Faibussowitsch this macro are placed before the `return` statement as-is. 18040a57284eSJacob Faibussowitsch 18050a57284eSJacob Faibussowitsch Any routine which returns via `PetscFunctionReturn()` must begin with a corresponding 18060a57284eSJacob Faibussowitsch `PetscFunctionBegin`. 18070a57284eSJacob Faibussowitsch 18080a57284eSJacob Faibussowitsch For routines which return `void` use `PetscFunctionReturnVoid()` instead. 18090a57284eSJacob Faibussowitsch 18100a57284eSJacob Faibussowitsch Example Usage: 181130de9b25SBarry Smith .vb 18120a57284eSJacob Faibussowitsch PetscErrorCode foo(int *x) 18130a57284eSJacob Faibussowitsch { 18140a57284eSJacob Faibussowitsch PetscFunctionBegin; // don't forget the begin! 18150a57284eSJacob Faibussowitsch *x = 10; 18163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 181730de9b25SBarry Smith } 181830de9b25SBarry Smith .ve 181930de9b25SBarry Smith 18200a57284eSJacob Faibussowitsch May return any arbitrary type\: 18210a57284eSJacob Faibussowitsch .vb 18220a57284eSJacob Faibussowitsch struct Foo 18230a57284eSJacob Faibussowitsch { 18240a57284eSJacob Faibussowitsch int x; 18250a57284eSJacob Faibussowitsch }; 18260a57284eSJacob Faibussowitsch 18270a57284eSJacob Faibussowitsch struct Foo make_foo(int value) 18280a57284eSJacob Faibussowitsch { 18290a57284eSJacob Faibussowitsch struct Foo f; 18300a57284eSJacob Faibussowitsch 18310a57284eSJacob Faibussowitsch PetscFunctionBegin; 18320a57284eSJacob Faibussowitsch f.x = value; 18330a57284eSJacob Faibussowitsch PetscFunctionReturn(f); 18340a57284eSJacob Faibussowitsch } 18350a57284eSJacob Faibussowitsch .ve 18360a57284eSJacob Faibussowitsch 18370a57284eSJacob Faibussowitsch .seealso: `PetscFunctionBegin`, `PetscFunctionBeginUser`, `PetscFunctionReturnVoid()`, 18380a57284eSJacob Faibussowitsch `PetscStackPopNoCheck()` 183930de9b25SBarry Smith M*/ 18405687f061SJacob Faibussowitsch #define PetscFunctionReturn(...) \ 18419371c9d4SSatish Balay do { \ 1842362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 18435687f061SJacob Faibussowitsch return __VA_ARGS__; \ 184427104ee2SJacob Faibussowitsch } while (0) 1845d64ed03dSBarry Smith 18460a57284eSJacob Faibussowitsch /*MC 18470a57284eSJacob Faibussowitsch PetscFunctionReturnVoid - Like `PetscFunctionReturn()` but returns `void` 18480a57284eSJacob Faibussowitsch 18490a57284eSJacob Faibussowitsch Synopsis: 18500a57284eSJacob Faibussowitsch #include <petscerror.h> 18510a57284eSJacob Faibussowitsch void PetscFunctionReturnVoid() 18520a57284eSJacob Faibussowitsch 18530a57284eSJacob Faibussowitsch Not Collective 18540a57284eSJacob Faibussowitsch 18555687f061SJacob Faibussowitsch Level: beginner 18565687f061SJacob Faibussowitsch 18575687f061SJacob Faibussowitsch Note: 18580a57284eSJacob Faibussowitsch Behaves identically to `PetscFunctionReturn()` except that it returns `void`. That is, this 18595687f061SJacob Faibussowitsch macro culminates with `return`. 18600a57284eSJacob Faibussowitsch 18610a57284eSJacob Faibussowitsch Example Usage: 18620a57284eSJacob Faibussowitsch .vb 18630a57284eSJacob Faibussowitsch void foo() 18640a57284eSJacob Faibussowitsch { 18650a57284eSJacob Faibussowitsch PetscFunctionBegin; // must start with PetscFunctionBegin! 18660a57284eSJacob Faibussowitsch bar(); 18670a57284eSJacob Faibussowitsch baz(); 18680a57284eSJacob Faibussowitsch PetscFunctionReturnVoid(); 18690a57284eSJacob Faibussowitsch } 18700a57284eSJacob Faibussowitsch .ve 18710a57284eSJacob Faibussowitsch 18720a57284eSJacob Faibussowitsch .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, PetscFunctionBeginUser` 18730a57284eSJacob Faibussowitsch M*/ 18749371c9d4SSatish Balay #define PetscFunctionReturnVoid() \ 18759371c9d4SSatish Balay do { \ 1876362febeeSStefano Zampini PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \ 187727104ee2SJacob Faibussowitsch return; \ 187827104ee2SJacob Faibussowitsch } while (0) 187927104ee2SJacob Faibussowitsch #else /* PETSC_USE_DEBUG */ 188027104ee2SJacob Faibussowitsch #define PetscStackPushNoCheck(funct, petsc_routine, hot) 188137154d10SBarry Smith #define PetscStackUpdateLine 1882792fecdfSBarry Smith #define PetscStackPushExternal(funct) 18833ba16761SJacob Faibussowitsch #define PetscStackPopNoCheck(...) 188427104ee2SJacob Faibussowitsch #define PetscStackClearTop 18853a40ed3dSBarry Smith #define PetscFunctionBegin 18860bdf7c52SPeter Brune #define PetscFunctionBeginUser 1887a2f94806SJed Brown #define PetscFunctionBeginHot 18880a57284eSJacob Faibussowitsch #define PetscFunctionReturn(...) return __VA_ARGS__ 18895665465eSBarry Smith #define PetscFunctionReturnVoid() return 1890812af9f3SBarry Smith #define PetscStackPop CHKMEMQ 1891812af9f3SBarry Smith #define PetscStackPush(f) CHKMEMQ 189227104ee2SJacob Faibussowitsch #endif /* PETSC_USE_DEBUG */ 18933a40ed3dSBarry Smith 18946d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 18953ba16761SJacob Faibussowitsch #define PetscStackCallExternalVoid(...) 18963ba16761SJacob Faibussowitsch template <typename F, typename... Args> 18973ba16761SJacob Faibussowitsch void PetscCallExternal(F, Args...); 18981fc6d13cSAlexander template <typename F, typename... Args> 18991fc6d13cSAlexander void PetscCallExternalAbort(F, Args...); 19006d210af2SJacob Faibussowitsch #else 1901586f9135SBarry Smith /*MC 1902e77caa6dSBarry Smith PetscStackCallExternalVoid - Calls an external library routine or user function after pushing the name of the routine on the stack. 1903eb6b5d47SBarry Smith 1904eb6b5d47SBarry Smith Input Parameters: 1905eb6b5d47SBarry Smith + name - string that gives the name of the function being called 1906586f9135SBarry Smith - routine - actual call to the routine, for example, functionname(a,b) 1907fd3f9acdSBarry Smith 1908586f9135SBarry Smith Level: developer 1909eb6b5d47SBarry Smith 191095bd0b28SBarry Smith Notes: 1911792fecdfSBarry Smith Often one should use `PetscCallExternal()` instead. This routine is intended for external library routines that DO NOT return error codes 1912eb6b5d47SBarry Smith 1913586f9135SBarry Smith In debug mode this also checks the memory for corruption at the end of the function call. 1914586f9135SBarry Smith 191595bd0b28SBarry Smith Certain external packages, such as BLAS/LAPACK may have their own macros, `PetscCallBLAS()` for managing the call, error checking, etc. 1916586f9135SBarry Smith 1917586f9135SBarry Smith Developer Note: 1918586f9135SBarry 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. 1919586f9135SBarry Smith 1920792fecdfSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscCallExternal()`, `PetscCallBLAS()` 1921586f9135SBarry Smith @*/ 19223ba16761SJacob Faibussowitsch #define PetscStackCallExternalVoid(name, ...) \ 19239371c9d4SSatish Balay do { \ 192428770333SStefano Zampini PetscStackPushExternal(name); \ 19253ba16761SJacob Faibussowitsch __VA_ARGS__; \ 19269371c9d4SSatish Balay PetscStackPop; \ 19279371c9d4SSatish Balay } while (0) 1928eb6b5d47SBarry Smith 1929586f9135SBarry Smith /*MC 1930792fecdfSBarry Smith PetscCallExternal - Calls an external library routine that returns an error code after pushing the name of the routine on the stack. 1931fd3f9acdSBarry Smith 1932fd3f9acdSBarry Smith Input Parameters: 1933fd3f9acdSBarry Smith + func - name of the routine 1934586f9135SBarry Smith - args - arguments to the routine 1935586f9135SBarry Smith 1936586f9135SBarry Smith Level: developer 1937fd3f9acdSBarry Smith 193895452b02SPatrick Sanan Notes: 1939e77caa6dSBarry Smith This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not. 1940dbf62e16SBarry Smith 1941586f9135SBarry Smith In debug mode this also checks the memory for corruption at the end of the function call. 1942fd3f9acdSBarry Smith 194387497f52SBarry Smith Assumes the error return code of the function is an integer and that a value of 0 indicates success 194487497f52SBarry Smith 1945586f9135SBarry Smith Developer Note: 1946d5b43468SJose E. Roman This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc. 1947586f9135SBarry Smith 19481fc6d13cSAlexander .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()`, `PetscCallExternalAbort()` 1949586f9135SBarry Smith M*/ 19509371c9d4SSatish Balay #define PetscCallExternal(func, ...) \ 19519371c9d4SSatish Balay do { \ 1952a74df02fSJacob Faibussowitsch PetscStackPush(PetscStringize(func)); \ 1953*9ad2cedaSPierre Jolivet int ierr_petsc_call_external_ = (int)func(__VA_ARGS__); \ 19541d4906efSStefano Zampini PetscStackPop; \ 19553ba16761SJacob 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_); \ 1956fd3f9acdSBarry Smith } while (0) 19571fc6d13cSAlexander 19581fc6d13cSAlexander /*MC 19591fc6d13cSAlexander PetscCallExternalAbort - Calls an external library routine that returns an error code after pushing the name of the routine on the stack. If the external library function return code indicates an error, this prints the error and aborts 19601fc6d13cSAlexander 19611fc6d13cSAlexander Input Parameters: 19621fc6d13cSAlexander + func - name of the routine 19631fc6d13cSAlexander - args - arguments to the routine 19641fc6d13cSAlexander 19651fc6d13cSAlexander Level: developer 19661fc6d13cSAlexander 19671fc6d13cSAlexander Notes: 19681fc6d13cSAlexander This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not. 19691fc6d13cSAlexander 19701fc6d13cSAlexander In debug mode this also checks the memory for corruption at the end of the function call. 19711fc6d13cSAlexander 19721fc6d13cSAlexander Assumes the error return code of the function is an integer and that a value of 0 indicates success 19731fc6d13cSAlexander 19741fc6d13cSAlexander Developer Note: 19751fc6d13cSAlexander This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc. 19761fc6d13cSAlexander 19771fc6d13cSAlexander .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()`, `PetscCallExternal()` 19781fc6d13cSAlexander M*/ 19791fc6d13cSAlexander #define PetscCallExternalAbort(func, ...) \ 19801fc6d13cSAlexander do { \ 19811fc6d13cSAlexander PetscStackUpdateLine; \ 19821fc6d13cSAlexander int ierr_petsc_call_external_ = func(__VA_ARGS__); \ 19831fc6d13cSAlexander if (PetscUnlikely(ierr_petsc_call_external_ != 0)) { \ 19841fc6d13cSAlexander (void)PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_LIB, PETSC_ERROR_INITIAL, "Error in %s(): error code %d", PetscStringize(func), ierr_petsc_call_external_); \ 19851fc6d13cSAlexander PETSCABORTWITHIERR_Private(PETSC_COMM_SELF, PETSC_ERR_LIB); \ 19861fc6d13cSAlexander } \ 19871fc6d13cSAlexander } while (0) 19886d210af2SJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */ 1989