xref: /petsc/include/petscerror.h (revision 139c8099082a12dfcf688dbf7e6afbfe11a33acd)
154a8ef01SBarry Smith /*
2f621e05eSBarry Smith     Contains all error handling interfaces for PETSc.
354a8ef01SBarry Smith */
4a4963045SJacob Faibussowitsch #pragma once
56c7e564aSBarry Smith 
6e1bf4ed2SJacob Faibussowitsch #include <petscmacros.h>
7e1bf4ed2SJacob Faibussowitsch #include <petscsystypes.h>
8e1bf4ed2SJacob Faibussowitsch 
9af75f258SJacob Faibussowitsch #if defined(__cplusplus)
10af75f258SJacob Faibussowitsch   #include <exception> // std::exception
11af75f258SJacob Faibussowitsch #endif
12af75f258SJacob Faibussowitsch 
13ac09b921SBarry Smith /* SUBMANSEC = Sys */
14ac09b921SBarry Smith 
15edd03b47SJacob Faibussowitsch #define SETERRQ1(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
16edd03b47SJacob Faibussowitsch #define SETERRQ2(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
17edd03b47SJacob Faibussowitsch #define SETERRQ3(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
18edd03b47SJacob Faibussowitsch #define SETERRQ4(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
19edd03b47SJacob Faibussowitsch #define SETERRQ5(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
20edd03b47SJacob Faibussowitsch #define SETERRQ6(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
21edd03b47SJacob Faibussowitsch #define SETERRQ7(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
22edd03b47SJacob Faibussowitsch #define SETERRQ8(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
23edd03b47SJacob Faibussowitsch #define SETERRQ9(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
2498921bdaSJacob Faibussowitsch 
2530de9b25SBarry Smith /*MC
261957e957SBarry Smith    SETERRQ - Macro to be called when an error has been detected,
2730de9b25SBarry Smith 
2830de9b25SBarry Smith    Synopsis:
29aaa7dc30SBarry Smith    #include <petscsys.h>
3098921bdaSJacob Faibussowitsch    PetscErrorCode SETERRQ(MPI_Comm comm, PetscErrorCode ierr, char *message, ...)
3130de9b25SBarry Smith 
32d083f849SBarry Smith    Collective
3330de9b25SBarry Smith 
3430de9b25SBarry Smith    Input Parameters:
3595bd0b28SBarry Smith +  comm    - An MPI communicator, use `PETSC_COMM_SELF` unless you know all ranks of another communicator will detect the error
363af045c5SBarry Smith .  ierr    - nonzero error code, see the list of standard error codes in include/petscerror.h
3730de9b25SBarry Smith -  message - error message
3830de9b25SBarry Smith 
3930de9b25SBarry Smith   Level: beginner
4030de9b25SBarry Smith 
4130de9b25SBarry Smith    Notes:
4287497f52SBarry Smith    This is rarely needed, one should use `PetscCheck()` and `PetscCall()` and friends to automatically handle error conditions.
4330de9b25SBarry Smith    Once the error handler is called the calling function is then returned from with the given error code.
4430de9b25SBarry Smith 
4587497f52SBarry Smith    Experienced users can set the error handler with `PetscPushErrorHandler()`.
4630de9b25SBarry Smith 
4716a05f60SBarry Smith    Fortran Note:
48989c49a2SBarry Smith    `SETERRQ()` may be called from Fortran subroutines but `SETERRA()` must be called from the
493b1008b8SBarry Smith    Fortran main program.
503b1008b8SBarry Smith 
51db781477SPatrick Sanan .seealso: `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
52af27ebaaSBarry Smith           `PetscError()`, `PetscCall()`, `CHKMEMQ`, `CHKERRA()`, `PetscCallMPI()`, `PetscErrorCode`
5330de9b25SBarry Smith M*/
54e809461dSJacob Faibussowitsch #define SETERRQ(comm, ierr, ...) \
55e809461dSJacob Faibussowitsch   do { \
56e809461dSJacob Faibussowitsch     PetscErrorCode ierr_seterrq_petsc_ = PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \
57e809461dSJacob Faibussowitsch     return ierr_seterrq_petsc_ ? ierr_seterrq_petsc_ : PETSC_ERR_RETURN; \
58e809461dSJacob Faibussowitsch   } while (0)
59986eef2eSBarry Smith 
60ffc4695bSBarry Smith /*
61ffc4695bSBarry Smith     Returned from PETSc functions that are called from MPI, such as related to attributes
62ffc4695bSBarry Smith       Do not confuse PETSC_MPI_ERROR_CODE and PETSC_ERR_MPI, the first is registered with MPI and returned to MPI as
63888a1fb3SBarry Smith       an error code, the latter is a regular PETSc error code passed within PETSc code indicating an error was detected in an MPI call.
64ffc4695bSBarry Smith */
65ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CLASS;
66ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CODE;
67ffc4695bSBarry Smith 
68986eef2eSBarry Smith /*MC
69986eef2eSBarry Smith    SETERRMPI - Macro to be called when an error has been detected within an MPI callback function
70986eef2eSBarry Smith 
71989c49a2SBarry Smith    No Fortran Support
72989c49a2SBarry Smith 
73986eef2eSBarry Smith    Synopsis:
74986eef2eSBarry Smith    #include <petscsys.h>
7598921bdaSJacob Faibussowitsch    PetscErrorCode SETERRMPI(MPI_Comm comm, PetscErrorCode ierr, char *message, ...)
76986eef2eSBarry Smith 
77d083f849SBarry Smith    Collective
78986eef2eSBarry Smith 
79986eef2eSBarry Smith    Input Parameters:
8095bd0b28SBarry Smith +  comm    - An MPI communicator, use `PETSC_COMM_SELF` unless you know all ranks of another communicator will detect the error
81986eef2eSBarry Smith .  ierr    - nonzero error code, see the list of standard error codes in include/petscerror.h
82986eef2eSBarry Smith -  message - error message
83986eef2eSBarry Smith 
84986eef2eSBarry Smith   Level: developer
85986eef2eSBarry Smith 
8695bd0b28SBarry Smith    Note:
8787497f52SBarry Smith    This macro is FOR USE IN MPI CALLBACK FUNCTIONS ONLY, such as those passed to `MPI_Comm_create_keyval()`. It always returns the error code `PETSC_MPI_ERROR_CODE`
8887497f52SBarry Smith   which is registered with `MPI_Add_error_code()` when PETSc is initialized.
89986eef2eSBarry Smith 
90af27ebaaSBarry Smith .seealso: `SETERRQ()`, `PetscCall()`, `PetscCallMPI()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `PetscErrorCode`
91986eef2eSBarry Smith M*/
923ba16761SJacob Faibussowitsch #define SETERRMPI(comm, ierr, ...) return ((void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__), PETSC_MPI_ERROR_CODE)
93ee8199e6SMatthew G. Knepley 
94ee8199e6SMatthew G. Knepley /*MC
95f388eb8bSPatrick Sanan    SETERRA - Fortran-only macro that can be called when an error has been detected from the main program
96f388eb8bSPatrick Sanan 
97f388eb8bSPatrick Sanan    Synopsis:
98f388eb8bSPatrick Sanan    #include <petscsys.h>
99f388eb8bSPatrick Sanan    PetscErrorCode SETERRA(MPI_Comm comm, PetscErrorCode ierr, char *message)
100f388eb8bSPatrick Sanan 
101f388eb8bSPatrick Sanan    Collective
102f388eb8bSPatrick Sanan 
103f388eb8bSPatrick Sanan    Input Parameters:
10495bd0b28SBarry Smith +  comm    - An MPI communicator, so that the error can be collective
105f388eb8bSPatrick Sanan .  ierr    - nonzero error code, see the list of standard error codes in include/petscerror.h
106f388eb8bSPatrick Sanan -  message - error message in the printf format
107f388eb8bSPatrick Sanan 
108f388eb8bSPatrick Sanan    Level: beginner
109f388eb8bSPatrick Sanan 
110f388eb8bSPatrick Sanan    Notes:
11187497f52SBarry Smith    This should only be used with Fortran. With C/C++, use `SETERRQ()`.
112f388eb8bSPatrick Sanan 
11387497f52SBarry Smith    `SETERRQ()` may be called from Fortran subroutines but `SETERRA()` must be called from the
114f388eb8bSPatrick Sanan     Fortran main program.
115f388eb8bSPatrick Sanan 
116af27ebaaSBarry Smith .seealso: `SETERRQ()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()`, `PetscErrorCode`
117f388eb8bSPatrick Sanan M*/
118f388eb8bSPatrick Sanan 
119f388eb8bSPatrick Sanan /*MC
120fa190f98SMatthew G. Knepley    SETERRABORT - Macro that can be called when an error has been detected,
121fa190f98SMatthew G. Knepley 
122fa190f98SMatthew G. Knepley    Synopsis:
123fa190f98SMatthew G. Knepley    #include <petscsys.h>
12498921bdaSJacob Faibussowitsch    PetscErrorCode SETERRABORT(MPI_Comm comm, PetscErrorCode ierr, char *message, ...)
125fa190f98SMatthew G. Knepley 
126d083f849SBarry Smith    Collective
127fa190f98SMatthew G. Knepley 
128fa190f98SMatthew G. Knepley    Input Parameters:
12995bd0b28SBarry Smith +  comm    - An MPI communicator, so that the error can be collective
1303af045c5SBarry Smith .  ierr    - nonzero error code, see the list of standard error codes in include/petscerror.h
131fa190f98SMatthew G. Knepley -  message - error message in the printf format
132fa190f98SMatthew G. Knepley 
133fa190f98SMatthew G. Knepley    Level: beginner
134fa190f98SMatthew G. Knepley 
135fa190f98SMatthew G. Knepley    Notes:
13687497f52SBarry Smith    This function just calls `MPI_Abort()`.
13787497f52SBarry Smith 
13887497f52SBarry Smith    This should only be called in routines that cannot return an error code, such as in C++ constructors.
139fa190f98SMatthew G. Knepley 
140989c49a2SBarry Smith    Fortran Note:
141989c49a2SBarry Smith    Use `SETERRA()` in Fortran main program and `SETERRQ()` in Fortran subroutines
142989c49a2SBarry Smith 
143989c49a2SBarry Smith    Developer Note:
144989c49a2SBarry Smith    In Fortran `SETERRA()` could be called `SETERRABORT()` since they serve the same purpose
145989c49a2SBarry Smith 
146af27ebaaSBarry Smith .seealso: `SETERRQ()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `PetscCall()`, `CHKMEMQ`, `PetscErrorCode`
147fa190f98SMatthew G. Knepley M*/
1489371c9d4SSatish Balay #define SETERRABORT(comm, ierr, ...) \
1499371c9d4SSatish Balay   do { \
1503ba16761SJacob Faibussowitsch     (void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \
15198921bdaSJacob Faibussowitsch     MPI_Abort(comm, ierr); \
15298921bdaSJacob Faibussowitsch   } while (0)
1539a00fa46SSatish Balay 
15430de9b25SBarry Smith /*MC
155*139c8099SSatish Balay   PetscCheck - Checks that a particular condition is true; if not true, then returns the provided error code
1562c71b3e2SJacob Faibussowitsch 
1572c71b3e2SJacob Faibussowitsch   Synopsis:
1582c71b3e2SJacob Faibussowitsch   #include <petscerror.h>
1592c71b3e2SJacob Faibussowitsch   void PetscCheck(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
1602c71b3e2SJacob Faibussowitsch 
1617cdbe19fSJose E. Roman   Collective; No Fortran Support
1622c71b3e2SJacob Faibussowitsch 
1632c71b3e2SJacob Faibussowitsch   Input Parameters:
1642c71b3e2SJacob Faibussowitsch + cond    - The boolean condition
1652c71b3e2SJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
1662c71b3e2SJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
1672c71b3e2SJacob Faibussowitsch - message - Error message in printf format
1682c71b3e2SJacob Faibussowitsch 
169667f096bSBarry Smith   Level: beginner
170667f096bSBarry Smith 
1712c71b3e2SJacob Faibussowitsch   Notes:
1722c71b3e2SJacob Faibussowitsch   Enabled in both optimized and debug builds.
1732c71b3e2SJacob Faibussowitsch 
174a2454668SBarry Smith   As a general rule, `PetscCheck()` is used to check "usage error" (for example, passing an incorrect value as a function argument),
175a2454668SBarry Smith   `PetscAssert()` is used to "check for bugs in PETSc" (for example, is a value in a PETSc data structure nonsensical).
176a2454668SBarry Smith   However, for functions that are called in a "hot spot", for example, thousands of times in a loop, `PetscAssert()` should be used instead
177a2454668SBarry Smith   of `PetscCheck()` since the former is compiled out in PETSc's optimization code.
178a2454668SBarry Smith 
17987497f52SBarry Smith   Calls `SETERRQ()` if the assertion fails, so can only be called from functions returning a
18087497f52SBarry Smith   `PetscErrorCode` (or equivalent type after conversion).
1812c71b3e2SJacob Faibussowitsch 
182af27ebaaSBarry Smith .seealso: `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheckAbort()`, `PetscErrorCode`
1839566063dSJacob Faibussowitsch M*/
1849371c9d4SSatish Balay #define PetscCheck(cond, comm, ierr, ...) \
1859371c9d4SSatish Balay   do { \
1869371c9d4SSatish Balay     if (PetscUnlikely(!(cond))) SETERRQ(comm, ierr, __VA_ARGS__); \
1879371c9d4SSatish Balay   } while (0)
1882c71b3e2SJacob Faibussowitsch 
1892c71b3e2SJacob Faibussowitsch /*MC
1904be741a6SBarry Smith   PetscCheckAbort - Check that a particular condition is true, otherwise prints error and aborts
1914be741a6SBarry Smith 
1924be741a6SBarry Smith   Synopsis:
1934be741a6SBarry Smith   #include <petscerror.h>
1944be741a6SBarry Smith   void PetscCheckAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
1954be741a6SBarry Smith 
1967cdbe19fSJose E. Roman   Collective; No Fortran Support
1974be741a6SBarry Smith 
1984be741a6SBarry Smith   Input Parameters:
1994be741a6SBarry Smith + cond    - The boolean condition
2004be741a6SBarry Smith . comm    - The communicator on which the check can be collective on
2014be741a6SBarry Smith . ierr    - A nonzero error code, see include/petscerror.h for the complete list
2024be741a6SBarry Smith - message - Error message in printf format
2034be741a6SBarry Smith 
204667f096bSBarry Smith   Level: developer
205667f096bSBarry Smith 
2064be741a6SBarry Smith   Notes:
2074be741a6SBarry Smith   Enabled in both optimized and debug builds.
2084be741a6SBarry Smith 
20987497f52SBarry Smith   Calls `SETERRABORT()` if the assertion fails, can be called from a function that does not return an
21087497f52SBarry Smith   error code, such as a C++ constructor. usually `PetscCheck()` should be used.
2114be741a6SBarry Smith 
212af27ebaaSBarry Smith .seealso: `PetscAssertAbort()`, `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheck()`, `SETERRABORT()`, `PetscErrorCode`
2134be741a6SBarry Smith M*/
2149371c9d4SSatish Balay #define PetscCheckAbort(cond, comm, ierr, ...) \
2150e6b6b59SJacob Faibussowitsch   do { \
2160e6b6b59SJacob Faibussowitsch     if (PetscUnlikely(!(cond))) SETERRABORT(comm, ierr, __VA_ARGS__); \
217c7481402SJacob Faibussowitsch   } while (0)
2184be741a6SBarry Smith 
2194be741a6SBarry Smith /*MC
2209ace16cdSJacob Faibussowitsch   PetscAssert - Assert that a particular condition is true
2219ace16cdSJacob Faibussowitsch 
2229ace16cdSJacob Faibussowitsch   Synopsis:
2239ace16cdSJacob Faibussowitsch   #include <petscerror.h>
2249ace16cdSJacob Faibussowitsch   void PetscAssert(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
2259ace16cdSJacob Faibussowitsch 
2267cdbe19fSJose E. Roman   Collective; No Fortran Support
2279ace16cdSJacob Faibussowitsch 
2289ace16cdSJacob Faibussowitsch   Input Parameters:
2299ace16cdSJacob Faibussowitsch + cond    - The boolean condition
2309ace16cdSJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
2319ace16cdSJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
232af27ebaaSBarry Smith - message - Error message in `printf()` format
2339ace16cdSJacob Faibussowitsch 
234667f096bSBarry Smith   Level: beginner
235667f096bSBarry Smith 
2369ace16cdSJacob Faibussowitsch   Notes:
237c7481402SJacob Faibussowitsch   Equivalent to `PetscCheck()` if debugging is enabled, and `PetscAssume(cond)` otherwise.
2382c71b3e2SJacob Faibussowitsch 
23987497f52SBarry Smith   See `PetscCheck()` for usage and behaviour.
24087497f52SBarry Smith 
24187497f52SBarry Smith   This is needed instead of simply using `assert()` because this correctly handles the collective nature of errors under MPI
2429ace16cdSJacob Faibussowitsch 
243af27ebaaSBarry Smith .seealso: `PetscCheck()`, `SETERRQ()`, `PetscError()`, `PetscAssertAbort()`, `PetscErrorCode`
2449566063dSJacob Faibussowitsch M*/
245c7481402SJacob Faibussowitsch #if PetscDefined(USE_DEBUG)
246c7481402SJacob Faibussowitsch   #define PetscAssert(cond, comm, ierr, ...) PetscCheck(cond, comm, ierr, __VA_ARGS__)
247c7481402SJacob Faibussowitsch #else
248c7481402SJacob Faibussowitsch   #define PetscAssert(cond, ...) PetscAssume(cond)
249c7481402SJacob Faibussowitsch #endif
2509ace16cdSJacob Faibussowitsch 
2519ace16cdSJacob Faibussowitsch /*MC
2520e6b6b59SJacob Faibussowitsch   PetscAssertAbort - Assert that a particular condition is true, otherwise prints error and aborts
2530e6b6b59SJacob Faibussowitsch 
2540e6b6b59SJacob Faibussowitsch   Synopsis:
2550e6b6b59SJacob Faibussowitsch   #include <petscerror.h>
2560e6b6b59SJacob Faibussowitsch   void PetscAssertAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
2570e6b6b59SJacob Faibussowitsch 
2587cdbe19fSJose E. Roman   Collective; No Fortran Support
2590e6b6b59SJacob Faibussowitsch 
2600e6b6b59SJacob Faibussowitsch   Input Parameters:
2610e6b6b59SJacob Faibussowitsch + cond    - The boolean condition
2620e6b6b59SJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
2630e6b6b59SJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
2640e6b6b59SJacob Faibussowitsch - message - Error message in printf format
2650e6b6b59SJacob Faibussowitsch 
266667f096bSBarry Smith   Level: beginner
267667f096bSBarry Smith 
26895bd0b28SBarry Smith   Note:
2690e6b6b59SJacob Faibussowitsch   Enabled only in debug builds. See `PetscCheckAbort()` for usage.
2700e6b6b59SJacob Faibussowitsch 
2710e6b6b59SJacob Faibussowitsch .seealso: `PetscCheckAbort()`, `PetscAssert()`, `PetscCheck()`, `SETERRABORT()`, `PetscError()`
2720e6b6b59SJacob Faibussowitsch M*/
2733ba16761SJacob Faibussowitsch #if PetscDefined(USE_DEBUG)
2743ba16761SJacob Faibussowitsch   #define PetscAssertAbort(cond, comm, ierr, ...) PetscCheckAbort(cond, comm, ierr, __VA_ARGS__)
2753ba16761SJacob Faibussowitsch #else
2763ba16761SJacob Faibussowitsch   #define PetscAssertAbort(cond, comm, ierr, ...) PetscAssume(cond)
2773ba16761SJacob Faibussowitsch #endif
2780e6b6b59SJacob Faibussowitsch 
2790e6b6b59SJacob Faibussowitsch /*MC
2803ba16761SJacob Faibussowitsch   PetscCall - Calls a PETSc function and then checks the resulting error code, if it is
2813ba16761SJacob Faibussowitsch   non-zero it calls the error handler and returns from the current function with the error
2823ba16761SJacob Faibussowitsch   code.
2839566063dSJacob Faibussowitsch 
2849566063dSJacob Faibussowitsch   Synopsis:
2859566063dSJacob Faibussowitsch   #include <petscerror.h>
28649c86fc7SBarry Smith   void PetscCall(PetscFunction(args))
2879566063dSJacob Faibussowitsch 
2889566063dSJacob Faibussowitsch   Not Collective
2899566063dSJacob Faibussowitsch 
2909566063dSJacob Faibussowitsch   Input Parameter:
29149c86fc7SBarry Smith . PetscFunction - any PETSc function that returns an error code
2929566063dSJacob Faibussowitsch 
293667f096bSBarry Smith   Level: beginner
294667f096bSBarry Smith 
2959566063dSJacob Faibussowitsch   Notes:
2969566063dSJacob Faibussowitsch   Once the error handler is called the calling function is then returned from with the given
29787497f52SBarry Smith   error code. Experienced users can set the error handler with `PetscPushErrorHandler()`.
2989566063dSJacob Faibussowitsch 
29987497f52SBarry Smith   `PetscCall()` cannot be used in functions returning a datatype not convertible to
3008af9f190SBarry Smith   `PetscErrorCode`. For example, `PetscCall()` may not be used in functions returning `void`, use
3013ba16761SJacob Faibussowitsch   `PetscCallAbort()` or `PetscCallVoid()` in this case.
3029566063dSJacob Faibussowitsch 
3039566063dSJacob Faibussowitsch   Example Usage:
3049566063dSJacob Faibussowitsch .vb
3059566063dSJacob Faibussowitsch   PetscCall(PetscInitiailize(...)); // OK to call even when PETSc is not yet initialized!
3069566063dSJacob Faibussowitsch 
3079566063dSJacob Faibussowitsch   struct my_struct
3089566063dSJacob Faibussowitsch   {
3099566063dSJacob Faibussowitsch     void *data;
3109566063dSJacob Faibussowitsch   } my_complex_type;
3119566063dSJacob Faibussowitsch 
3129566063dSJacob Faibussowitsch   struct my_struct bar(void)
3139566063dSJacob Faibussowitsch   {
3146aad120cSJose E. Roman     PetscCall(foo(15)); // ERROR PetscErrorCode not convertible to struct my_struct!
3159566063dSJacob Faibussowitsch   }
3169566063dSJacob Faibussowitsch 
3176aad120cSJose E. Roman   PetscCall(bar()) // ERROR input not convertible to PetscErrorCode
3189566063dSJacob Faibussowitsch .ve
3199566063dSJacob Faibussowitsch 
32087497f52SBarry Smith   It is also possible to call this directly on a `PetscErrorCode` variable
32149c86fc7SBarry Smith .vb
32249c86fc7SBarry Smith   PetscCall(ierr);  // check if ierr is nonzero
32349c86fc7SBarry Smith .ve
32449c86fc7SBarry Smith 
325792fecdfSBarry Smith   Should not be used to call callback functions provided by users, `PetscCallBack()` should be used in that situation.
326ef1023bdSBarry Smith 
3276a8be23eSBarry Smith   `PetscUseTypeMethod()` or `PetscTryTypeMethod()` should be used when calling functions pointers contained in a PETSc object's `ops` array
3286a8be23eSBarry Smith 
32949c86fc7SBarry Smith   Fortran Notes:
33098d50953SPierre Jolivet     The Fortran function in which this is used must declare a `PetscErrorCode` variable necessarily named `ierr`, and `ierr` must be
33187497f52SBarry Smith     the final argument to the PETSc function being called.
33249c86fc7SBarry Smith 
33398d50953SPierre Jolivet     In the main program and in Fortran subroutines that do not have `ierr` as the final return parameter, one
33487497f52SBarry Smith     should use `PetscCallA()`
33549c86fc7SBarry Smith 
33649c86fc7SBarry Smith   Example Fortran Usage:
33749c86fc7SBarry Smith .vb
33849c86fc7SBarry Smith   PetscErrorCode ierr
33949c86fc7SBarry Smith   Vec v
34049c86fc7SBarry Smith 
34149c86fc7SBarry Smith   ...
34249c86fc7SBarry Smith   PetscCall(VecShift(v, 1.0, ierr))
34349c86fc7SBarry Smith   PetscCallA(VecShift(v, 1.0, ierr))
34449c86fc7SBarry Smith .ve
34549c86fc7SBarry Smith 
346667f096bSBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`,
347667f096bSBarry Smith           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`,
3483ba16761SJacob Faibussowitsch           `CHKERRMPI()`, `PetscCallBack()`, `PetscCallAbort()`, `PetscCallVoid()`
3499566063dSJacob Faibussowitsch M*/
350ef1023bdSBarry Smith 
351ef1023bdSBarry Smith /*MC
35298d50953SPierre Jolivet    PetscCallA - Fortran-only macro that should be used in the main program and subroutines that do not have `ierr` as the final return parameter, to call PETSc functions instead of using
35398d50953SPierre Jolivet    `PetscCall()` which should be used in other Fortran subroutines
354989c49a2SBarry Smith 
355989c49a2SBarry Smith    Synopsis:
356989c49a2SBarry Smith    #include <petscsys.h>
357989c49a2SBarry Smith    PetscErrorCode PetscCallA(PetscFunction(arguments, ierr))
358989c49a2SBarry Smith 
359989c49a2SBarry Smith    Collective
360989c49a2SBarry Smith 
361989c49a2SBarry Smith    Input Parameter:
362989c49a2SBarry Smith .  PetscFunction(arguments,ierr) - the call to the function
363989c49a2SBarry Smith 
364989c49a2SBarry Smith   Level: beginner
365989c49a2SBarry Smith 
366989c49a2SBarry Smith    Notes:
367989c49a2SBarry Smith    This should only be used with Fortran. With C/C++, use `PetscCall()` always.
368989c49a2SBarry Smith 
36998d50953SPierre Jolivet    The Fortran function in which this is used must declare a `PetscErrorCode` variable necessarily named `ierr`
370989c49a2SBarry Smith    Use `SETERRA()` to set an error in a Fortran main program and `SETERRQ()` in Fortran subroutines
371989c49a2SBarry Smith 
372989c49a2SBarry Smith .seealso: `SETERRQ()`, `SETERRA()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()`
373989c49a2SBarry Smith M*/
374989c49a2SBarry Smith 
375989c49a2SBarry Smith /*MC
376792fecdfSBarry Smith   PetscCallBack - Calls a user provided PETSc callback function and then checks the resulting error code, if it is non-zero it calls the error
377ef1023bdSBarry Smith   handler and returns from the current function with the error code.
378ef1023bdSBarry Smith 
379ef1023bdSBarry Smith   Synopsis:
380ef1023bdSBarry Smith   #include <petscerror.h>
381792fecdfSBarry Smith   void PetscCallBack(const char *functionname, PetscFunction(args))
382ef1023bdSBarry Smith 
3837cdbe19fSJose E. Roman   Not Collective; No Fortran Support
384ef1023bdSBarry Smith 
385ef1023bdSBarry Smith   Input Parameters:
386ef1023bdSBarry Smith + functionname - the name of the function being called, this can be a string with spaces that describes the meaning of the callback
38787497f52SBarry Smith - PetscFunction - user provided callback function that returns an error code
388ef1023bdSBarry Smith 
389ef1023bdSBarry Smith   Example Usage:
390ef1023bdSBarry Smith .vb
391792fecdfSBarry Smith   PetscCallBack("XXX callback to do something", a->callback(...));
392ef1023bdSBarry Smith .ve
393ef1023bdSBarry Smith 
394ef1023bdSBarry Smith   Level: developer
395ef1023bdSBarry Smith 
396667f096bSBarry Smith   Notes:
397667f096bSBarry Smith   Once the error handler is called the calling function is then returned from with the given
398667f096bSBarry Smith   error code. Experienced users can set the error handler with `PetscPushErrorHandler()`.
399667f096bSBarry Smith 
400667f096bSBarry Smith   `PetscCallBack()` should only be called in PETSc when a call is being made to a user provided call-back routine.
401667f096bSBarry Smith 
40287497f52SBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscCall()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`
403ef1023bdSBarry Smith           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, `CHKERRMPI()`, `PetscCall()`
404ef1023bdSBarry Smith M*/
405ef1023bdSBarry Smith 
4063ba16761SJacob Faibussowitsch /*MC
4078af9f190SBarry Smith   PetscCallVoid - Like `PetscCall()` but for use in functions that return `void`
4083ba16761SJacob Faibussowitsch 
4093ba16761SJacob Faibussowitsch   Synopsis:
4103ba16761SJacob Faibussowitsch   #include <petscerror.h>
4118af9f190SBarry Smith   void PetscCallVoid(PetscFunction(args))
4123ba16761SJacob Faibussowitsch 
4137cdbe19fSJose E. Roman   Not Collective; No Fortran Support
4143ba16761SJacob Faibussowitsch 
4153ba16761SJacob Faibussowitsch   Input Parameter:
4163ba16761SJacob Faibussowitsch . PetscFunction - any PETSc function that returns an error code
4173ba16761SJacob Faibussowitsch 
4183ba16761SJacob Faibussowitsch   Example Usage:
4193ba16761SJacob Faibussowitsch .vb
4203ba16761SJacob Faibussowitsch   void foo()
4213ba16761SJacob Faibussowitsch   {
4223ba16761SJacob Faibussowitsch     KSP ksp;
4233ba16761SJacob Faibussowitsch 
4243ba16761SJacob Faibussowitsch     PetscFunctionBeginUser;
4253ba16761SJacob Faibussowitsch     // OK, properly handles PETSc error codes
4263ba16761SJacob Faibussowitsch     PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp));
4278af9f190SBarry Smith     PetscFunctionReturnVoid();
4283ba16761SJacob Faibussowitsch   }
4293ba16761SJacob Faibussowitsch 
4303ba16761SJacob Faibussowitsch   PetscErrorCode bar()
4313ba16761SJacob Faibussowitsch   {
4323ba16761SJacob Faibussowitsch     KSP ksp;
4333ba16761SJacob Faibussowitsch 
4343ba16761SJacob Faibussowitsch     PetscFunctionBeginUser;
4353ba16761SJacob Faibussowitsch     // ERROR, Non-void function 'bar' should return a value
4363ba16761SJacob Faibussowitsch     PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp));
4373ba16761SJacob Faibussowitsch     // OK, returning PetscErrorCode
4383ba16761SJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
4393ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
4403ba16761SJacob Faibussowitsch   }
441667f096bSBarry Smith .ve
4423ba16761SJacob Faibussowitsch 
4433ba16761SJacob Faibussowitsch   Level: beginner
4443ba16761SJacob Faibussowitsch 
445667f096bSBarry Smith   Notes:
446667f096bSBarry Smith   Has identical usage to `PetscCall()`, except that it returns `void` on error instead of a
447667f096bSBarry Smith   `PetscErrorCode`. See `PetscCall()` for more detailed discussion.
448667f096bSBarry Smith 
449667f096bSBarry Smith   Note that users should prefer `PetscCallAbort()` to this routine. While this routine does
450667f096bSBarry Smith   "handle" errors by returning from the enclosing function, it effectively gobbles the
451667f096bSBarry Smith   error. Since the enclosing function itself returns `void`, its callers have no way of knowing
452667f096bSBarry Smith   that the routine returned early due to an error. `PetscCallAbort()` at least ensures that the
453667f096bSBarry Smith   program crashes gracefully.
454667f096bSBarry Smith 
4558af9f190SBarry Smith .seealso: `PetscCall()`, `PetscErrorCode`, `PetscCallAbort()`
4563ba16761SJacob Faibussowitsch M*/
4579566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
4589566063dSJacob Faibussowitsch void PetscCall(PetscErrorCode);
459792fecdfSBarry Smith void PetscCallBack(const char *, PetscErrorCode);
4609566063dSJacob Faibussowitsch void PetscCallVoid(PetscErrorCode);
4619566063dSJacob Faibussowitsch #else
4629371c9d4SSatish Balay   #define PetscCall(...) \
4639371c9d4SSatish Balay     do { \
4643ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_q_; \
46537154d10SBarry Smith       PetscStackUpdateLine; \
4663ba16761SJacob Faibussowitsch       ierr_petsc_call_q_ = __VA_ARGS__; \
4673ba16761SJacob 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, " "); \
4689566063dSJacob Faibussowitsch     } while (0)
4699371c9d4SSatish Balay   #define PetscCallBack(function, ...) \
4709371c9d4SSatish Balay     do { \
4713ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_q_; \
472e33ced7fSLisandro Dalcin       PetscStackUpdateLine; \
473e33ced7fSLisandro Dalcin       PetscStackPushExternal(function); \
4743ba16761SJacob Faibussowitsch       ierr_petsc_call_q_ = __VA_ARGS__; \
475ef1023bdSBarry Smith       PetscStackPop; \
4763ba16761SJacob 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, " "); \
477ef1023bdSBarry Smith     } while (0)
4789371c9d4SSatish Balay   #define PetscCallVoid(...) \
4799371c9d4SSatish Balay     do { \
4803ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_void_; \
481e33ced7fSLisandro Dalcin       PetscStackUpdateLine; \
4823ba16761SJacob Faibussowitsch       ierr_petsc_call_void_ = __VA_ARGS__; \
4833ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_void_ != PETSC_SUCCESS)) { \
4843ba16761SJacob Faibussowitsch         ierr_petsc_call_void_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_void_, PETSC_ERROR_REPEAT, " "); \
4853ba16761SJacob Faibussowitsch         (void)ierr_petsc_call_void_; \
4869371c9d4SSatish Balay         return; \
4879371c9d4SSatish Balay       } \
4889566063dSJacob Faibussowitsch     } while (0)
4899566063dSJacob Faibussowitsch #endif
4909566063dSJacob Faibussowitsch 
4919566063dSJacob Faibussowitsch /*MC
4929566063dSJacob Faibussowitsch   CHKERRQ - Checks error code returned from PETSc function
49330de9b25SBarry Smith 
49430de9b25SBarry Smith   Synopsis:
495aaa7dc30SBarry Smith   #include <petscsys.h>
4969566063dSJacob Faibussowitsch   void CHKERRQ(PetscErrorCode ierr)
4979566063dSJacob Faibussowitsch 
4989566063dSJacob Faibussowitsch   Not Collective
4999566063dSJacob Faibussowitsch 
5002fe279fdSBarry Smith   Input Parameter:
5019566063dSJacob Faibussowitsch . ierr - nonzero error code
5029566063dSJacob Faibussowitsch 
5039566063dSJacob Faibussowitsch   Level: deprecated
5049566063dSJacob Faibussowitsch 
505667f096bSBarry Smith   Note:
506667f096bSBarry Smith   Deprecated in favor of `PetscCall()`. This routine behaves identically to it.
507667f096bSBarry Smith 
508db781477SPatrick Sanan .seealso: `PetscCall()`
5099566063dSJacob Faibussowitsch M*/
5109566063dSJacob Faibussowitsch #define CHKERRQ(...) PetscCall(__VA_ARGS__)
5119566063dSJacob Faibussowitsch #define CHKERRV(...) PetscCallVoid(__VA_ARGS__)
5129566063dSJacob Faibussowitsch 
513db9cea48SBarry Smith PETSC_EXTERN void PetscMPIErrorString(PetscMPIInt, char *);
514db9cea48SBarry Smith 
5159566063dSJacob Faibussowitsch /*MC
5169566063dSJacob Faibussowitsch   PetscCallMPI - Checks error code returned from MPI calls, if non-zero it calls the error
5179566063dSJacob Faibussowitsch   handler and then returns
5189566063dSJacob Faibussowitsch 
5199566063dSJacob Faibussowitsch   Synopsis:
5209566063dSJacob Faibussowitsch   #include <petscerror.h>
52149c86fc7SBarry Smith   void PetscCallMPI(MPI_Function(args))
52230de9b25SBarry Smith 
523eca87e8dSBarry Smith   Not Collective
52430de9b25SBarry Smith 
5252fe279fdSBarry Smith   Input Parameter:
52649c86fc7SBarry Smith . MPI_Function - an MPI function that returns an MPI error code
52730de9b25SBarry Smith 
528667f096bSBarry Smith   Level: beginner
529667f096bSBarry Smith 
5309566063dSJacob Faibussowitsch   Notes:
53187497f52SBarry Smith   Always returns the error code `PETSC_ERR_MPI`; the MPI error code and string are embedded in
5329566063dSJacob Faibussowitsch   the string error message. Do not use this to call any other routines (for example PETSc
5333ba16761SJacob Faibussowitsch   routines), it should only be used for direct MPI calls. The user may configure PETSc with the
5343ba16761SJacob Faibussowitsch   `--with-strict-petscerrorcode` option to check this at compile-time, otherwise they must
5359566063dSJacob Faibussowitsch   check this themselves.
5369566063dSJacob Faibussowitsch 
537aaa8cc7dSPierre Jolivet   This routine can only be used in functions returning `PetscErrorCode` themselves. If the
5383ba16761SJacob Faibussowitsch   calling function returns a different type, use `PetscCallMPIAbort()` instead.
5393ba16761SJacob Faibussowitsch 
5409566063dSJacob Faibussowitsch   Example Usage:
5419566063dSJacob Faibussowitsch .vb
5429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(...)); // OK, calling MPI function
5439566063dSJacob Faibussowitsch 
5449566063dSJacob Faibussowitsch   PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead!
5459566063dSJacob Faibussowitsch .ve
5469566063dSJacob Faibussowitsch 
54749c86fc7SBarry Smith   Fortran Notes:
54887497f52SBarry Smith     The Fortran function from which this is used must declare a variable `PetscErrorCode` ierr and ierr must be
54949c86fc7SBarry Smith     the final argument to the MPI function being called.
55049c86fc7SBarry Smith 
55149c86fc7SBarry Smith     In the main program and in Fortran subroutines that do not have ierr as the final return parameter one
55287497f52SBarry Smith     should use `PetscCallMPIA()`
55349c86fc7SBarry Smith 
55449c86fc7SBarry Smith   Fortran Usage:
55549c86fc7SBarry Smith .vb
55649c86fc7SBarry Smith   PetscErrorCode ierr or integer ierr
55749c86fc7SBarry Smith   ...
55849c86fc7SBarry Smith   PetscCallMPI(MPI_Comm_size(...,ierr))
55949c86fc7SBarry Smith   PetscCallMPIA(MPI_Comm_size(...,ierr)) ! Will abort after calling error handler
56049c86fc7SBarry Smith 
56149c86fc7SBarry Smith   PetscCallMPI(MPI_Comm_size(...,eflag)) ! ERROR, final argument must be ierr
56249c86fc7SBarry Smith .ve
56349c86fc7SBarry Smith 
564db781477SPatrick Sanan .seealso: `SETERRMPI()`, `PetscCall()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`,
5653ba16761SJacob Faibussowitsch           `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
5663ba16761SJacob Faibussowitsch           `PetscError()`, `CHKMEMQ`
5673ba16761SJacob Faibussowitsch M*/
5683ba16761SJacob Faibussowitsch 
5693ba16761SJacob Faibussowitsch /*MC
5703ba16761SJacob Faibussowitsch   PetscCallMPIAbort - Like `PetscCallMPI()` but calls `MPI_Abort()` on error
5713ba16761SJacob Faibussowitsch 
5723ba16761SJacob Faibussowitsch   Synopsis:
5733ba16761SJacob Faibussowitsch   #include <petscerror.h>
5743ba16761SJacob Faibussowitsch   void PetscCallMPIAbort(MPI_Comm comm, MPI_Function(args))
5753ba16761SJacob Faibussowitsch 
5763ba16761SJacob Faibussowitsch   Not Collective
5773ba16761SJacob Faibussowitsch 
5783ba16761SJacob Faibussowitsch   Input Parameters:
5793ba16761SJacob Faibussowitsch + comm         - the MPI communicator to abort on
5803ba16761SJacob Faibussowitsch - MPI_Function - an MPI function that returns an MPI error code
5813ba16761SJacob Faibussowitsch 
582667f096bSBarry Smith   Level: beginner
583667f096bSBarry Smith 
5843ba16761SJacob Faibussowitsch   Notes:
5853ba16761SJacob Faibussowitsch   Usage is identical to `PetscCallMPI()`. See `PetscCallMPI()` for detailed discussion.
5863ba16761SJacob Faibussowitsch 
5873ba16761SJacob Faibussowitsch   This routine may be used in functions returning `void` or other non-`PetscErrorCode` types.
5883ba16761SJacob Faibussowitsch 
589989c49a2SBarry Smith   Fortran Note:
590989c49a2SBarry Smith   In Fortran this is called `PetscCallMPIA()` and is intended to be used in the main program while `PetscCallMPI()` is
591989c49a2SBarry Smith   used in Fortran subroutines.
592989c49a2SBarry Smith 
593989c49a2SBarry Smith   Developer Note:
594989c49a2SBarry Smith   This should have the same name in Fortran.
595989c49a2SBarry Smith 
5963ba16761SJacob Faibussowitsch .seealso: `PetscCallMPI()`, `PetscCallAbort()`, `SETERRABORT()`
59730de9b25SBarry Smith M*/
5983fcd9f07SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
59963a3b9bcSJacob Faibussowitsch void PetscCallMPI(PetscMPIInt);
6003ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm, PetscMPIInt);
601064a246eSJacob Faibussowitsch #else
6023ba16761SJacob Faibussowitsch   #define PetscCallMPI_Private(__PETSC_STACK_POP_FUNC__, __SETERR_FUNC__, __COMM__, ...) \
6039371c9d4SSatish Balay     do { \
6043ba16761SJacob Faibussowitsch       PetscMPIInt ierr_petsc_call_mpi_; \
605ef1023bdSBarry Smith       PetscStackUpdateLine; \
606792fecdfSBarry Smith       PetscStackPushExternal("MPI function"); \
607d71ae5a4SJacob Faibussowitsch       { \
6083ba16761SJacob Faibussowitsch         ierr_petsc_call_mpi_ = __VA_ARGS__; \
609d71ae5a4SJacob Faibussowitsch       } \
6103ba16761SJacob Faibussowitsch       __PETSC_STACK_POP_FUNC__; \
6113ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_mpi_ != MPI_SUCCESS)) { \
6123ba16761SJacob Faibussowitsch         char petsc_mpi_7_errorstring[2 * MPI_MAX_ERROR_STRING]; \
6133ba16761SJacob Faibussowitsch         PetscMPIErrorString(ierr_petsc_call_mpi_, (char *)petsc_mpi_7_errorstring); \
6143ba16761SJacob Faibussowitsch         __SETERR_FUNC__(__COMM__, PETSC_ERR_MPI, "MPI error %d %s", (int)ierr_petsc_call_mpi_, petsc_mpi_7_errorstring); \
6153fcd9f07SJacob Faibussowitsch       } \
6163fcd9f07SJacob Faibussowitsch     } while (0)
6173ba16761SJacob Faibussowitsch 
6183ba16761SJacob Faibussowitsch   #define PetscCallMPI(...)            PetscCallMPI_Private(PetscStackPop, SETERRQ, PETSC_COMM_SELF, __VA_ARGS__)
6193ba16761SJacob Faibussowitsch   #define PetscCallMPIAbort(comm, ...) PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRABORT, comm, __VA_ARGS__)
6209566063dSJacob Faibussowitsch #endif
621064a246eSJacob Faibussowitsch 
6227037b0edSPatrick Sanan /*MC
6239566063dSJacob Faibussowitsch   CHKERRMPI - Checks error code returned from MPI calls, if non-zero it calls the error
6249566063dSJacob Faibussowitsch   handler and then returns
6259566063dSJacob Faibussowitsch 
6269566063dSJacob Faibussowitsch   Synopsis:
6279566063dSJacob Faibussowitsch   #include <petscerror.h>
6289566063dSJacob Faibussowitsch   void CHKERRMPI(PetscErrorCode ierr)
6299566063dSJacob Faibussowitsch 
6309566063dSJacob Faibussowitsch   Not Collective
6319566063dSJacob Faibussowitsch 
6329566063dSJacob Faibussowitsch   Input Parameter:
6339566063dSJacob Faibussowitsch . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
6349566063dSJacob Faibussowitsch 
6359566063dSJacob Faibussowitsch   Level: deprecated
6369566063dSJacob Faibussowitsch 
637667f096bSBarry Smith   Note:
638667f096bSBarry Smith   Deprecated in favor of `PetscCallMPI()`. This routine behaves identically to it.
639667f096bSBarry Smith 
640db781477SPatrick Sanan .seealso: `PetscCallMPI()`
6419566063dSJacob Faibussowitsch M*/
6429566063dSJacob Faibussowitsch #define CHKERRMPI(...) PetscCallMPI(__VA_ARGS__)
6439566063dSJacob Faibussowitsch 
6449566063dSJacob Faibussowitsch /*MC
645989c49a2SBarry Smith   PetscCallAbort - Checks error code returned from PETSc function, if non-zero it aborts immediately by calling `MPI_Abort()`
6469566063dSJacob Faibussowitsch 
6479566063dSJacob Faibussowitsch   Synopsis:
6489566063dSJacob Faibussowitsch   #include <petscerror.h>
6499566063dSJacob Faibussowitsch   void PetscCallAbort(MPI_Comm comm, PetscErrorCode ierr)
6509566063dSJacob Faibussowitsch 
651c3339decSBarry Smith   Collective
6529566063dSJacob Faibussowitsch 
6539566063dSJacob Faibussowitsch   Input Parameters:
6549566063dSJacob Faibussowitsch + comm - the MPI communicator on which to abort
6559566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
6569566063dSJacob Faibussowitsch 
657667f096bSBarry Smith   Level: intermediate
658667f096bSBarry Smith 
6599566063dSJacob Faibussowitsch   Notes:
66087497f52SBarry Smith   This macro has identical type and usage semantics to `PetscCall()` with the important caveat
6619566063dSJacob Faibussowitsch   that this macro does not return. Instead, if ierr is nonzero it calls the PETSc error handler
66287497f52SBarry Smith   and then immediately calls `MPI_Abort()`. It can therefore be used anywhere.
6639566063dSJacob Faibussowitsch 
66487497f52SBarry Smith   As per `MPI_Abort()` semantics the communicator passed must be valid, although there is currently
66587497f52SBarry Smith   no attempt made at handling any potential errors from `MPI_Abort()`. Note that while
66687497f52SBarry Smith   `MPI_Abort()` is required to terminate only those processes which reside on comm, it is often
66787497f52SBarry Smith   the case that `MPI_Abort()` terminates *all* processes.
6689566063dSJacob Faibussowitsch 
6699566063dSJacob Faibussowitsch   Example Usage:
6709566063dSJacob Faibussowitsch .vb
6719566063dSJacob Faibussowitsch   PetscErrorCode boom(void) { return PETSC_ERR_MEM; }
6729566063dSJacob Faibussowitsch 
6739566063dSJacob Faibussowitsch   void foo(void)
6749566063dSJacob Faibussowitsch   {
6759566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type
6769566063dSJacob Faibussowitsch   }
6779566063dSJacob Faibussowitsch 
6789566063dSJacob Faibussowitsch   double bar(void)
6799566063dSJacob Faibussowitsch   {
6809566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type
6819566063dSJacob Faibussowitsch   }
6829566063dSJacob Faibussowitsch 
6839566063dSJacob Faibussowitsch   PetscCallAbort(MPI_COMM_NULL,boom()); // ERROR, communicator should be valid
6849566063dSJacob Faibussowitsch 
6859566063dSJacob Faibussowitsch   struct baz
6869566063dSJacob Faibussowitsch   {
6879566063dSJacob Faibussowitsch     baz()
6889566063dSJacob Faibussowitsch     {
6899566063dSJacob Faibussowitsch       PetscCallAbort(PETSC_COMM_SELF,boom()); // OK
6909566063dSJacob Faibussowitsch     }
6919566063dSJacob Faibussowitsch 
6929566063dSJacob Faibussowitsch     ~baz()
6939566063dSJacob Faibussowitsch     {
6949566063dSJacob Faibussowitsch       PetscCallAbort(PETSC_COMM_SELF,boom()); // OK (in fact the only way to handle PETSc errors)
6959566063dSJacob Faibussowitsch     }
6969566063dSJacob Faibussowitsch   };
6979566063dSJacob Faibussowitsch .ve
6989566063dSJacob Faibussowitsch 
699989c49a2SBarry Smith   Fortran Note:
700989c49a2SBarry Smith   Use `PetscCallA()`.
701989c49a2SBarry Smith 
702989c49a2SBarry Smith   Developer Note:
703989c49a2SBarry Smith   This should have the same name in Fortran as in C.
704989c49a2SBarry Smith 
705db781477SPatrick Sanan .seealso: `SETERRABORT()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`,
7065687f061SJacob Faibussowitsch           `SETERRQ()`, `CHKMEMQ`, `PetscCallMPI()`, `PetscCallCXXAbort()`
7079566063dSJacob Faibussowitsch M*/
7089566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
7099566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm, PetscErrorCode);
7109566063dSJacob Faibussowitsch void PetscCallContinue(PetscErrorCode);
7119566063dSJacob Faibussowitsch #else
7129371c9d4SSatish Balay   #define PetscCallAbort(comm, ...) \
7139371c9d4SSatish Balay     do { \
7147da6b3ccSLisandro Dalcin       PetscErrorCode ierr_petsc_call_abort_; \
7153ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
7167da6b3ccSLisandro Dalcin       ierr_petsc_call_abort_ = __VA_ARGS__; \
7173ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_abort_ != PETSC_SUCCESS)) { \
7183ba16761SJacob Faibussowitsch         ierr_petsc_call_abort_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_abort_, PETSC_ERROR_REPEAT, " "); \
7193ba16761SJacob Faibussowitsch         (void)MPI_Abort(comm, (PetscMPIInt)ierr_petsc_call_abort_); \
7209566063dSJacob Faibussowitsch       } \
7219566063dSJacob Faibussowitsch     } while (0)
7229371c9d4SSatish Balay   #define PetscCallContinue(...) \
7239371c9d4SSatish Balay     do { \
7247da6b3ccSLisandro Dalcin       PetscErrorCode ierr_petsc_call_continue_; \
7253ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
7267da6b3ccSLisandro Dalcin       ierr_petsc_call_continue_ = __VA_ARGS__; \
7273ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_continue_ != PETSC_SUCCESS)) { \
7283ba16761SJacob Faibussowitsch         ierr_petsc_call_continue_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_continue_, PETSC_ERROR_REPEAT, " "); \
7293ba16761SJacob Faibussowitsch         (void)ierr_petsc_call_continue_; \
7303ba16761SJacob Faibussowitsch       } \
7319566063dSJacob Faibussowitsch     } while (0)
7329566063dSJacob Faibussowitsch #endif
7339566063dSJacob Faibussowitsch 
7349566063dSJacob Faibussowitsch /*MC
7359566063dSJacob Faibussowitsch   CHKERRABORT - Checks error code returned from PETSc function. If non-zero it aborts immediately.
7369566063dSJacob Faibussowitsch 
7379566063dSJacob Faibussowitsch   Synopsis:
7389566063dSJacob Faibussowitsch   #include <petscerror.h>
7399566063dSJacob Faibussowitsch   void CHKERRABORT(MPI_Comm comm, PetscErrorCode ierr)
7409566063dSJacob Faibussowitsch 
7419566063dSJacob Faibussowitsch   Not Collective
7429566063dSJacob Faibussowitsch 
7439566063dSJacob Faibussowitsch   Input Parameters:
7449566063dSJacob Faibussowitsch + comm - the MPI communicator
7459566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
7469566063dSJacob Faibussowitsch 
7479566063dSJacob Faibussowitsch   Level: deprecated
7489566063dSJacob Faibussowitsch 
749667f096bSBarry Smith   Note:
750667f096bSBarry Smith   Deprecated in favor of `PetscCallAbort()`. This routine behaves identically to it.
751667f096bSBarry Smith 
752af27ebaaSBarry Smith .seealso: `PetscCallAbort()`, `PetscErrorCode`
7539566063dSJacob Faibussowitsch M*/
7549566063dSJacob Faibussowitsch #define CHKERRABORT(comm, ...) PetscCallAbort(comm, __VA_ARGS__)
7559566063dSJacob Faibussowitsch #define CHKERRCONTINUE(...)    PetscCallContinue(__VA_ARGS__)
7569566063dSJacob Faibussowitsch 
7579566063dSJacob Faibussowitsch /*MC
75887497f52SBarry Smith    CHKERRA - Fortran-only replacement for use of `CHKERRQ()` in the main program, which aborts immediately
759f388eb8bSPatrick Sanan 
760f388eb8bSPatrick Sanan    Synopsis:
761f388eb8bSPatrick Sanan    #include <petscsys.h>
762f388eb8bSPatrick Sanan    PetscErrorCode CHKERRA(PetscErrorCode ierr)
763f388eb8bSPatrick Sanan 
764f388eb8bSPatrick Sanan    Not Collective
765f388eb8bSPatrick Sanan 
7662fe279fdSBarry Smith    Input Parameter:
767f388eb8bSPatrick Sanan .  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
768f388eb8bSPatrick Sanan 
76987497f52SBarry Smith    Level: deprecated
770f388eb8bSPatrick Sanan 
77187497f52SBarry Smith    Note:
77287497f52SBarry Smith    This macro is rarely needed, normal usage is `PetscCallA()` in the main Fortran program.
773f388eb8bSPatrick Sanan 
774989c49a2SBarry Smith    Developer Note:
775989c49a2SBarry Smith    Why isn't this named `CHKERRABORT()` in Fortran?
776989c49a2SBarry Smith 
77787497f52SBarry Smith .seealso: `PetscCall()`, `PetscCallA()`, `PetscCallAbort()`, `CHKERRQ()`, `SETERRA()`, `SETERRQ()`, `SETERRABORT()`
778f388eb8bSPatrick Sanan M*/
779f388eb8bSPatrick Sanan 
7801e4f893aSSatish Balay PETSC_EXTERN PetscBool petscwaitonerrorflg;
7811e4f893aSSatish Balay PETSC_EXTERN PetscBool petscindebugger;
78235f00c14SToby Isaac PETSC_EXTERN PetscBool petscabortmpifinalize;
7837c66cc67SJunchao Zhang 
7847c66cc67SJunchao Zhang /*MC
785667f096bSBarry Smith    PETSCABORT - Call `MPI_Abort()` with an informative error code
7867c66cc67SJunchao Zhang 
7877c66cc67SJunchao Zhang    Synopsis:
7887c66cc67SJunchao Zhang    #include <petscsys.h>
7897c66cc67SJunchao Zhang    PETSCABORT(MPI_Comm comm, PetscErrorCode ierr)
7907c66cc67SJunchao Zhang 
7917cdbe19fSJose E. Roman    Collective; No Fortran Support
7927c66cc67SJunchao Zhang 
7937c66cc67SJunchao Zhang    Input Parameters:
79495bd0b28SBarry Smith +  comm - An MPI communicator, so that the error can be collective
7957c66cc67SJunchao Zhang -  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
7967c66cc67SJunchao Zhang 
797bf4d2887SBarry Smith    Level: advanced
7987c66cc67SJunchao Zhang 
799bf4d2887SBarry Smith    Notes:
800989c49a2SBarry Smith    If the option `-start_in_debugger` was used then this calls `abort()` to stop the program in the debugger.
801bf4d2887SBarry Smith 
802989c49a2SBarry Smith    if `PetscCIEnabledPortableErrorOutput` is set, which means the code is running in the PETSc test harness (make test),
803989c49a2SBarry Smith    and `comm` is `MPI_COMM_WORLD` it strives to exit cleanly without calling `MPI_Abort()` and instead calling `MPI_Finalize()`.
804660278c0SBarry Smith 
805989c49a2SBarry Smith    This is currently only used when an error propagates up to the C `main()` program and is detected by a `PetscCall()`, `PetscCallMPI()`,
806989c49a2SBarry Smith    or is set in `main()` with `SETERRQ()`. Abort calls such as `SETERRABORT()`,
807989c49a2SBarry Smith    `PetscCheckAbort()`, `PetscCallMPIAbort()`, and `PetscCallAbort()` always call `MPI_Abort()` and do not have any special
808989c49a2SBarry Smith    handling for the test harness.
809989c49a2SBarry Smith 
810989c49a2SBarry Smith    Developer Note:
811989c49a2SBarry Smith    Should the other abort calls also pass through this call instead of calling `MPI_Abort()` directly?
812989c49a2SBarry Smith 
813989c49a2SBarry Smith .seealso: `PetscError()`, `PetscCall()`, `SETERRABORT()`, `PetscCheckAbort()`, `PetscCallMPIAbort()`, `PetscCall()`, `PetscCallMPI()`,
814af27ebaaSBarry Smith           `PetscCallAbort()`, `MPI_Abort()`, `PetscErrorCode`
8157c66cc67SJunchao Zhang M*/
81629f88feaSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
81729f88feaSJacob Faibussowitsch void PETSCABORT(MPI_Comm, PetscErrorCode);
81829f88feaSJacob Faibussowitsch #else
8199371c9d4SSatish Balay   #define PETSCABORT(comm, ...) \
8209371c9d4SSatish Balay     do { \
8213ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_abort_; \
8223ba16761SJacob Faibussowitsch       if (petscwaitonerrorflg) { ierr_petsc_abort_ = PetscSleep(1000); } \
8233ba16761SJacob Faibussowitsch       if (petscindebugger) { \
8243ba16761SJacob Faibussowitsch         abort(); \
8253ba16761SJacob Faibussowitsch       } else { \
8267da6b3ccSLisandro Dalcin         PetscMPIInt size_; \
8273ba16761SJacob Faibussowitsch         ierr_petsc_abort_ = __VA_ARGS__; \
8287da6b3ccSLisandro Dalcin         MPI_Comm_size(comm, &size_); \
82935f00c14SToby Isaac         if (PetscCIEnabledPortableErrorOutput && (size_ == PetscGlobalSize || petscabortmpifinalize) && ierr_petsc_abort_ != PETSC_ERR_SIG) { \
8309371c9d4SSatish Balay           MPI_Finalize(); \
8319371c9d4SSatish Balay           exit(0); \
832660278c0SBarry Smith         } else if (PetscCIEnabledPortableErrorOutput && PetscGlobalSize == 1) { \
833660278c0SBarry Smith           exit(0); \
834660278c0SBarry Smith         } else { \
835660278c0SBarry Smith           MPI_Abort(comm, (PetscMPIInt)ierr_petsc_abort_); \
836660278c0SBarry Smith         } \
8373fcd9f07SJacob Faibussowitsch       } \
8387c66cc67SJunchao Zhang     } while (0)
83929f88feaSJacob Faibussowitsch #endif
840986eef2eSBarry Smith 
8419566063dSJacob Faibussowitsch #ifdef PETSC_CLANGUAGE_CXX
842986eef2eSBarry Smith   /*MC
8439566063dSJacob Faibussowitsch   PetscCallThrow - Checks error code, if non-zero it calls the C++ error handler which throws
8449566063dSJacob Faibussowitsch   an exception
845986eef2eSBarry Smith 
846986eef2eSBarry Smith   Synopsis:
8479566063dSJacob Faibussowitsch   #include <petscerror.h>
8489566063dSJacob Faibussowitsch   void PetscCallThrow(PetscErrorCode ierr)
849986eef2eSBarry Smith 
850986eef2eSBarry Smith   Not Collective
851986eef2eSBarry Smith 
8529566063dSJacob Faibussowitsch   Input Parameter:
853986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
854986eef2eSBarry Smith 
855667f096bSBarry Smith   Level: beginner
856667f096bSBarry Smith 
857986eef2eSBarry Smith   Notes:
858af27ebaaSBarry Smith   Requires PETSc to be configured with clanguage of c++. Throws a std::runtime_error() on error.
859986eef2eSBarry Smith 
86087497f52SBarry Smith   Once the error handler throws the exception you can use `PetscCallVoid()` which returns without
86187497f52SBarry Smith   an error code (bad idea since the error is ignored) or `PetscCallAbort()` to have `MPI_Abort()`
8629566063dSJacob Faibussowitsch   called immediately.
8639566063dSJacob Faibussowitsch 
864db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscCall()`, `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`,
865db781477SPatrick Sanan           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`
866986eef2eSBarry Smith M*/
8679371c9d4SSatish Balay   #define PetscCallThrow(...) \
8689371c9d4SSatish Balay     do { \
8693ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
8703ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_throw_ = __VA_ARGS__; \
8713ba16761SJacob 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); \
872ffc4695bSBarry Smith     } while (0)
87385614651SBarry Smith 
874cc26af49SMatthew Knepley   /*MC
875cc26af49SMatthew Knepley   CHKERRXX - Checks error code, if non-zero it calls the C++ error handler which throws an exception
876cc26af49SMatthew Knepley 
877cc26af49SMatthew Knepley   Synopsis:
8789566063dSJacob Faibussowitsch   #include <petscerror.h>
8793af045c5SBarry Smith   void CHKERRXX(PetscErrorCode ierr)
880cc26af49SMatthew Knepley 
881eca87e8dSBarry Smith   Not Collective
882cc26af49SMatthew Knepley 
8839566063dSJacob Faibussowitsch   Input Parameter:
8843af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
885cc26af49SMatthew Knepley 
8869566063dSJacob Faibussowitsch   Level: deprecated
887cc26af49SMatthew Knepley 
888667f096bSBarry Smith   Note:
889667f096bSBarry Smith   Deprecated in favor of `PetscCallThrow()`. This routine behaves identically to it.
890667f096bSBarry Smith 
891db781477SPatrick Sanan .seealso: `PetscCallThrow()`
892cc26af49SMatthew Knepley M*/
8939566063dSJacob Faibussowitsch   #define CHKERRXX(...) PetscCallThrow(__VA_ARGS__)
894fd705b32SMatthew Knepley #endif
895fd705b32SMatthew Knepley 
8963ba16761SJacob Faibussowitsch #define PetscCallCXX_Private(__SETERR_FUNC__, __COMM__, ...) \
8973ba16761SJacob Faibussowitsch   do { \
8983ba16761SJacob Faibussowitsch     PetscStackUpdateLine; \
8993ba16761SJacob Faibussowitsch     try { \
9003ba16761SJacob Faibussowitsch       __VA_ARGS__; \
9013ba16761SJacob Faibussowitsch     } catch (const std::exception &e) { \
9023ba16761SJacob Faibussowitsch       __SETERR_FUNC__(__COMM__, PETSC_ERR_LIB, "%s", e.what()); \
9033ba16761SJacob Faibussowitsch     } \
9043ba16761SJacob Faibussowitsch   } while (0)
9053ba16761SJacob Faibussowitsch 
9063f520e80SJunchao Zhang /*MC
9079566063dSJacob Faibussowitsch   PetscCallCXX - Checks C++ function calls and if they throw an exception, catch it and then
9089566063dSJacob Faibussowitsch   return a PETSc error code
9093f520e80SJunchao Zhang 
9103f520e80SJunchao Zhang   Synopsis:
9119566063dSJacob Faibussowitsch   #include <petscerror.h>
9125687f061SJacob Faibussowitsch   void PetscCallCXX(...) noexcept;
9133f520e80SJunchao Zhang 
9143f520e80SJunchao Zhang   Not Collective
9153f520e80SJunchao Zhang 
9169566063dSJacob Faibussowitsch   Input Parameter:
9175687f061SJacob Faibussowitsch . __VA_ARGS__ - An arbitrary expression
9185687f061SJacob Faibussowitsch 
9195687f061SJacob Faibussowitsch   Level: beginner
9203f520e80SJunchao Zhang 
9213f520e80SJunchao Zhang   Notes:
9225687f061SJacob Faibussowitsch   `PetscCallCXX(...)` is a macro replacement for
9239566063dSJacob Faibussowitsch .vb
9249566063dSJacob Faibussowitsch   try {
9255687f061SJacob Faibussowitsch     __VA_ARGS__;
9269566063dSJacob Faibussowitsch   } catch (const std::exception& e) {
9279566063dSJacob Faibussowitsch     return ConvertToPetscErrorCode(e);
9289566063dSJacob Faibussowitsch   }
9299566063dSJacob Faibussowitsch .ve
9309566063dSJacob Faibussowitsch   Due to the fact that it catches any (reasonable) exception, it is essentially noexcept.
9313f520e80SJunchao Zhang 
9325687f061SJacob Faibussowitsch   If you cannot return a `PetscErrorCode` use `PetscCallCXXAbort()` instead.
9335687f061SJacob Faibussowitsch 
9349566063dSJacob Faibussowitsch   Example Usage:
9359566063dSJacob Faibussowitsch .vb
9369566063dSJacob Faibussowitsch   void foo(void) { throw std::runtime_error("error"); }
9373f520e80SJunchao Zhang 
9389566063dSJacob Faibussowitsch   void bar()
9399566063dSJacob Faibussowitsch   {
940e8952933SJacob Faibussowitsch     PetscCallCXX(foo()); // ERROR bar() does not return PetscErrorCode
9419566063dSJacob Faibussowitsch   }
9429566063dSJacob Faibussowitsch 
9439566063dSJacob Faibussowitsch   PetscErrorCode baz()
9449566063dSJacob Faibussowitsch   {
9459566063dSJacob Faibussowitsch     PetscCallCXX(foo()); // OK
9469566063dSJacob Faibussowitsch 
9479566063dSJacob Faibussowitsch     PetscCallCXX(
9489566063dSJacob Faibussowitsch       bar();
949d5b43468SJose E. Roman       foo(); // OK multiple statements allowed
9509566063dSJacob Faibussowitsch     );
9519566063dSJacob Faibussowitsch   }
9529566063dSJacob Faibussowitsch 
9539566063dSJacob Faibussowitsch   struct bop
9549566063dSJacob Faibussowitsch   {
9559566063dSJacob Faibussowitsch     bop()
9569566063dSJacob Faibussowitsch     {
957e8952933SJacob Faibussowitsch       PetscCallCXX(foo()); // ERROR returns PetscErrorCode, cannot be used in constructors
9589566063dSJacob Faibussowitsch     }
9599566063dSJacob Faibussowitsch   };
9609566063dSJacob Faibussowitsch 
961e8952933SJacob Faibussowitsch   // ERROR contains do-while, cannot be used as function-try block
9629566063dSJacob Faibussowitsch   PetscErrorCode qux() PetscCallCXX(
9639566063dSJacob Faibussowitsch     bar();
9649566063dSJacob Faibussowitsch     baz();
9659566063dSJacob Faibussowitsch     foo();
9669566063dSJacob Faibussowitsch     return 0;
9679566063dSJacob Faibussowitsch   )
9689566063dSJacob Faibussowitsch .ve
9699566063dSJacob Faibussowitsch 
9705687f061SJacob Faibussowitsch .seealso: `PetscCallCXXAbort()`, `PetscCallThrow()`, `SETERRQ()`, `PetscCall()`,
9715687f061SJacob Faibussowitsch           `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
9725687f061SJacob Faibussowitsch           `PetscError()`, `CHKMEMQ`
9733f520e80SJunchao Zhang M*/
9743ba16761SJacob Faibussowitsch #define PetscCallCXX(...) PetscCallCXX_Private(SETERRQ, PETSC_COMM_SELF, __VA_ARGS__)
9755687f061SJacob Faibussowitsch 
9765687f061SJacob Faibussowitsch /*MC
9775687f061SJacob Faibussowitsch   PetscCallCXXAbort - Like `PetscCallCXX()` but calls `MPI_Abort()` instead of returning an
9785687f061SJacob Faibussowitsch   error-code
9795687f061SJacob Faibussowitsch 
9805687f061SJacob Faibussowitsch   Synopsis:
9815687f061SJacob Faibussowitsch   #include <petscerror.h>
9825687f061SJacob Faibussowitsch   void PetscCallCXXAbort(MPI_Comm comm, ...) noexcept;
9835687f061SJacob Faibussowitsch 
98420f4b53cSBarry Smith   Collective; No Fortran Support
9855687f061SJacob Faibussowitsch 
9865687f061SJacob Faibussowitsch   Input Parameters:
9875687f061SJacob Faibussowitsch + comm        - The MPI communicator to abort on
9885687f061SJacob Faibussowitsch - __VA_ARGS__ - An arbitrary expression
9895687f061SJacob Faibussowitsch 
9905687f061SJacob Faibussowitsch   Level: beginner
9915687f061SJacob Faibussowitsch 
9925687f061SJacob Faibussowitsch   Notes:
9935687f061SJacob Faibussowitsch   This macro may be used to check C++ expressions for exceptions in cases where you cannot
9945687f061SJacob Faibussowitsch   return an error code. This includes constructors, destructors, copy/move assignment functions
9955687f061SJacob Faibussowitsch   or constructors among others.
9965687f061SJacob Faibussowitsch 
9975687f061SJacob Faibussowitsch   If an exception is caught, the macro calls `SETERRABORT()` on `comm`. The exception must
9985687f061SJacob Faibussowitsch   derive from `std::exception` in order to be caught.
9995687f061SJacob Faibussowitsch 
10005687f061SJacob Faibussowitsch   If the routine _can_ return an error-code it is highly advised to use `PetscCallCXX()`
10015687f061SJacob Faibussowitsch   instead.
10025687f061SJacob Faibussowitsch 
10035687f061SJacob Faibussowitsch   See `PetscCallCXX()` for additional discussion.
10045687f061SJacob Faibussowitsch 
10055687f061SJacob Faibussowitsch   Example Usage:
10065687f061SJacob Faibussowitsch .vb
10075687f061SJacob Faibussowitsch   class Foo
10085687f061SJacob Faibussowitsch   {
10095687f061SJacob Faibussowitsch     std::vector<int> data_;
10105687f061SJacob Faibussowitsch 
10115687f061SJacob Faibussowitsch   public:
10125687f061SJacob Faibussowitsch     // normally std::vector::reserve() may raise an exception, but since we handle it with
10135687f061SJacob Faibussowitsch     // PetscCallCXXAbort() we may mark this routine as noexcept!
10145687f061SJacob Faibussowitsch     Foo() noexcept
10155687f061SJacob Faibussowitsch     {
10165687f061SJacob Faibussowitsch       PetscCallCXXAbort(PETSC_COMM_SELF, data_.reserve(10));
10175687f061SJacob Faibussowitsch     }
10185687f061SJacob Faibussowitsch   };
10195687f061SJacob Faibussowitsch 
10205687f061SJacob Faibussowitsch   std::vector<int> bar()
10215687f061SJacob Faibussowitsch   {
10225687f061SJacob Faibussowitsch     std::vector<int> v;
10235687f061SJacob Faibussowitsch 
10245687f061SJacob Faibussowitsch     PetscFunctionBegin;
10255687f061SJacob Faibussowitsch     // OK!
10265687f061SJacob Faibussowitsch     PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1));
10275687f061SJacob Faibussowitsch     PetscFunctionReturn(v);
10285687f061SJacob Faibussowitsch   }
10295687f061SJacob Faibussowitsch 
10305687f061SJacob Faibussowitsch   PetscErrorCode baz()
10315687f061SJacob Faibussowitsch   {
10325687f061SJacob Faibussowitsch     std::vector<int> v;
10335687f061SJacob Faibussowitsch 
10345687f061SJacob Faibussowitsch     PetscFunctionBegin;
10355687f061SJacob Faibussowitsch     // WRONG! baz() returns a PetscErrorCode, prefer PetscCallCXX() instead
10365687f061SJacob Faibussowitsch     PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1));
10373ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
10385687f061SJacob Faibussowitsch   }
10395687f061SJacob Faibussowitsch .ve
10405687f061SJacob Faibussowitsch 
10415687f061SJacob Faibussowitsch .seealso: `PetscCallCXX()`, `SETERRABORT()`, `PetscCallAbort()`
10425687f061SJacob Faibussowitsch M*/
10433ba16761SJacob Faibussowitsch #define PetscCallCXXAbort(comm, ...) PetscCallCXX_Private(SETERRABORT, comm, __VA_ARGS__)
10443f520e80SJunchao Zhang 
104530de9b25SBarry Smith /*MC
10469566063dSJacob Faibussowitsch   CHKERRCXX - Checks C++ function calls and if they throw an exception, catch it and then
10479566063dSJacob Faibussowitsch   return a PETSc error code
10489566063dSJacob Faibussowitsch 
10499566063dSJacob Faibussowitsch   Synopsis:
10509566063dSJacob Faibussowitsch   #include <petscerror.h>
10519566063dSJacob Faibussowitsch   void CHKERRCXX(func) noexcept;
10529566063dSJacob Faibussowitsch 
10539566063dSJacob Faibussowitsch   Not Collective
10549566063dSJacob Faibussowitsch 
10559566063dSJacob Faibussowitsch   Input Parameter:
10569566063dSJacob Faibussowitsch . func - C++ function calls
10579566063dSJacob Faibussowitsch 
10589566063dSJacob Faibussowitsch   Level: deprecated
10599566063dSJacob Faibussowitsch 
1060667f096bSBarry Smith   Note:
1061667f096bSBarry Smith   Deprecated in favor of `PetscCallCXX()`. This routine behaves identically to it.
1062667f096bSBarry Smith 
1063db781477SPatrick Sanan .seealso: `PetscCallCXX()`
10649566063dSJacob Faibussowitsch M*/
10659566063dSJacob Faibussowitsch #define CHKERRCXX(...) PetscCallCXX(__VA_ARGS__)
10669566063dSJacob Faibussowitsch 
10679566063dSJacob Faibussowitsch /*MC
106830de9b25SBarry Smith    CHKMEMQ - Checks the memory for corruption, calls error handler if any is detected
106930de9b25SBarry Smith 
107030de9b25SBarry Smith    Synopsis:
1071aaa7dc30SBarry Smith    #include <petscsys.h>
107291d3bdf4SKris Buschelman    CHKMEMQ;
107330de9b25SBarry Smith 
1074eca87e8dSBarry Smith    Not Collective
1075eca87e8dSBarry Smith 
107630de9b25SBarry Smith    Level: beginner
107730de9b25SBarry Smith 
107830de9b25SBarry Smith    Notes:
1079af27ebaaSBarry Smith    We recommend using Valgrind <https://petsc.org/release/faq/#valgrind> or for NVIDIA CUDA systems
1080af27ebaaSBarry Smith    <https://docs.nvidia.com/cuda/cuda-memcheck/index.html> for finding memory problems. The ``CHKMEMQ`` macro is useful on systems that
10815ed36255SBarry Smith    do not have valgrind, but is not as good as valgrind or cuda-memcheck.
10821957e957SBarry Smith 
1083667f096bSBarry Smith    Must run with the option `-malloc_debug` (`-malloc_test` in debug mode; or if `PetscMallocSetDebug()` called) to enable this option
108430de9b25SBarry Smith 
108530de9b25SBarry Smith    Once the error handler is called the calling function is then returned from with the given error code.
108630de9b25SBarry Smith 
108730de9b25SBarry Smith    By defaults prints location where memory that is corrupted was allocated.
108830de9b25SBarry Smith 
10898af9f190SBarry Smith    Use `CHKMEMA` for functions that return `void`
1090f621e05eSBarry Smith 
1091db781477SPatrick Sanan .seealso: `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `SETERRQ()`, `PetscMallocValidate()`
109230de9b25SBarry Smith M*/
10936d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
1094064a246eSJacob Faibussowitsch   #define CHKMEMQ
1095064a246eSJacob Faibussowitsch   #define CHKMEMA
10966d210af2SJacob Faibussowitsch #else
10979371c9d4SSatish Balay   #define CHKMEMQ \
10989371c9d4SSatish Balay     do { \
10993ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_memq_ = PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__); \
11003ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_memq_ != PETSC_SUCCESS)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_memq_, PETSC_ERROR_REPEAT, " "); \
110186d09637SLisandro Dalcin     } while (0)
11026d210af2SJacob Faibussowitsch   #define CHKMEMA PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__)
1103064a246eSJacob Faibussowitsch #endif
11049566063dSJacob Faibussowitsch 
1105668f157eSBarry Smith /*E
1106668f157eSBarry Smith   PetscErrorType - passed to the PETSc error handling routines indicating if this is the first or a later call to the error handlers
1107668f157eSBarry Smith 
1108668f157eSBarry Smith   Level: advanced
1109668f157eSBarry Smith 
1110667f096bSBarry Smith   Note:
111187497f52SBarry Smith   `PETSC_ERROR_IN_CXX` indicates the error was detected in C++ and an exception should be generated
1112d736bfebSBarry Smith 
1113667f096bSBarry Smith   Developer Note:
1114667f096bSBarry Smith   This is currently used to decide when to print the detailed information about the run in `PetscTraceBackErrorHandler()`
1115668f157eSBarry Smith 
111687497f52SBarry Smith .seealso: `PetscError()`, `SETERRQ()`
1117668f157eSBarry Smith E*/
11189371c9d4SSatish Balay typedef enum {
11199371c9d4SSatish Balay   PETSC_ERROR_INITIAL = 0,
11209371c9d4SSatish Balay   PETSC_ERROR_REPEAT  = 1,
11219371c9d4SSatish Balay   PETSC_ERROR_IN_CXX  = 2
11229371c9d4SSatish Balay } PetscErrorType;
11234b209cf6SBarry Smith 
1124eb9e708aSLisandro Dalcin #if defined(__clang_analyzer__)
1125eb9e708aSLisandro Dalcin __attribute__((analyzer_noreturn))
1126eb9e708aSLisandro Dalcin #endif
11279371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode
11289371c9d4SSatish Balay PetscError(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, ...) PETSC_ATTRIBUTE_COLD PETSC_ATTRIBUTE_FORMAT(7, 8);
1129eb9e708aSLisandro Dalcin 
1130014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscErrorPrintfInitialize(void);
11313ba16761SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscErrorMessage(PetscErrorCode, const char *[], char **);
1132d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscTraceBackErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1133d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscIgnoreErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1134d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscEmacsClientErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1135d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscMPIAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1136d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1137d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAttachDebuggerErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1138d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscReturnErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1139efca3c55SSatish Balay PETSC_EXTERN PetscErrorCode PetscPushErrorHandler(PetscErrorCode (*handler)(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *), void *);
1140014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopErrorHandler(void);
11418d359177SBarry Smith PETSC_EXTERN PetscErrorCode PetscSignalHandlerDefault(int, void *);
1142014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPushSignalHandler(PetscErrorCode (*)(int, void *), void *);
1143014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopSignalHandler(void);
114428559dc8SJed Brown PETSC_EXTERN PetscErrorCode PetscCheckPointerSetIntensity(PetscInt);
1145c2a741eeSJunchao Zhang PETSC_EXTERN void           PetscSignalSegvCheckPointerOrMpi(void);
1146edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 13, 0, "PetscSignalSegvCheckPointerOrMpi()", ) static inline void PetscSignalSegvCheckPointer(void)
1147d71ae5a4SJacob Faibussowitsch {
11489371c9d4SSatish Balay   PetscSignalSegvCheckPointerOrMpi();
11499371c9d4SSatish Balay }
1150329f5518SBarry Smith 
1151639ff905SBarry Smith /*MC
1152639ff905SBarry Smith     PetscErrorPrintf - Prints error messages.
1153639ff905SBarry Smith 
1154639ff905SBarry Smith    Synopsis:
1155aaa7dc30SBarry Smith     #include <petscsys.h>
1156639ff905SBarry Smith      PetscErrorCode (*PetscErrorPrintf)(const char format[], ...);
1157639ff905SBarry Smith 
11587cdbe19fSJose E. Roman     Not Collective; No Fortran Support
11597cdbe19fSJose E. Roman 
1160f899ff85SJose E. Roman     Input Parameter:
1161cd05f99aSJacob Faibussowitsch .   format - the usual `printf()` format string
1162639ff905SBarry Smith 
1163639ff905SBarry Smith    Options Database Keys:
11641957e957SBarry Smith +  -error_output_stdout - cause error messages to be printed to stdout instead of the (default) stderr
1165e1bc860dSBarry Smith -  -error_output_none   - to turn off all printing of error messages (does not change the way the error is handled.)
1166639ff905SBarry Smith 
1167cf53795eSBarry Smith    Level: developer
1168cf53795eSBarry Smith 
116995452b02SPatrick Sanan    Notes:
117095452b02SPatrick Sanan    Use
1171667f096bSBarry Smith .vb
117295bd0b28SBarry Smith      PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the error is handled) and
1173667f096bSBarry Smith      PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on or you can use your own function
1174667f096bSBarry Smith .ve
1175639ff905SBarry Smith    Use
1176667f096bSBarry Smith .vb
117787497f52SBarry Smith      `PETSC_STDERR` = FILE* obtained from a file open etc. to have stderr printed to the file.
117887497f52SBarry Smith      `PETSC_STDOUT` = FILE* obtained from a file open etc. to have stdout printed to the file.
1179667f096bSBarry Smith .ve
1180639ff905SBarry Smith    Use
1181af27ebaaSBarry Smith .vb
118287497f52SBarry Smith       `PetscPushErrorHandler()` to provide your own error handler that determines what kind of messages to print
1183af27ebaaSBarry Smith .ve
1184639ff905SBarry Smith 
1185db781477SPatrick Sanan .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscHelpPrintf()`, `PetscPrintf()`, `PetscPushErrorHandler()`, `PetscVFPrintf()`, `PetscHelpPrintf()`
1186639ff905SBarry Smith M*/
11873ca90d2dSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode (*PetscErrorPrintf)(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1188639ff905SBarry Smith 
1189cf0818bdSBarry Smith /*E
1190cf0818bdSBarry Smith      PetscFPTrap - types of floating point exceptions that may be trapped
1191cf0818bdSBarry Smith 
119287497f52SBarry Smith      Currently only `PETSC_FP_TRAP_OFF` and `PETSC_FP_TRAP_ON` are handled. All others are treated as `PETSC_FP_TRAP_ON`.
1193cf0818bdSBarry Smith 
1194cf0818bdSBarry Smith      Level: intermediate
1195cf0818bdSBarry Smith 
11967de69702SBarry Smith .seealso: `PetscSetFPTrap()`, `PetscFPTrapPush()`
1197cf0818bdSBarry Smith  E*/
11989371c9d4SSatish Balay typedef enum {
11999371c9d4SSatish Balay   PETSC_FP_TRAP_OFF      = 0,
12009371c9d4SSatish Balay   PETSC_FP_TRAP_INDIV    = 1,
12019371c9d4SSatish Balay   PETSC_FP_TRAP_FLTOPERR = 2,
12029371c9d4SSatish Balay   PETSC_FP_TRAP_FLTOVF   = 4,
12039371c9d4SSatish Balay   PETSC_FP_TRAP_FLTUND   = 8,
12049371c9d4SSatish Balay   PETSC_FP_TRAP_FLTDIV   = 16,
12059371c9d4SSatish Balay   PETSC_FP_TRAP_FLTINEX  = 32
12069371c9d4SSatish Balay } PetscFPTrap;
1207bd2b07b1SBarry Smith #define PETSC_FP_TRAP_ON (PetscFPTrap)(PETSC_FP_TRAP_INDIV | PETSC_FP_TRAP_FLTOPERR | PETSC_FP_TRAP_FLTOVF | PETSC_FP_TRAP_FLTDIV | PETSC_FP_TRAP_FLTINEX)
1208014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSetFPTrap(PetscFPTrap);
1209014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPush(PetscFPTrap);
1210014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPop(void);
1211aba4c478SBarry Smith PETSC_EXTERN PetscErrorCode PetscDetermineInitialFPTrap(void);
121254a8ef01SBarry Smith 
12133a40ed3dSBarry Smith /*
12143a40ed3dSBarry Smith       Allows the code to build a stack frame as it runs
12153a40ed3dSBarry Smith */
12163a40ed3dSBarry Smith 
121799cd645aSJed Brown #define PETSCSTACKSIZE 64
12183a40ed3dSBarry Smith typedef struct {
12190e33f6ddSBarry Smith   const char *function[PETSCSTACKSIZE];
12200e33f6ddSBarry Smith   const char *file[PETSCSTACKSIZE];
1221184914b5SBarry Smith   int         line[PETSCSTACKSIZE];
1222362febeeSStefano Zampini   int         petscroutine[PETSCSTACKSIZE]; /* 0 external called from petsc, 1 petsc functions, 2 petsc user functions */
1223184914b5SBarry Smith   int         currentsize;
1224a2f94806SJed Brown   int         hotdepth;
12254be741a6SBarry Smith   PetscBool   check; /* option to check for correct Push/Pop semantics, true for default petscstack but not other stacks */
12263a40ed3dSBarry Smith } PetscStack;
1227dfb7d7afSStefano Zampini #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
122827104ee2SJacob Faibussowitsch PETSC_EXTERN PetscStack petscstack;
122927104ee2SJacob Faibussowitsch #endif
12303a40ed3dSBarry Smith 
12315d12eec7SSatish Balay #if defined(PETSC_SERIALIZE_FUNCTIONS)
12325d12eec7SSatish Balay   #include <petsc/private/petscfptimpl.h>
12335d12eec7SSatish Balay   /*
12345d12eec7SSatish Balay    Registers the current function into the global function pointer to function name table
12355d12eec7SSatish Balay 
12365d12eec7SSatish Balay    Have to fix this to handle errors but cannot return error since used in PETSC_VIEWER_DRAW_() etc
12375d12eec7SSatish Balay */
12389371c9d4SSatish Balay   #define PetscRegister__FUNCT__() \
12399371c9d4SSatish Balay     do { \
12405d12eec7SSatish Balay       static PetscBool __chked = PETSC_FALSE; \
12415d12eec7SSatish Balay       if (!__chked) { \
12429371c9d4SSatish Balay         void *ptr; \
12433ba16761SJacob Faibussowitsch         PetscCallAbort(PETSC_COMM_SELF, PetscDLSym(NULL, PETSC_FUNCTION_NAME, &ptr)); \
12445d12eec7SSatish Balay         __chked = PETSC_TRUE; \
12459371c9d4SSatish Balay       } \
12469371c9d4SSatish Balay     } while (0)
12475d12eec7SSatish Balay #else
12485d12eec7SSatish Balay   #define PetscRegister__FUNCT__()
12495d12eec7SSatish Balay #endif
12505d12eec7SSatish Balay 
1251eae3dc7dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) || defined(__clang_analyzer__)
12526d210af2SJacob Faibussowitsch   #define PetscStackPushNoCheck(funct, petsc_routine, hot)
125337154d10SBarry Smith   #define PetscStackUpdateLine
1254792fecdfSBarry Smith   #define PetscStackPushExternal(funct)
12556d210af2SJacob Faibussowitsch   #define PetscStackPopNoCheck
12566d210af2SJacob Faibussowitsch   #define PetscStackClearTop
12576d210af2SJacob Faibussowitsch   #define PetscFunctionBegin
12586d210af2SJacob Faibussowitsch   #define PetscFunctionBeginUser
12596d210af2SJacob Faibussowitsch   #define PetscFunctionBeginHot
12600a57284eSJacob Faibussowitsch   #define PetscFunctionReturn(...)  return __VA_ARGS__
12616d210af2SJacob Faibussowitsch   #define PetscFunctionReturnVoid() return
12626d210af2SJacob Faibussowitsch   #define PetscStackPop
12636d210af2SJacob Faibussowitsch   #define PetscStackPush(f)
1264dfb7d7afSStefano Zampini #elif defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1265660278c0SBarry Smith 
12669371c9d4SSatish Balay   #define PetscStackPush_Private(stack__, file__, func__, line__, petsc_routine__, hot__) \
12679371c9d4SSatish Balay     do { \
12685a96b57dSJacob Faibussowitsch       if (stack__.currentsize < PETSCSTACKSIZE) { \
12695a96b57dSJacob Faibussowitsch         stack__.function[stack__.currentsize] = func__; \
1270ef1023bdSBarry Smith         if (petsc_routine__) { \
1271ef1023bdSBarry Smith           stack__.file[stack__.currentsize] = file__; \
12725a96b57dSJacob Faibussowitsch           stack__.line[stack__.currentsize] = line__; \
1273ef1023bdSBarry Smith         } else { \
1274648bc8c4SBarry Smith           stack__.file[stack__.currentsize] = PETSC_NULLPTR; \
1275ef1023bdSBarry Smith           stack__.line[stack__.currentsize] = 0; \
1276ef1023bdSBarry Smith         } \
12775a96b57dSJacob Faibussowitsch         stack__.petscroutine[stack__.currentsize] = petsc_routine__; \
12785a96b57dSJacob Faibussowitsch       } \
12795a96b57dSJacob Faibussowitsch       ++stack__.currentsize; \
12805a96b57dSJacob Faibussowitsch       stack__.hotdepth += (hot__ || stack__.hotdepth); \
12815a96b57dSJacob Faibussowitsch     } while (0)
12825a96b57dSJacob Faibussowitsch 
12834be741a6SBarry Smith   /* uses PetscCheckAbort() because may be used in a function that does not return an error code */
12849371c9d4SSatish Balay   #define PetscStackPop_Private(stack__, func__) \
12859371c9d4SSatish Balay     do { \
12864be741a6SBarry 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__); \
12875a96b57dSJacob Faibussowitsch       if (--stack__.currentsize < PETSCSTACKSIZE) { \
12889371c9d4SSatish 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", \
12899371c9d4SSatish Balay                         stack__.function[stack__.currentsize], stack__.file[stack__.currentsize], stack__.line[stack__.currentsize], func__, __FILE__, __LINE__); \
12905a96b57dSJacob Faibussowitsch         stack__.function[stack__.currentsize]     = PETSC_NULLPTR; \
12915a96b57dSJacob Faibussowitsch         stack__.file[stack__.currentsize]         = PETSC_NULLPTR; \
12925a96b57dSJacob Faibussowitsch         stack__.line[stack__.currentsize]         = 0; \
12935a96b57dSJacob Faibussowitsch         stack__.petscroutine[stack__.currentsize] = 0; \
12945a96b57dSJacob Faibussowitsch       } \
12955a96b57dSJacob Faibussowitsch       stack__.hotdepth = PetscMax(stack__.hotdepth - 1, 0); \
12965a96b57dSJacob Faibussowitsch     } while (0)
12975a96b57dSJacob Faibussowitsch 
1298586f9135SBarry Smith   /*MC
1299586f9135SBarry Smith    PetscStackPushNoCheck - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is
1300586f9135SBarry Smith    currently in the source code.
1301586f9135SBarry Smith 
1302586f9135SBarry Smith    Synopsis:
1303586f9135SBarry Smith    #include <petscsys.h>
1304586f9135SBarry Smith    void PetscStackPushNoCheck(char *funct,int petsc_routine,PetscBool hot);
1305586f9135SBarry Smith 
13067cdbe19fSJose E. Roman    Not Collective
13077cdbe19fSJose E. Roman 
1308586f9135SBarry Smith    Input Parameters:
1309586f9135SBarry Smith +  funct - the function name
1310586f9135SBarry Smith .  petsc_routine - 2 user function, 1 PETSc function, 0 some other function
1311586f9135SBarry Smith -  hot - indicates that the function may be called often so expensive error checking should be turned off inside the function
1312586f9135SBarry Smith 
1313586f9135SBarry Smith    Level: developer
1314586f9135SBarry Smith 
1315586f9135SBarry Smith    Notes:
1316586f9135SBarry 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
131787497f52SBarry 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
1318586f9135SBarry Smith    help debug the problem.
1319586f9135SBarry Smith 
1320ef1023bdSBarry Smith    This version does not check the memory corruption (an expensive operation), use `PetscStackPush()` to check the memory.
1321ef1023bdSBarry Smith 
1322792fecdfSBarry Smith    Use `PetscStackPushExternal()` for a function call that is about to be made to a non-PETSc or user function (such as BLAS etc).
1323ef1023bdSBarry Smith 
1324586f9135SBarry Smith    The default stack is a global variable called `petscstack`.
1325586f9135SBarry Smith 
1326586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1327ef1023bdSBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPush()`, `PetscStackPop`,
1328792fecdfSBarry Smith           `PetscStackPushExternal()`
1329586f9135SBarry Smith M*/
13309371c9d4SSatish Balay   #define PetscStackPushNoCheck(funct, petsc_routine, hot) \
13319371c9d4SSatish Balay     do { \
1332e04113cfSBarry Smith       PetscStackSAWsTakeAccess(); \
13335a96b57dSJacob Faibussowitsch       PetscStackPush_Private(petscstack, __FILE__, funct, __LINE__, petsc_routine, hot); \
1334e04113cfSBarry Smith       PetscStackSAWsGrantAccess(); \
1335441dd030SJed Brown     } while (0)
1336441dd030SJed Brown 
1337586f9135SBarry Smith   /*MC
133887497f52SBarry Smith    PetscStackUpdateLine - in a function that has a `PetscFunctionBegin` or `PetscFunctionBeginUser` updates the stack line number to the
133937154d10SBarry Smith    current line number.
134037154d10SBarry Smith 
134137154d10SBarry Smith    Synopsis:
134237154d10SBarry Smith    #include <petscsys.h>
134337154d10SBarry Smith    void PetscStackUpdateLine
134437154d10SBarry Smith 
13457cdbe19fSJose E. Roman    Not Collective
13467cdbe19fSJose E. Roman 
134737154d10SBarry Smith    Level: developer
134837154d10SBarry Smith 
134937154d10SBarry Smith    Notes:
135087497f52SBarry Smith    Using `PetscCall()` and friends automatically handles this process
135187497f52SBarry Smith 
135237154d10SBarry 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
135337154d10SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
135437154d10SBarry Smith    help debug the problem.
135537154d10SBarry Smith 
135695bd0b28SBarry Smith    The default stack is a global variable called `petscstack`.
135737154d10SBarry Smith 
135837154d10SBarry Smith    This is used by `PetscCall()` and is otherwise not like to be needed
135937154d10SBarry Smith 
136037154d10SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`, `PetscCall()`
136137154d10SBarry Smith M*/
136237154d10SBarry Smith   #define PetscStackUpdateLine \
13633ba16761SJacob Faibussowitsch     do { \
1364f1e99121SPierre Jolivet       if (petscstack.currentsize > 0 && petscstack.currentsize < PETSCSTACKSIZE && petscstack.function[petscstack.currentsize - 1] == PETSC_FUNCTION_NAME) { petscstack.line[petscstack.currentsize - 1] = __LINE__; } \
13653ba16761SJacob Faibussowitsch     } while (0)
136637154d10SBarry Smith 
136737154d10SBarry Smith   /*MC
1368792fecdfSBarry Smith    PetscStackPushExternal - Pushes a new function name onto the PETSc default stack that tracks where the running program is
1369ef1023bdSBarry Smith    currently in the source code. Does not include the filename or line number since this is called by the calling routine
1370ef1023bdSBarry Smith    for non-PETSc or user functions.
1371ef1023bdSBarry Smith 
1372ef1023bdSBarry Smith    Synopsis:
1373ef1023bdSBarry Smith    #include <petscsys.h>
1374792fecdfSBarry Smith    void PetscStackPushExternal(char *funct);
1375ef1023bdSBarry Smith 
13767cdbe19fSJose E. Roman    Not Collective
13777cdbe19fSJose E. Roman 
13782fe279fdSBarry Smith    Input Parameter:
1379ef1023bdSBarry Smith .  funct - the function name
1380ef1023bdSBarry Smith 
1381ef1023bdSBarry Smith    Level: developer
1382ef1023bdSBarry Smith 
1383ef1023bdSBarry Smith    Notes:
138487497f52SBarry Smith    Using `PetscCallExternal()` and friends automatically handles this process
138587497f52SBarry Smith 
1386ef1023bdSBarry 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
1387ef1023bdSBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1388ef1023bdSBarry Smith    help debug the problem.
1389ef1023bdSBarry Smith 
1390ef1023bdSBarry Smith    The default stack is a global variable called `petscstack`.
1391ef1023bdSBarry Smith 
1392ef1023bdSBarry Smith    This is to be used when calling an external package function such as a BLAS function.
1393ef1023bdSBarry Smith 
1394ef1023bdSBarry Smith    This also updates the stack line number for the current stack function.
1395ef1023bdSBarry Smith 
1396ef1023bdSBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1397ef1023bdSBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1398ef1023bdSBarry Smith M*/
13999371c9d4SSatish Balay   #define PetscStackPushExternal(funct) \
14009371c9d4SSatish Balay     do { \
14019371c9d4SSatish Balay       PetscStackUpdateLine; \
14029371c9d4SSatish Balay       PetscStackPushNoCheck(funct, 0, PETSC_TRUE); \
1403a8f51744SPierre Jolivet     } while (0)
1404ef1023bdSBarry Smith 
1405ef1023bdSBarry Smith   /*MC
1406586f9135SBarry Smith    PetscStackPopNoCheck - Pops a function name from the PETSc default stack that tracks where the running program is
1407586f9135SBarry Smith    currently in the source code.
1408586f9135SBarry Smith 
1409586f9135SBarry Smith    Synopsis:
1410586f9135SBarry Smith    #include <petscsys.h>
1411586f9135SBarry Smith    void PetscStackPopNoCheck(char *funct);
1412586f9135SBarry Smith 
14137cdbe19fSJose E. Roman    Not Collective
14147cdbe19fSJose E. Roman 
1415586f9135SBarry Smith    Input Parameter:
1416586f9135SBarry Smith .   funct - the function name
1417586f9135SBarry Smith 
1418586f9135SBarry Smith    Level: developer
1419586f9135SBarry Smith 
1420586f9135SBarry Smith    Notes:
142187497f52SBarry Smith    Using `PetscCall()`, `PetscCallExternal()`, `PetscCallBack()` and friends negates the need to call this
142287497f52SBarry Smith 
1423586f9135SBarry 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
1424586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1425586f9135SBarry Smith    help debug the problem.
1426586f9135SBarry Smith 
142795bd0b28SBarry Smith    The default stack is a global variable called `petscstack`.
1428586f9135SBarry Smith 
1429586f9135SBarry Smith    Developer Note:
1430586f9135SBarry Smith    `PetscStackPopNoCheck()` takes a function argument while  `PetscStackPop` does not, this difference is likely just historical.
1431586f9135SBarry Smith 
1432586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1433586f9135SBarry Smith M*/
14349371c9d4SSatish Balay   #define PetscStackPopNoCheck(funct) \
14359371c9d4SSatish Balay     do { \
1436362febeeSStefano Zampini       PetscStackSAWsTakeAccess(); \
14375a96b57dSJacob Faibussowitsch       PetscStackPop_Private(petscstack, funct); \
1438362febeeSStefano Zampini       PetscStackSAWsGrantAccess(); \
1439362febeeSStefano Zampini     } while (0)
1440362febeeSStefano Zampini 
14419371c9d4SSatish Balay   #define PetscStackClearTop \
14429371c9d4SSatish Balay     do { \
1443e04113cfSBarry Smith       PetscStackSAWsTakeAccess(); \
14449371c9d4SSatish Balay       if (petscstack.currentsize > 0 && --petscstack.currentsize < PETSCSTACKSIZE) { \
144527104ee2SJacob Faibussowitsch         petscstack.function[petscstack.currentsize]     = PETSC_NULLPTR; \
144627104ee2SJacob Faibussowitsch         petscstack.file[petscstack.currentsize]         = PETSC_NULLPTR; \
144727104ee2SJacob Faibussowitsch         petscstack.line[petscstack.currentsize]         = 0; \
144827104ee2SJacob Faibussowitsch         petscstack.petscroutine[petscstack.currentsize] = 0; \
1449441dd030SJed Brown       } \
145027104ee2SJacob Faibussowitsch       petscstack.hotdepth = PetscMax(petscstack.hotdepth - 1, 0); \
1451e04113cfSBarry Smith       PetscStackSAWsGrantAccess(); \
1452441dd030SJed Brown     } while (0)
1453441dd030SJed Brown 
145430de9b25SBarry Smith   /*MC
14551957e957SBarry Smith    PetscFunctionBegin - First executable line of each PETSc function,  used for error handling. Final
145687497f52SBarry Smith    line of PETSc functions should be `PetscFunctionReturn`(0);
145730de9b25SBarry Smith 
145830de9b25SBarry Smith    Synopsis:
1459aaa7dc30SBarry Smith    #include <petscsys.h>
146030de9b25SBarry Smith    void PetscFunctionBegin;
146130de9b25SBarry Smith 
146220f4b53cSBarry Smith    Not Collective; No Fortran Support
1463eca87e8dSBarry Smith 
146430de9b25SBarry Smith    Usage:
146530de9b25SBarry Smith .vb
146630de9b25SBarry Smith      int something;
146730de9b25SBarry Smith 
146830de9b25SBarry Smith      PetscFunctionBegin;
146930de9b25SBarry Smith .ve
147030de9b25SBarry Smith 
147130de9b25SBarry Smith    Level: developer
147230de9b25SBarry Smith 
147320f4b53cSBarry Smith    Note:
147420f4b53cSBarry Smith    Use `PetscFunctionBeginUser` for application codes.
147520f4b53cSBarry Smith 
1476586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`
147730de9b25SBarry Smith 
147830de9b25SBarry Smith M*/
14799371c9d4SSatish Balay   #define PetscFunctionBegin \
14809371c9d4SSatish Balay     do { \
1481362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_FALSE); \
1482a2f94806SJed Brown       PetscRegister__FUNCT__(); \
1483a2f94806SJed Brown     } while (0)
1484a2f94806SJed Brown 
1485a2f94806SJed Brown   /*MC
148687497f52SBarry Smith    PetscFunctionBeginHot - Substitute for `PetscFunctionBegin` to be used in functions that are called in
1487a2f94806SJed Brown    performance-critical circumstances.  Use of this function allows for lighter profiling by default.
1488a2f94806SJed Brown 
1489a2f94806SJed Brown    Synopsis:
1490aaa7dc30SBarry Smith    #include <petscsys.h>
1491a2f94806SJed Brown    void PetscFunctionBeginHot;
1492a2f94806SJed Brown 
149320f4b53cSBarry Smith    Not Collective; No Fortran Support
1494a2f94806SJed Brown 
1495a2f94806SJed Brown    Usage:
1496a2f94806SJed Brown .vb
1497a2f94806SJed Brown      int something;
1498a2f94806SJed Brown 
1499a2f94806SJed Brown      PetscFunctionBeginHot;
1500a2f94806SJed Brown .ve
1501a2f94806SJed Brown 
1502a2f94806SJed Brown    Level: developer
1503a2f94806SJed Brown 
1504586f9135SBarry Smith .seealso: `PetscFunctionBegin`, `PetscFunctionReturn()`, `PetscStackPushNoCheck()`
1505a2f94806SJed Brown 
1506a2f94806SJed Brown M*/
15079371c9d4SSatish Balay   #define PetscFunctionBeginHot \
15089371c9d4SSatish Balay     do { \
1509362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_TRUE); \
15102d53ad75SBarry Smith       PetscRegister__FUNCT__(); \
151153c77d0aSJed Brown     } while (0)
151253c77d0aSJed Brown 
1513a8d2bbe5SBarry Smith   /*MC
1514530d85c1SBarry Smith    PetscFunctionBeginUser - First executable line of user provided routines
1515a8d2bbe5SBarry Smith 
1516a8d2bbe5SBarry Smith    Synopsis:
1517aaa7dc30SBarry Smith    #include <petscsys.h>
1518a8d2bbe5SBarry Smith    void PetscFunctionBeginUser;
1519a8d2bbe5SBarry Smith 
152020f4b53cSBarry Smith    Not Collective; No Fortran Support
1521a8d2bbe5SBarry Smith 
1522a8d2bbe5SBarry Smith    Usage:
1523a8d2bbe5SBarry Smith .vb
1524a8d2bbe5SBarry Smith      int something;
1525a8d2bbe5SBarry Smith 
1526ac285190SBarry Smith      PetscFunctionBeginUser;
1527a8d2bbe5SBarry Smith .ve
1528a8d2bbe5SBarry Smith 
152920f4b53cSBarry Smith    Level: intermediate
153020f4b53cSBarry Smith 
1531a8d2bbe5SBarry Smith    Notes:
1532530d85c1SBarry Smith    Functions that incorporate this must call `PetscFunctionReturn()` instead of return except for main().
1533530d85c1SBarry Smith 
1534530d85c1SBarry Smith    May be used before `PetscInitialize()`
15351957e957SBarry Smith 
1536530d85c1SBarry Smith    This is identical to `PetscFunctionBegin` except it labels the routine as a user
1537ac285190SBarry Smith    routine instead of as a PETSc library routine.
1538ac285190SBarry Smith 
1539586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, `PetscFunctionBeginHot`, `PetscStackPushNoCheck()`
1540a8d2bbe5SBarry Smith M*/
15419371c9d4SSatish Balay   #define PetscFunctionBeginUser \
15429371c9d4SSatish Balay     do { \
1543362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 2, PETSC_FALSE); \
1544a8d2bbe5SBarry Smith       PetscRegister__FUNCT__(); \
1545a8d2bbe5SBarry Smith     } while (0)
1546a8d2bbe5SBarry Smith 
1547586f9135SBarry Smith   /*MC
1548586f9135SBarry Smith    PetscStackPush - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is
1549586f9135SBarry Smith    currently in the source code and verifies the memory is not corrupted.
1550586f9135SBarry Smith 
1551586f9135SBarry Smith    Synopsis:
1552586f9135SBarry Smith    #include <petscsys.h>
1553586f9135SBarry Smith    void PetscStackPush(char *funct)
1554586f9135SBarry Smith 
15557cdbe19fSJose E. Roman    Not Collective
15567cdbe19fSJose E. Roman 
1557586f9135SBarry Smith    Input Parameter:
1558586f9135SBarry Smith .  funct - the function name
1559586f9135SBarry Smith 
1560586f9135SBarry Smith    Level: developer
1561586f9135SBarry Smith 
1562586f9135SBarry Smith    Notes:
1563586f9135SBarry 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
1564586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1565586f9135SBarry Smith    help debug the problem.
1566586f9135SBarry Smith 
156795bd0b28SBarry Smith    The default stack is a global variable called `petscstack`.
1568586f9135SBarry Smith 
1569586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1570586f9135SBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1571586f9135SBarry Smith M*/
15729371c9d4SSatish Balay   #define PetscStackPush(n) \
15739371c9d4SSatish Balay     do { \
1574362febeeSStefano Zampini       PetscStackPushNoCheck(n, 0, PETSC_FALSE); \
157515681b3cSBarry Smith       CHKMEMQ; \
157615681b3cSBarry Smith     } while (0)
15773a40ed3dSBarry Smith 
1578586f9135SBarry Smith   /*MC
1579586f9135SBarry Smith    PetscStackPop - Pops a function name from the PETSc default stack that tracks where the running program is
1580586f9135SBarry Smith    currently in the source code and verifies the memory is not corrupted.
1581586f9135SBarry Smith 
1582586f9135SBarry Smith    Synopsis:
1583586f9135SBarry Smith    #include <petscsys.h>
1584586f9135SBarry Smith    void PetscStackPop
1585586f9135SBarry Smith 
15867cdbe19fSJose E. Roman    Not Collective
15877cdbe19fSJose E. Roman 
1588586f9135SBarry Smith    Level: developer
1589586f9135SBarry Smith 
1590586f9135SBarry Smith    Notes:
1591586f9135SBarry 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
1592586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1593586f9135SBarry Smith    help debug the problem.
1594586f9135SBarry Smith 
159595bd0b28SBarry Smith    The default stack is a global variable called `petscstack`.
1596586f9135SBarry Smith 
1597586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPopNoCheck()`, `PetscStackPush()`
1598586f9135SBarry Smith M*/
15999371c9d4SSatish Balay   #define PetscStackPop \
16009371c9d4SSatish Balay     do { \
1601441dd030SJed Brown       CHKMEMQ; \
1602362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
160315681b3cSBarry Smith     } while (0)
1604d64ed03dSBarry Smith 
160530de9b25SBarry Smith   /*MC
16060a57284eSJacob Faibussowitsch    PetscFunctionReturn - Last executable line of each PETSc function used for error
16070a57284eSJacob Faibussowitsch    handling. Replaces `return()`.
160830de9b25SBarry Smith 
160930de9b25SBarry Smith    Synopsis:
16100a57284eSJacob Faibussowitsch    #include <petscerror.h>
16110a57284eSJacob Faibussowitsch    void PetscFunctionReturn(...)
161230de9b25SBarry Smith 
1613667f096bSBarry Smith    Not Collective; No Fortran Support
1614eca87e8dSBarry Smith 
16155687f061SJacob Faibussowitsch    Level: beginner
16165687f061SJacob Faibussowitsch 
16170a57284eSJacob Faibussowitsch    Notes:
16180a57284eSJacob Faibussowitsch    This routine is a macro, so while it does not "return" anything itself, it does return from
16190a57284eSJacob Faibussowitsch    the function in the literal sense.
16200a57284eSJacob Faibussowitsch 
16210a57284eSJacob Faibussowitsch    Usually the return value is the integer literal `0` (for example in any function returning
16220a57284eSJacob Faibussowitsch    `PetscErrorCode`), however it is possible to return any arbitrary type. The arguments of
16230a57284eSJacob Faibussowitsch    this macro are placed before the `return` statement as-is.
16240a57284eSJacob Faibussowitsch 
16250a57284eSJacob Faibussowitsch    Any routine which returns via `PetscFunctionReturn()` must begin with a corresponding
16260a57284eSJacob Faibussowitsch    `PetscFunctionBegin`.
16270a57284eSJacob Faibussowitsch 
16280a57284eSJacob Faibussowitsch    For routines which return `void` use `PetscFunctionReturnVoid()` instead.
16290a57284eSJacob Faibussowitsch 
16300a57284eSJacob Faibussowitsch    Example Usage:
163130de9b25SBarry Smith .vb
16320a57284eSJacob Faibussowitsch    PetscErrorCode foo(int *x)
16330a57284eSJacob Faibussowitsch    {
16340a57284eSJacob Faibussowitsch      PetscFunctionBegin; // don't forget the begin!
16350a57284eSJacob Faibussowitsch      *x = 10;
16363ba16761SJacob Faibussowitsch      PetscFunctionReturn(PETSC_SUCCESS);
163730de9b25SBarry Smith    }
163830de9b25SBarry Smith .ve
163930de9b25SBarry Smith 
16400a57284eSJacob Faibussowitsch    May return any arbitrary type\:
16410a57284eSJacob Faibussowitsch .vb
16420a57284eSJacob Faibussowitsch   struct Foo
16430a57284eSJacob Faibussowitsch   {
16440a57284eSJacob Faibussowitsch     int x;
16450a57284eSJacob Faibussowitsch   };
16460a57284eSJacob Faibussowitsch 
16470a57284eSJacob Faibussowitsch   struct Foo make_foo(int value)
16480a57284eSJacob Faibussowitsch   {
16490a57284eSJacob Faibussowitsch     struct Foo f;
16500a57284eSJacob Faibussowitsch 
16510a57284eSJacob Faibussowitsch     PetscFunctionBegin;
16520a57284eSJacob Faibussowitsch     f.x = value;
16530a57284eSJacob Faibussowitsch     PetscFunctionReturn(f);
16540a57284eSJacob Faibussowitsch   }
16550a57284eSJacob Faibussowitsch .ve
16560a57284eSJacob Faibussowitsch 
16570a57284eSJacob Faibussowitsch .seealso: `PetscFunctionBegin`, `PetscFunctionBeginUser`, `PetscFunctionReturnVoid()`,
16580a57284eSJacob Faibussowitsch           `PetscStackPopNoCheck()`
165930de9b25SBarry Smith M*/
16605687f061SJacob Faibussowitsch   #define PetscFunctionReturn(...) \
16619371c9d4SSatish Balay     do { \
1662362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
16635687f061SJacob Faibussowitsch       return __VA_ARGS__; \
166427104ee2SJacob Faibussowitsch     } while (0)
1665d64ed03dSBarry Smith 
16660a57284eSJacob Faibussowitsch   /*MC
16670a57284eSJacob Faibussowitsch   PetscFunctionReturnVoid - Like `PetscFunctionReturn()` but returns `void`
16680a57284eSJacob Faibussowitsch 
16690a57284eSJacob Faibussowitsch   Synopsis:
16700a57284eSJacob Faibussowitsch   #include <petscerror.h>
16710a57284eSJacob Faibussowitsch   void PetscFunctionReturnVoid()
16720a57284eSJacob Faibussowitsch 
16730a57284eSJacob Faibussowitsch   Not Collective
16740a57284eSJacob Faibussowitsch 
16755687f061SJacob Faibussowitsch   Level: beginner
16765687f061SJacob Faibussowitsch 
16775687f061SJacob Faibussowitsch   Note:
16780a57284eSJacob Faibussowitsch   Behaves identically to `PetscFunctionReturn()` except that it returns `void`. That is, this
16795687f061SJacob Faibussowitsch   macro culminates with `return`.
16800a57284eSJacob Faibussowitsch 
16810a57284eSJacob Faibussowitsch   Example Usage:
16820a57284eSJacob Faibussowitsch .vb
16830a57284eSJacob Faibussowitsch   void foo()
16840a57284eSJacob Faibussowitsch   {
16850a57284eSJacob Faibussowitsch     PetscFunctionBegin; // must start with PetscFunctionBegin!
16860a57284eSJacob Faibussowitsch     bar();
16870a57284eSJacob Faibussowitsch     baz();
16880a57284eSJacob Faibussowitsch     PetscFunctionReturnVoid();
16890a57284eSJacob Faibussowitsch   }
16900a57284eSJacob Faibussowitsch .ve
16910a57284eSJacob Faibussowitsch 
16920a57284eSJacob Faibussowitsch .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, PetscFunctionBeginUser`
16930a57284eSJacob Faibussowitsch M*/
16949371c9d4SSatish Balay   #define PetscFunctionReturnVoid() \
16959371c9d4SSatish Balay     do { \
1696362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
169727104ee2SJacob Faibussowitsch       return; \
169827104ee2SJacob Faibussowitsch     } while (0)
169927104ee2SJacob Faibussowitsch #else /* PETSC_USE_DEBUG */
170027104ee2SJacob Faibussowitsch   #define PetscStackPushNoCheck(funct, petsc_routine, hot)
170137154d10SBarry Smith   #define PetscStackUpdateLine
1702792fecdfSBarry Smith   #define PetscStackPushExternal(funct)
17033ba16761SJacob Faibussowitsch   #define PetscStackPopNoCheck(...)
170427104ee2SJacob Faibussowitsch   #define PetscStackClearTop
17053a40ed3dSBarry Smith   #define PetscFunctionBegin
17060bdf7c52SPeter Brune   #define PetscFunctionBeginUser
1707a2f94806SJed Brown   #define PetscFunctionBeginHot
17080a57284eSJacob Faibussowitsch   #define PetscFunctionReturn(...)  return __VA_ARGS__
17095665465eSBarry Smith   #define PetscFunctionReturnVoid() return
1710812af9f3SBarry Smith   #define PetscStackPop             CHKMEMQ
1711812af9f3SBarry Smith   #define PetscStackPush(f)         CHKMEMQ
171227104ee2SJacob Faibussowitsch #endif /* PETSC_USE_DEBUG */
17133a40ed3dSBarry Smith 
17146d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
17153ba16761SJacob Faibussowitsch   #define PetscStackCallExternalVoid(...)
17163ba16761SJacob Faibussowitsch template <typename F, typename... Args>
17173ba16761SJacob Faibussowitsch void PetscCallExternal(F, Args...);
17186d210af2SJacob Faibussowitsch #else
1719586f9135SBarry Smith   /*MC
1720e77caa6dSBarry Smith     PetscStackCallExternalVoid - Calls an external library routine or user function after pushing the name of the routine on the stack.
1721eb6b5d47SBarry Smith 
1722eb6b5d47SBarry Smith    Input Parameters:
1723eb6b5d47SBarry Smith +   name    - string that gives the name of the function being called
1724586f9135SBarry Smith -   routine - actual call to the routine, for example, functionname(a,b)
1725fd3f9acdSBarry Smith 
1726586f9135SBarry Smith    Level: developer
1727eb6b5d47SBarry Smith 
172895bd0b28SBarry Smith    Notes:
1729792fecdfSBarry Smith    Often one should use `PetscCallExternal()` instead. This routine is intended for external library routines that DO NOT return error codes
1730eb6b5d47SBarry Smith 
1731586f9135SBarry Smith    In debug mode this also checks the memory for corruption at the end of the function call.
1732586f9135SBarry Smith 
173395bd0b28SBarry Smith    Certain external packages, such as BLAS/LAPACK may have their own macros, `PetscCallBLAS()` for managing the call, error checking, etc.
1734586f9135SBarry Smith 
1735586f9135SBarry Smith    Developer Note:
1736586f9135SBarry 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.
1737586f9135SBarry Smith 
1738792fecdfSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscCallExternal()`, `PetscCallBLAS()`
1739586f9135SBarry Smith @*/
17403ba16761SJacob Faibussowitsch   #define PetscStackCallExternalVoid(name, ...) \
17419371c9d4SSatish Balay     do { \
174228770333SStefano Zampini       PetscStackPushExternal(name); \
17433ba16761SJacob Faibussowitsch       __VA_ARGS__; \
17449371c9d4SSatish Balay       PetscStackPop; \
17459371c9d4SSatish Balay     } while (0)
1746eb6b5d47SBarry Smith 
1747586f9135SBarry Smith   /*MC
1748792fecdfSBarry Smith     PetscCallExternal - Calls an external library routine that returns an error code after pushing the name of the routine on the stack.
1749fd3f9acdSBarry Smith 
1750fd3f9acdSBarry Smith    Input Parameters:
1751fd3f9acdSBarry Smith +  func - name of the routine
1752586f9135SBarry Smith -  args - arguments to the routine
1753586f9135SBarry Smith 
1754586f9135SBarry Smith    Level: developer
1755fd3f9acdSBarry Smith 
175695452b02SPatrick Sanan    Notes:
1757e77caa6dSBarry Smith    This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not.
1758dbf62e16SBarry Smith 
1759586f9135SBarry Smith    In debug mode this also checks the memory for corruption at the end of the function call.
1760fd3f9acdSBarry Smith 
176187497f52SBarry Smith    Assumes the error return code of the function is an integer and that a value of 0 indicates success
176287497f52SBarry Smith 
1763586f9135SBarry Smith    Developer Note:
1764d5b43468SJose E. Roman    This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc.
1765586f9135SBarry Smith 
1766e77caa6dSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()`
1767586f9135SBarry Smith M*/
17689371c9d4SSatish Balay   #define PetscCallExternal(func, ...) \
17699371c9d4SSatish Balay     do { \
1770a74df02fSJacob Faibussowitsch       PetscStackPush(PetscStringize(func)); \
17713ba16761SJacob Faibussowitsch       int ierr_petsc_call_external_ = func(__VA_ARGS__); \
17721d4906efSStefano Zampini       PetscStackPop; \
17733ba16761SJacob 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_); \
1774fd3f9acdSBarry Smith     } while (0)
17756d210af2SJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */
1776