xref: /petsc/include/petscerror.h (revision a496304597bacff3545e802853d69e8765312868)
154a8ef01SBarry Smith /*
2f621e05eSBarry Smith     Contains all error handling interfaces for PETSc.
354a8ef01SBarry Smith */
4*a4963045SJacob 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:
3587497f52SBarry Smith +  comm - A 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()`,
52db781477SPatrick Sanan           `PetscError()`, `PetscCall()`, `CHKMEMQ`, `CHKERRA()`, `PetscCallMPI()`
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:
8087497f52SBarry Smith +  comm - A 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 
86986eef2eSBarry Smith    Notes:
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 
90db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscCall()`, `PetscCallMPI()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`
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:
104f388eb8bSPatrick Sanan +  comm - A 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 
116db781477SPatrick Sanan .seealso: `SETERRQ()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()`
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:
129fa190f98SMatthew G. Knepley +  comm - A 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 
146db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `PetscCall()`, `CHKMEMQ`
147fa190f98SMatthew G. Knepley M*/
1489371c9d4SSatish Balay #define SETERRABORT(comm, ierr, ...) \
1499371c9d4SSatish Balay   do { \
1503ba16761SJacob Faibussowitsch     (void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \
15198921bdaSJacob Faibussowitsch     MPI_Abort(comm, ierr); \
15298921bdaSJacob Faibussowitsch   } while (0)
1539a00fa46SSatish Balay 
15430de9b25SBarry Smith /*MC
1552c71b3e2SJacob Faibussowitsch   PetscCheck - Check that a particular condition is true
1562c71b3e2SJacob Faibussowitsch 
1572c71b3e2SJacob Faibussowitsch   Synopsis:
1582c71b3e2SJacob Faibussowitsch   #include <petscerror.h>
1592c71b3e2SJacob Faibussowitsch   void PetscCheck(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
1602c71b3e2SJacob Faibussowitsch 
1617cdbe19fSJose E. Roman   Collective; No Fortran Support
1622c71b3e2SJacob Faibussowitsch 
1632c71b3e2SJacob Faibussowitsch   Input Parameters:
1642c71b3e2SJacob Faibussowitsch + cond    - The boolean condition
1652c71b3e2SJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
1662c71b3e2SJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
1672c71b3e2SJacob Faibussowitsch - message - Error message in printf format
1682c71b3e2SJacob Faibussowitsch 
169667f096bSBarry Smith   Level: beginner
170667f096bSBarry Smith 
1712c71b3e2SJacob Faibussowitsch   Notes:
1722c71b3e2SJacob Faibussowitsch   Enabled in both optimized and debug builds.
1732c71b3e2SJacob Faibussowitsch 
17487497f52SBarry Smith   Calls `SETERRQ()` if the assertion fails, so can only be called from functions returning a
17587497f52SBarry Smith   `PetscErrorCode` (or equivalent type after conversion).
1762c71b3e2SJacob Faibussowitsch 
1774be741a6SBarry Smith  .seealso: `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheckAbort()`
1789566063dSJacob Faibussowitsch M*/
1799371c9d4SSatish Balay #define PetscCheck(cond, comm, ierr, ...) \
1809371c9d4SSatish Balay   do { \
1819371c9d4SSatish Balay     if (PetscUnlikely(!(cond))) SETERRQ(comm, ierr, __VA_ARGS__); \
1829371c9d4SSatish Balay   } while (0)
1832c71b3e2SJacob Faibussowitsch 
1842c71b3e2SJacob Faibussowitsch /*MC
1854be741a6SBarry Smith   PetscCheckAbort - Check that a particular condition is true, otherwise prints error and aborts
1864be741a6SBarry Smith 
1874be741a6SBarry Smith   Synopsis:
1884be741a6SBarry Smith   #include <petscerror.h>
1894be741a6SBarry Smith   void PetscCheckAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
1904be741a6SBarry Smith 
1917cdbe19fSJose E. Roman   Collective; No Fortran Support
1924be741a6SBarry Smith 
1934be741a6SBarry Smith   Input Parameters:
1944be741a6SBarry Smith + cond    - The boolean condition
1954be741a6SBarry Smith . comm    - The communicator on which the check can be collective on
1964be741a6SBarry Smith . ierr    - A nonzero error code, see include/petscerror.h for the complete list
1974be741a6SBarry Smith - message - Error message in printf format
1984be741a6SBarry Smith 
199667f096bSBarry Smith   Level: developer
200667f096bSBarry Smith 
2014be741a6SBarry Smith   Notes:
2024be741a6SBarry Smith   Enabled in both optimized and debug builds.
2034be741a6SBarry Smith 
20487497f52SBarry Smith   Calls `SETERRABORT()` if the assertion fails, can be called from a function that does not return an
20587497f52SBarry Smith   error code, such as a C++ constructor. usually `PetscCheck()` should be used.
2064be741a6SBarry Smith 
207989c49a2SBarry Smith .seealso: `PetscAssertAbort()`, `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheck()`, `SETERRABORT()`
2084be741a6SBarry Smith M*/
2099371c9d4SSatish Balay #define PetscCheckAbort(cond, comm, ierr, ...) \
2100e6b6b59SJacob Faibussowitsch   do { \
2110e6b6b59SJacob Faibussowitsch     if (PetscUnlikely(!(cond))) SETERRABORT(comm, ierr, __VA_ARGS__); \
212c7481402SJacob Faibussowitsch   } while (0)
2134be741a6SBarry Smith 
2144be741a6SBarry Smith /*MC
2159ace16cdSJacob Faibussowitsch   PetscAssert - Assert that a particular condition is true
2169ace16cdSJacob Faibussowitsch 
2179ace16cdSJacob Faibussowitsch   Synopsis:
2189ace16cdSJacob Faibussowitsch   #include <petscerror.h>
2199ace16cdSJacob Faibussowitsch   void PetscAssert(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
2209ace16cdSJacob Faibussowitsch 
2217cdbe19fSJose E. Roman   Collective; No Fortran Support
2229ace16cdSJacob Faibussowitsch 
2239ace16cdSJacob Faibussowitsch   Input Parameters:
2249ace16cdSJacob Faibussowitsch + cond    - The boolean condition
2259ace16cdSJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
2269ace16cdSJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
2279ace16cdSJacob Faibussowitsch - message - Error message in printf format
2289ace16cdSJacob Faibussowitsch 
229667f096bSBarry Smith   Level: beginner
230667f096bSBarry Smith 
2319ace16cdSJacob Faibussowitsch   Notes:
232c7481402SJacob Faibussowitsch   Equivalent to `PetscCheck()` if debugging is enabled, and `PetscAssume(cond)` otherwise.
2332c71b3e2SJacob Faibussowitsch 
23487497f52SBarry Smith   See `PetscCheck()` for usage and behaviour.
23587497f52SBarry Smith 
23687497f52SBarry Smith   This is needed instead of simply using `assert()` because this correctly handles the collective nature of errors under MPI
2379ace16cdSJacob Faibussowitsch 
2380e6b6b59SJacob Faibussowitsch .seealso: `PetscCheck()`, `SETERRQ()`, `PetscError()`, `PetscAssertAbort()`
2399566063dSJacob Faibussowitsch M*/
240c7481402SJacob Faibussowitsch #if PetscDefined(USE_DEBUG)
241c7481402SJacob Faibussowitsch   #define PetscAssert(cond, comm, ierr, ...) PetscCheck(cond, comm, ierr, __VA_ARGS__)
242c7481402SJacob Faibussowitsch #else
243c7481402SJacob Faibussowitsch   #define PetscAssert(cond, ...) PetscAssume(cond)
244c7481402SJacob Faibussowitsch #endif
2459ace16cdSJacob Faibussowitsch 
2469ace16cdSJacob Faibussowitsch /*MC
2470e6b6b59SJacob Faibussowitsch   PetscAssertAbort - Assert that a particular condition is true, otherwise prints error and aborts
2480e6b6b59SJacob Faibussowitsch 
2490e6b6b59SJacob Faibussowitsch   Synopsis:
2500e6b6b59SJacob Faibussowitsch   #include <petscerror.h>
2510e6b6b59SJacob Faibussowitsch   void PetscAssertAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
2520e6b6b59SJacob Faibussowitsch 
2537cdbe19fSJose E. Roman   Collective; No Fortran Support
2540e6b6b59SJacob Faibussowitsch 
2550e6b6b59SJacob Faibussowitsch   Input Parameters:
2560e6b6b59SJacob Faibussowitsch + cond    - The boolean condition
2570e6b6b59SJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
2580e6b6b59SJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
2590e6b6b59SJacob Faibussowitsch - message - Error message in printf format
2600e6b6b59SJacob Faibussowitsch 
261667f096bSBarry Smith   Level: beginner
262667f096bSBarry Smith 
2630e6b6b59SJacob Faibussowitsch   Notes:
2640e6b6b59SJacob Faibussowitsch   Enabled only in debug builds. See `PetscCheckAbort()` for usage.
2650e6b6b59SJacob Faibussowitsch 
2660e6b6b59SJacob Faibussowitsch .seealso: `PetscCheckAbort()`, `PetscAssert()`, `PetscCheck()`, `SETERRABORT()`, `PetscError()`
2670e6b6b59SJacob Faibussowitsch M*/
2683ba16761SJacob Faibussowitsch #if PetscDefined(USE_DEBUG)
2693ba16761SJacob Faibussowitsch   #define PetscAssertAbort(cond, comm, ierr, ...) PetscCheckAbort(cond, comm, ierr, __VA_ARGS__)
2703ba16761SJacob Faibussowitsch #else
2713ba16761SJacob Faibussowitsch   #define PetscAssertAbort(cond, comm, ierr, ...) PetscAssume(cond)
2723ba16761SJacob Faibussowitsch #endif
2730e6b6b59SJacob Faibussowitsch 
2740e6b6b59SJacob Faibussowitsch /*MC
2753ba16761SJacob Faibussowitsch   PetscCall - Calls a PETSc function and then checks the resulting error code, if it is
2763ba16761SJacob Faibussowitsch   non-zero it calls the error handler and returns from the current function with the error
2773ba16761SJacob Faibussowitsch   code.
2789566063dSJacob Faibussowitsch 
2799566063dSJacob Faibussowitsch   Synopsis:
2809566063dSJacob Faibussowitsch   #include <petscerror.h>
28149c86fc7SBarry Smith   void PetscCall(PetscFunction(args))
2829566063dSJacob Faibussowitsch 
2839566063dSJacob Faibussowitsch   Not Collective
2849566063dSJacob Faibussowitsch 
2859566063dSJacob Faibussowitsch   Input Parameter:
28649c86fc7SBarry Smith . PetscFunction - any PETSc function that returns an error code
2879566063dSJacob Faibussowitsch 
288667f096bSBarry Smith   Level: beginner
289667f096bSBarry Smith 
2909566063dSJacob Faibussowitsch   Notes:
2919566063dSJacob Faibussowitsch   Once the error handler is called the calling function is then returned from with the given
29287497f52SBarry Smith   error code. Experienced users can set the error handler with `PetscPushErrorHandler()`.
2939566063dSJacob Faibussowitsch 
29487497f52SBarry Smith   `PetscCall()` cannot be used in functions returning a datatype not convertible to
29587497f52SBarry Smith   `PetscErrorCode`. For example, `PetscCall()` may not be used in functions returning void, use
2963ba16761SJacob Faibussowitsch   `PetscCallAbort()` or `PetscCallVoid()` in this case.
2979566063dSJacob Faibussowitsch 
2989566063dSJacob Faibussowitsch   Example Usage:
2999566063dSJacob Faibussowitsch .vb
3009566063dSJacob Faibussowitsch   PetscCall(PetscInitiailize(...)); // OK to call even when PETSc is not yet initialized!
3019566063dSJacob Faibussowitsch 
3029566063dSJacob Faibussowitsch   struct my_struct
3039566063dSJacob Faibussowitsch   {
3049566063dSJacob Faibussowitsch     void *data;
3059566063dSJacob Faibussowitsch   } my_complex_type;
3069566063dSJacob Faibussowitsch 
3079566063dSJacob Faibussowitsch   struct my_struct bar(void)
3089566063dSJacob Faibussowitsch   {
3096aad120cSJose E. Roman     PetscCall(foo(15)); // ERROR PetscErrorCode not convertible to struct my_struct!
3109566063dSJacob Faibussowitsch   }
3119566063dSJacob Faibussowitsch 
3126aad120cSJose E. Roman   PetscCall(bar()) // ERROR input not convertible to PetscErrorCode
3139566063dSJacob Faibussowitsch .ve
3149566063dSJacob Faibussowitsch 
31587497f52SBarry Smith   It is also possible to call this directly on a `PetscErrorCode` variable
31649c86fc7SBarry Smith .vb
31749c86fc7SBarry Smith   PetscCall(ierr);  // check if ierr is nonzero
31849c86fc7SBarry Smith .ve
31949c86fc7SBarry Smith 
320792fecdfSBarry Smith   Should not be used to call callback functions provided by users, `PetscCallBack()` should be used in that situation.
321ef1023bdSBarry Smith 
3226a8be23eSBarry Smith   `PetscUseTypeMethod()` or `PetscTryTypeMethod()` should be used when calling functions pointers contained in a PETSc object's `ops` array
3236a8be23eSBarry Smith 
32449c86fc7SBarry Smith   Fortran Notes:
32549c86fc7SBarry Smith     The Fortran function from which this is used must declare a variable PetscErrorCode ierr and ierr must be
32687497f52SBarry Smith     the final argument to the PETSc function being called.
32749c86fc7SBarry Smith 
32849c86fc7SBarry Smith     In the main program and in Fortran subroutines that do not have ierr as the final return parameter one
32987497f52SBarry Smith     should use `PetscCallA()`
33049c86fc7SBarry Smith 
33149c86fc7SBarry Smith   Example Fortran Usage:
33249c86fc7SBarry Smith .vb
33349c86fc7SBarry Smith   PetscErrorCode ierr
33449c86fc7SBarry Smith   Vec v
33549c86fc7SBarry Smith 
33649c86fc7SBarry Smith   ...
33749c86fc7SBarry Smith   PetscCall(VecShift(v,1.0,ierr))
33849c86fc7SBarry Smith   PetscCallA(VecShift(v,1.0,ierr))
33949c86fc7SBarry Smith .ve
34049c86fc7SBarry Smith 
341667f096bSBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`,
342667f096bSBarry Smith           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`,
3433ba16761SJacob Faibussowitsch           `CHKERRMPI()`, `PetscCallBack()`, `PetscCallAbort()`, `PetscCallVoid()`
3449566063dSJacob Faibussowitsch M*/
345ef1023bdSBarry Smith 
346ef1023bdSBarry Smith /*MC
347989c49a2SBarry Smith    PetscCallA - Fortran-only macro that should be used in the main program to call PETSc functions instead of using
348989c49a2SBarry Smith    PetscCall() which should be used in other Fortran subroutines
349989c49a2SBarry Smith 
350989c49a2SBarry Smith    Synopsis:
351989c49a2SBarry Smith    #include <petscsys.h>
352989c49a2SBarry Smith    PetscErrorCode PetscCallA(PetscFunction(arguments,ierr))
353989c49a2SBarry Smith 
354989c49a2SBarry Smith    Collective
355989c49a2SBarry Smith 
356989c49a2SBarry Smith    Input Parameter:
357989c49a2SBarry Smith .  PetscFunction(arguments,ierr) - the call to the function
358989c49a2SBarry Smith 
359989c49a2SBarry Smith   Level: beginner
360989c49a2SBarry Smith 
361989c49a2SBarry Smith    Notes:
362989c49a2SBarry Smith    This should only be used with Fortran. With C/C++, use `PetscCall()` always.
363989c49a2SBarry Smith 
364989c49a2SBarry Smith    Use `SETERRA()` to set an error in a Fortran main program and `SETERRQ()` in Fortran subroutines
365989c49a2SBarry Smith 
366989c49a2SBarry Smith .seealso: `SETERRQ()`, `SETERRA()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()`
367989c49a2SBarry Smith M*/
368989c49a2SBarry Smith 
369989c49a2SBarry Smith /*MC
370792fecdfSBarry 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
371ef1023bdSBarry Smith   handler and returns from the current function with the error code.
372ef1023bdSBarry Smith 
373ef1023bdSBarry Smith   Synopsis:
374ef1023bdSBarry Smith   #include <petscerror.h>
375792fecdfSBarry Smith   void PetscCallBack(const char *functionname,PetscFunction(args))
376ef1023bdSBarry Smith 
3777cdbe19fSJose E. Roman   Not Collective; No Fortran Support
378ef1023bdSBarry Smith 
379ef1023bdSBarry Smith   Input Parameters:
380ef1023bdSBarry Smith + functionname - the name of the function being called, this can be a string with spaces that describes the meaning of the callback
38187497f52SBarry Smith - PetscFunction - user provided callback function that returns an error code
382ef1023bdSBarry Smith 
383ef1023bdSBarry Smith   Example Usage:
384ef1023bdSBarry Smith .vb
385792fecdfSBarry Smith   PetscCallBack("XXX callback to do something",a->callback(...));
386ef1023bdSBarry Smith .ve
387ef1023bdSBarry Smith 
388ef1023bdSBarry Smith   Level: developer
389ef1023bdSBarry Smith 
390667f096bSBarry Smith   Notes:
391667f096bSBarry Smith   Once the error handler is called the calling function is then returned from with the given
392667f096bSBarry Smith   error code. Experienced users can set the error handler with `PetscPushErrorHandler()`.
393667f096bSBarry Smith 
394667f096bSBarry Smith   `PetscCallBack()` should only be called in PETSc when a call is being made to a user provided call-back routine.
395667f096bSBarry Smith 
39687497f52SBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscCall()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`
397ef1023bdSBarry Smith           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, `CHKERRMPI()`, `PetscCall()`
398ef1023bdSBarry Smith M*/
399ef1023bdSBarry Smith 
4003ba16761SJacob Faibussowitsch /*MC
4013ba16761SJacob Faibussowitsch   PetscCallVoid - Like `PetscCall()` but for functions returning `void`
4023ba16761SJacob Faibussowitsch 
4033ba16761SJacob Faibussowitsch   Synopsis:
4043ba16761SJacob Faibussowitsch   #include <petscerror.h>
4053ba16761SJacob Faibussowitsch   void PetscCall(PetscFunction(args))
4063ba16761SJacob Faibussowitsch 
4077cdbe19fSJose E. Roman   Not Collective; No Fortran Support
4083ba16761SJacob Faibussowitsch 
4093ba16761SJacob Faibussowitsch   Input Parameter:
4103ba16761SJacob Faibussowitsch . PetscFunction - any PETSc function that returns an error code
4113ba16761SJacob Faibussowitsch 
4123ba16761SJacob Faibussowitsch   Example Usage:
4133ba16761SJacob Faibussowitsch .vb
4143ba16761SJacob Faibussowitsch   void foo()
4153ba16761SJacob Faibussowitsch   {
4163ba16761SJacob Faibussowitsch     KSP ksp;
4173ba16761SJacob Faibussowitsch 
4183ba16761SJacob Faibussowitsch     PetscFunctionBeginUser;
4193ba16761SJacob Faibussowitsch     // OK, properly handles PETSc error codes
4203ba16761SJacob Faibussowitsch     PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp));
4213ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
4223ba16761SJacob Faibussowitsch   }
4233ba16761SJacob Faibussowitsch 
4243ba16761SJacob Faibussowitsch   PetscErrorCode bar()
4253ba16761SJacob Faibussowitsch   {
4263ba16761SJacob Faibussowitsch     KSP ksp;
4273ba16761SJacob Faibussowitsch 
4283ba16761SJacob Faibussowitsch     PetscFunctionBeginUser;
4293ba16761SJacob Faibussowitsch     // ERROR, Non-void function 'bar' should return a value
4303ba16761SJacob Faibussowitsch     PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp));
4313ba16761SJacob Faibussowitsch     // OK, returning PetscErrorCode
4323ba16761SJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
4333ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
4343ba16761SJacob Faibussowitsch   }
435667f096bSBarry Smith .ve
4363ba16761SJacob Faibussowitsch 
4373ba16761SJacob Faibussowitsch   Level: beginner
4383ba16761SJacob Faibussowitsch 
439667f096bSBarry Smith   Notes:
440667f096bSBarry Smith   Has identical usage to `PetscCall()`, except that it returns `void` on error instead of a
441667f096bSBarry Smith   `PetscErrorCode`. See `PetscCall()` for more detailed discussion.
442667f096bSBarry Smith 
443667f096bSBarry Smith   Note that users should prefer `PetscCallAbort()` to this routine. While this routine does
444667f096bSBarry Smith   "handle" errors by returning from the enclosing function, it effectively gobbles the
445667f096bSBarry Smith   error. Since the enclosing function itself returns `void`, its callers have no way of knowing
446667f096bSBarry Smith   that the routine returned early due to an error. `PetscCallAbort()` at least ensures that the
447667f096bSBarry Smith   program crashes gracefully.
448667f096bSBarry Smith 
4493ba16761SJacob Faibussowitsch .seealso: `PetscCall()`, `PetscErrorCode`
4503ba16761SJacob Faibussowitsch M*/
4519566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
4529566063dSJacob Faibussowitsch void PetscCall(PetscErrorCode);
453792fecdfSBarry Smith void PetscCallBack(const char *, PetscErrorCode);
4549566063dSJacob Faibussowitsch void PetscCallVoid(PetscErrorCode);
4559566063dSJacob Faibussowitsch #else
4569371c9d4SSatish Balay   #define PetscCall(...) \
4579371c9d4SSatish Balay     do { \
4583ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_q_; \
45937154d10SBarry Smith       PetscStackUpdateLine; \
4603ba16761SJacob Faibussowitsch       ierr_petsc_call_q_ = __VA_ARGS__; \
4613ba16761SJacob 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, " "); \
4629566063dSJacob Faibussowitsch     } while (0)
4639371c9d4SSatish Balay   #define PetscCallBack(function, ...) \
4649371c9d4SSatish Balay     do { \
4653ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_q_; \
466e33ced7fSLisandro Dalcin       PetscStackUpdateLine; \
467e33ced7fSLisandro Dalcin       PetscStackPushExternal(function); \
4683ba16761SJacob Faibussowitsch       ierr_petsc_call_q_ = __VA_ARGS__; \
469ef1023bdSBarry Smith       PetscStackPop; \
4703ba16761SJacob 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, " "); \
471ef1023bdSBarry Smith     } while (0)
4729371c9d4SSatish Balay   #define PetscCallVoid(...) \
4739371c9d4SSatish Balay     do { \
4743ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_void_; \
475e33ced7fSLisandro Dalcin       PetscStackUpdateLine; \
4763ba16761SJacob Faibussowitsch       ierr_petsc_call_void_ = __VA_ARGS__; \
4773ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_void_ != PETSC_SUCCESS)) { \
4783ba16761SJacob Faibussowitsch         ierr_petsc_call_void_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_void_, PETSC_ERROR_REPEAT, " "); \
4793ba16761SJacob Faibussowitsch         (void)ierr_petsc_call_void_; \
4809371c9d4SSatish Balay         return; \
4819371c9d4SSatish Balay       } \
4829566063dSJacob Faibussowitsch     } while (0)
4839566063dSJacob Faibussowitsch #endif
4849566063dSJacob Faibussowitsch 
4859566063dSJacob Faibussowitsch /*MC
4869566063dSJacob Faibussowitsch   CHKERRQ - Checks error code returned from PETSc function
48730de9b25SBarry Smith 
48830de9b25SBarry Smith   Synopsis:
489aaa7dc30SBarry Smith   #include <petscsys.h>
4909566063dSJacob Faibussowitsch   void CHKERRQ(PetscErrorCode ierr)
4919566063dSJacob Faibussowitsch 
4929566063dSJacob Faibussowitsch   Not Collective
4939566063dSJacob Faibussowitsch 
4942fe279fdSBarry Smith   Input Parameter:
4959566063dSJacob Faibussowitsch . ierr - nonzero error code
4969566063dSJacob Faibussowitsch 
4979566063dSJacob Faibussowitsch   Level: deprecated
4989566063dSJacob Faibussowitsch 
499667f096bSBarry Smith   Note:
500667f096bSBarry Smith   Deprecated in favor of `PetscCall()`. This routine behaves identically to it.
501667f096bSBarry Smith 
502db781477SPatrick Sanan .seealso: `PetscCall()`
5039566063dSJacob Faibussowitsch M*/
5049566063dSJacob Faibussowitsch #define CHKERRQ(...) PetscCall(__VA_ARGS__)
5059566063dSJacob Faibussowitsch #define CHKERRV(...) PetscCallVoid(__VA_ARGS__)
5069566063dSJacob Faibussowitsch 
507db9cea48SBarry Smith PETSC_EXTERN void PetscMPIErrorString(PetscMPIInt, char *);
508db9cea48SBarry Smith 
5099566063dSJacob Faibussowitsch /*MC
5109566063dSJacob Faibussowitsch   PetscCallMPI - Checks error code returned from MPI calls, if non-zero it calls the error
5119566063dSJacob Faibussowitsch   handler and then returns
5129566063dSJacob Faibussowitsch 
5139566063dSJacob Faibussowitsch   Synopsis:
5149566063dSJacob Faibussowitsch   #include <petscerror.h>
51549c86fc7SBarry Smith   void PetscCallMPI(MPI_Function(args))
51630de9b25SBarry Smith 
517eca87e8dSBarry Smith   Not Collective
51830de9b25SBarry Smith 
5192fe279fdSBarry Smith   Input Parameter:
52049c86fc7SBarry Smith . MPI_Function - an MPI function that returns an MPI error code
52130de9b25SBarry Smith 
522667f096bSBarry Smith   Level: beginner
523667f096bSBarry Smith 
5249566063dSJacob Faibussowitsch   Notes:
52587497f52SBarry Smith   Always returns the error code `PETSC_ERR_MPI`; the MPI error code and string are embedded in
5269566063dSJacob Faibussowitsch   the string error message. Do not use this to call any other routines (for example PETSc
5273ba16761SJacob Faibussowitsch   routines), it should only be used for direct MPI calls. The user may configure PETSc with the
5283ba16761SJacob Faibussowitsch   `--with-strict-petscerrorcode` option to check this at compile-time, otherwise they must
5299566063dSJacob Faibussowitsch   check this themselves.
5309566063dSJacob Faibussowitsch 
531aaa8cc7dSPierre Jolivet   This routine can only be used in functions returning `PetscErrorCode` themselves. If the
5323ba16761SJacob Faibussowitsch   calling function returns a different type, use `PetscCallMPIAbort()` instead.
5333ba16761SJacob Faibussowitsch 
5349566063dSJacob Faibussowitsch   Example Usage:
5359566063dSJacob Faibussowitsch .vb
5369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(...)); // OK, calling MPI function
5379566063dSJacob Faibussowitsch 
5389566063dSJacob Faibussowitsch   PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead!
5399566063dSJacob Faibussowitsch .ve
5409566063dSJacob Faibussowitsch 
54149c86fc7SBarry Smith   Fortran Notes:
54287497f52SBarry Smith     The Fortran function from which this is used must declare a variable `PetscErrorCode` ierr and ierr must be
54349c86fc7SBarry Smith     the final argument to the MPI function being called.
54449c86fc7SBarry Smith 
54549c86fc7SBarry Smith     In the main program and in Fortran subroutines that do not have ierr as the final return parameter one
54687497f52SBarry Smith     should use `PetscCallMPIA()`
54749c86fc7SBarry Smith 
54849c86fc7SBarry Smith   Fortran Usage:
54949c86fc7SBarry Smith .vb
55049c86fc7SBarry Smith   PetscErrorCode ierr or integer ierr
55149c86fc7SBarry Smith   ...
55249c86fc7SBarry Smith   PetscCallMPI(MPI_Comm_size(...,ierr))
55349c86fc7SBarry Smith   PetscCallMPIA(MPI_Comm_size(...,ierr)) ! Will abort after calling error handler
55449c86fc7SBarry Smith 
55549c86fc7SBarry Smith   PetscCallMPI(MPI_Comm_size(...,eflag)) ! ERROR, final argument must be ierr
55649c86fc7SBarry Smith .ve
55749c86fc7SBarry Smith 
558db781477SPatrick Sanan .seealso: `SETERRMPI()`, `PetscCall()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`,
5593ba16761SJacob Faibussowitsch           `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
5603ba16761SJacob Faibussowitsch           `PetscError()`, `CHKMEMQ`
5613ba16761SJacob Faibussowitsch M*/
5623ba16761SJacob Faibussowitsch 
5633ba16761SJacob Faibussowitsch /*MC
5643ba16761SJacob Faibussowitsch   PetscCallMPIAbort - Like `PetscCallMPI()` but calls `MPI_Abort()` on error
5653ba16761SJacob Faibussowitsch 
5663ba16761SJacob Faibussowitsch   Synopsis:
5673ba16761SJacob Faibussowitsch   #include <petscerror.h>
5683ba16761SJacob Faibussowitsch   void PetscCallMPIAbort(MPI_Comm comm, MPI_Function(args))
5693ba16761SJacob Faibussowitsch 
5703ba16761SJacob Faibussowitsch   Not Collective
5713ba16761SJacob Faibussowitsch 
5723ba16761SJacob Faibussowitsch   Input Parameters:
5733ba16761SJacob Faibussowitsch + comm         - the MPI communicator to abort on
5743ba16761SJacob Faibussowitsch - MPI_Function - an MPI function that returns an MPI error code
5753ba16761SJacob Faibussowitsch 
576667f096bSBarry Smith   Level: beginner
577667f096bSBarry Smith 
5783ba16761SJacob Faibussowitsch   Notes:
5793ba16761SJacob Faibussowitsch   Usage is identical to `PetscCallMPI()`. See `PetscCallMPI()` for detailed discussion.
5803ba16761SJacob Faibussowitsch 
5813ba16761SJacob Faibussowitsch   This routine may be used in functions returning `void` or other non-`PetscErrorCode` types.
5823ba16761SJacob Faibussowitsch 
583989c49a2SBarry Smith   Fortran Note:
584989c49a2SBarry Smith   In Fortran this is called `PetscCallMPIA()` and is intended to be used in the main program while `PetscCallMPI()` is
585989c49a2SBarry Smith   used in Fortran subroutines.
586989c49a2SBarry Smith 
587989c49a2SBarry Smith   Developer Note:
588989c49a2SBarry Smith   This should have the same name in Fortran.
589989c49a2SBarry Smith 
5903ba16761SJacob Faibussowitsch .seealso: `PetscCallMPI()`, `PetscCallAbort()`, `SETERRABORT()`
59130de9b25SBarry Smith M*/
5923fcd9f07SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
59363a3b9bcSJacob Faibussowitsch void PetscCallMPI(PetscMPIInt);
5943ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm, PetscMPIInt);
595064a246eSJacob Faibussowitsch #else
5963ba16761SJacob Faibussowitsch   #define PetscCallMPI_Private(__PETSC_STACK_POP_FUNC__, __SETERR_FUNC__, __COMM__, ...) \
5979371c9d4SSatish Balay     do { \
5983ba16761SJacob Faibussowitsch       PetscMPIInt ierr_petsc_call_mpi_; \
599ef1023bdSBarry Smith       PetscStackUpdateLine; \
600792fecdfSBarry Smith       PetscStackPushExternal("MPI function"); \
601d71ae5a4SJacob Faibussowitsch       { \
6023ba16761SJacob Faibussowitsch         ierr_petsc_call_mpi_ = __VA_ARGS__; \
603d71ae5a4SJacob Faibussowitsch       } \
6043ba16761SJacob Faibussowitsch       __PETSC_STACK_POP_FUNC__; \
6053ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_mpi_ != MPI_SUCCESS)) { \
6063ba16761SJacob Faibussowitsch         char petsc_mpi_7_errorstring[2 * MPI_MAX_ERROR_STRING]; \
6073ba16761SJacob Faibussowitsch         PetscMPIErrorString(ierr_petsc_call_mpi_, (char *)petsc_mpi_7_errorstring); \
6083ba16761SJacob Faibussowitsch         __SETERR_FUNC__(__COMM__, PETSC_ERR_MPI, "MPI error %d %s", (int)ierr_petsc_call_mpi_, petsc_mpi_7_errorstring); \
6093fcd9f07SJacob Faibussowitsch       } \
6103fcd9f07SJacob Faibussowitsch     } while (0)
6113ba16761SJacob Faibussowitsch 
6123ba16761SJacob Faibussowitsch   #define PetscCallMPI(...)            PetscCallMPI_Private(PetscStackPop, SETERRQ, PETSC_COMM_SELF, __VA_ARGS__)
6133ba16761SJacob Faibussowitsch   #define PetscCallMPIAbort(comm, ...) PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRABORT, comm, __VA_ARGS__)
6149566063dSJacob Faibussowitsch #endif
615064a246eSJacob Faibussowitsch 
6167037b0edSPatrick Sanan /*MC
6179566063dSJacob Faibussowitsch   CHKERRMPI - Checks error code returned from MPI calls, if non-zero it calls the error
6189566063dSJacob Faibussowitsch   handler and then returns
6199566063dSJacob Faibussowitsch 
6209566063dSJacob Faibussowitsch   Synopsis:
6219566063dSJacob Faibussowitsch   #include <petscerror.h>
6229566063dSJacob Faibussowitsch   void CHKERRMPI(PetscErrorCode ierr)
6239566063dSJacob Faibussowitsch 
6249566063dSJacob Faibussowitsch   Not Collective
6259566063dSJacob Faibussowitsch 
6269566063dSJacob Faibussowitsch   Input Parameter:
6279566063dSJacob Faibussowitsch . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
6289566063dSJacob Faibussowitsch 
6299566063dSJacob Faibussowitsch   Level: deprecated
6309566063dSJacob Faibussowitsch 
631667f096bSBarry Smith   Note:
632667f096bSBarry Smith   Deprecated in favor of `PetscCallMPI()`. This routine behaves identically to it.
633667f096bSBarry Smith 
634db781477SPatrick Sanan .seealso: `PetscCallMPI()`
6359566063dSJacob Faibussowitsch M*/
6369566063dSJacob Faibussowitsch #define CHKERRMPI(...) PetscCallMPI(__VA_ARGS__)
6379566063dSJacob Faibussowitsch 
6389566063dSJacob Faibussowitsch /*MC
639989c49a2SBarry Smith   PetscCallAbort - Checks error code returned from PETSc function, if non-zero it aborts immediately by calling `MPI_Abort()`
6409566063dSJacob Faibussowitsch 
6419566063dSJacob Faibussowitsch   Synopsis:
6429566063dSJacob Faibussowitsch   #include <petscerror.h>
6439566063dSJacob Faibussowitsch   void PetscCallAbort(MPI_Comm comm, PetscErrorCode ierr)
6449566063dSJacob Faibussowitsch 
645c3339decSBarry Smith   Collective
6469566063dSJacob Faibussowitsch 
6479566063dSJacob Faibussowitsch   Input Parameters:
6489566063dSJacob Faibussowitsch + comm - the MPI communicator on which to abort
6499566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
6509566063dSJacob Faibussowitsch 
651667f096bSBarry Smith   Level: intermediate
652667f096bSBarry Smith 
6539566063dSJacob Faibussowitsch   Notes:
65487497f52SBarry Smith   This macro has identical type and usage semantics to `PetscCall()` with the important caveat
6559566063dSJacob Faibussowitsch   that this macro does not return. Instead, if ierr is nonzero it calls the PETSc error handler
65687497f52SBarry Smith   and then immediately calls `MPI_Abort()`. It can therefore be used anywhere.
6579566063dSJacob Faibussowitsch 
65887497f52SBarry Smith   As per `MPI_Abort()` semantics the communicator passed must be valid, although there is currently
65987497f52SBarry Smith   no attempt made at handling any potential errors from `MPI_Abort()`. Note that while
66087497f52SBarry Smith   `MPI_Abort()` is required to terminate only those processes which reside on comm, it is often
66187497f52SBarry Smith   the case that `MPI_Abort()` terminates *all* processes.
6629566063dSJacob Faibussowitsch 
6639566063dSJacob Faibussowitsch   Example Usage:
6649566063dSJacob Faibussowitsch .vb
6659566063dSJacob Faibussowitsch   PetscErrorCode boom(void) { return PETSC_ERR_MEM; }
6669566063dSJacob Faibussowitsch 
6679566063dSJacob Faibussowitsch   void foo(void)
6689566063dSJacob Faibussowitsch   {
6699566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type
6709566063dSJacob Faibussowitsch   }
6719566063dSJacob Faibussowitsch 
6729566063dSJacob Faibussowitsch   double bar(void)
6739566063dSJacob Faibussowitsch   {
6749566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type
6759566063dSJacob Faibussowitsch   }
6769566063dSJacob Faibussowitsch 
6779566063dSJacob Faibussowitsch   PetscCallAbort(MPI_COMM_NULL,boom()); // ERROR, communicator should be valid
6789566063dSJacob Faibussowitsch 
6799566063dSJacob Faibussowitsch   struct baz
6809566063dSJacob Faibussowitsch   {
6819566063dSJacob Faibussowitsch     baz()
6829566063dSJacob Faibussowitsch     {
6839566063dSJacob Faibussowitsch       PetscCallAbort(PETSC_COMM_SELF,boom()); // OK
6849566063dSJacob Faibussowitsch     }
6859566063dSJacob Faibussowitsch 
6869566063dSJacob Faibussowitsch     ~baz()
6879566063dSJacob Faibussowitsch     {
6889566063dSJacob Faibussowitsch       PetscCallAbort(PETSC_COMM_SELF,boom()); // OK (in fact the only way to handle PETSc errors)
6899566063dSJacob Faibussowitsch     }
6909566063dSJacob Faibussowitsch   };
6919566063dSJacob Faibussowitsch .ve
6929566063dSJacob Faibussowitsch 
693989c49a2SBarry Smith   Fortran Note:
694989c49a2SBarry Smith   Use `PetscCallA()`.
695989c49a2SBarry Smith 
696989c49a2SBarry Smith   Developer Note:
697989c49a2SBarry Smith   This should have the same name in Fortran as in C.
698989c49a2SBarry Smith 
699db781477SPatrick Sanan .seealso: `SETERRABORT()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`,
7005687f061SJacob Faibussowitsch           `SETERRQ()`, `CHKMEMQ`, `PetscCallMPI()`, `PetscCallCXXAbort()`
7019566063dSJacob Faibussowitsch M*/
7029566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
7039566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm, PetscErrorCode);
7049566063dSJacob Faibussowitsch void PetscCallContinue(PetscErrorCode);
7059566063dSJacob Faibussowitsch #else
7069371c9d4SSatish Balay   #define PetscCallAbort(comm, ...) \
7079371c9d4SSatish Balay     do { \
7087da6b3ccSLisandro Dalcin       PetscErrorCode ierr_petsc_call_abort_; \
7093ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
7107da6b3ccSLisandro Dalcin       ierr_petsc_call_abort_ = __VA_ARGS__; \
7113ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_abort_ != PETSC_SUCCESS)) { \
7123ba16761SJacob Faibussowitsch         ierr_petsc_call_abort_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_abort_, PETSC_ERROR_REPEAT, " "); \
7133ba16761SJacob Faibussowitsch         (void)MPI_Abort(comm, (PetscMPIInt)ierr_petsc_call_abort_); \
7149566063dSJacob Faibussowitsch       } \
7159566063dSJacob Faibussowitsch     } while (0)
7169371c9d4SSatish Balay   #define PetscCallContinue(...) \
7179371c9d4SSatish Balay     do { \
7187da6b3ccSLisandro Dalcin       PetscErrorCode ierr_petsc_call_continue_; \
7193ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
7207da6b3ccSLisandro Dalcin       ierr_petsc_call_continue_ = __VA_ARGS__; \
7213ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_continue_ != PETSC_SUCCESS)) { \
7223ba16761SJacob Faibussowitsch         ierr_petsc_call_continue_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_continue_, PETSC_ERROR_REPEAT, " "); \
7233ba16761SJacob Faibussowitsch         (void)ierr_petsc_call_continue_; \
7243ba16761SJacob Faibussowitsch       } \
7259566063dSJacob Faibussowitsch     } while (0)
7269566063dSJacob Faibussowitsch #endif
7279566063dSJacob Faibussowitsch 
7289566063dSJacob Faibussowitsch /*MC
7299566063dSJacob Faibussowitsch   CHKERRABORT - Checks error code returned from PETSc function. If non-zero it aborts immediately.
7309566063dSJacob Faibussowitsch 
7319566063dSJacob Faibussowitsch   Synopsis:
7329566063dSJacob Faibussowitsch   #include <petscerror.h>
7339566063dSJacob Faibussowitsch   void CHKERRABORT(MPI_Comm comm, PetscErrorCode ierr)
7349566063dSJacob Faibussowitsch 
7359566063dSJacob Faibussowitsch   Not Collective
7369566063dSJacob Faibussowitsch 
7379566063dSJacob Faibussowitsch   Input Parameters:
7389566063dSJacob Faibussowitsch + comm - the MPI communicator
7399566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
7409566063dSJacob Faibussowitsch 
7419566063dSJacob Faibussowitsch   Level: deprecated
7429566063dSJacob Faibussowitsch 
743667f096bSBarry Smith   Note:
744667f096bSBarry Smith   Deprecated in favor of `PetscCallAbort()`. This routine behaves identically to it.
745667f096bSBarry Smith 
746db781477SPatrick Sanan .seealso: `PetscCallAbort()`
7479566063dSJacob Faibussowitsch M*/
7489566063dSJacob Faibussowitsch #define CHKERRABORT(comm, ...) PetscCallAbort(comm, __VA_ARGS__)
7499566063dSJacob Faibussowitsch #define CHKERRCONTINUE(...)    PetscCallContinue(__VA_ARGS__)
7509566063dSJacob Faibussowitsch 
7519566063dSJacob Faibussowitsch /*MC
75287497f52SBarry Smith    CHKERRA - Fortran-only replacement for use of `CHKERRQ()` in the main program, which aborts immediately
753f388eb8bSPatrick Sanan 
754f388eb8bSPatrick Sanan    Synopsis:
755f388eb8bSPatrick Sanan    #include <petscsys.h>
756f388eb8bSPatrick Sanan    PetscErrorCode CHKERRA(PetscErrorCode ierr)
757f388eb8bSPatrick Sanan 
758f388eb8bSPatrick Sanan    Not Collective
759f388eb8bSPatrick Sanan 
7602fe279fdSBarry Smith    Input Parameter:
761f388eb8bSPatrick Sanan .  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
762f388eb8bSPatrick Sanan 
76387497f52SBarry Smith   Level: deprecated
764f388eb8bSPatrick Sanan 
76587497f52SBarry Smith    Note:
76687497f52SBarry Smith    This macro is rarely needed, normal usage is `PetscCallA()` in the main Fortran program.
767f388eb8bSPatrick Sanan 
768989c49a2SBarry Smith    Developer Note:
769989c49a2SBarry Smith    Why isn't this named `CHKERRABORT()` in Fortran?
770989c49a2SBarry Smith 
77187497f52SBarry Smith .seealso: `PetscCall()`, `PetscCallA()`, `PetscCallAbort()`, `CHKERRQ()`, `SETERRA()`, `SETERRQ()`, `SETERRABORT()`
772f388eb8bSPatrick Sanan M*/
773f388eb8bSPatrick Sanan 
7741e4f893aSSatish Balay PETSC_EXTERN PetscBool petscwaitonerrorflg;
7751e4f893aSSatish Balay PETSC_EXTERN PetscBool petscindebugger;
7767c66cc67SJunchao Zhang 
7777c66cc67SJunchao Zhang /*MC
778667f096bSBarry Smith    PETSCABORT - Call `MPI_Abort()` with an informative error code
7797c66cc67SJunchao Zhang 
7807c66cc67SJunchao Zhang    Synopsis:
7817c66cc67SJunchao Zhang    #include <petscsys.h>
7827c66cc67SJunchao Zhang    PETSCABORT(MPI_Comm comm, PetscErrorCode ierr)
7837c66cc67SJunchao Zhang 
7847cdbe19fSJose E. Roman    Collective; No Fortran Support
7857c66cc67SJunchao Zhang 
7867c66cc67SJunchao Zhang    Input Parameters:
7877c66cc67SJunchao Zhang +  comm - A communicator, so that the error can be collective
7887c66cc67SJunchao Zhang -  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
7897c66cc67SJunchao Zhang 
790bf4d2887SBarry Smith    Level: advanced
7917c66cc67SJunchao Zhang 
792bf4d2887SBarry Smith    Notes:
793989c49a2SBarry Smith    If the option `-start_in_debugger` was used then this calls `abort()` to stop the program in the debugger.
794bf4d2887SBarry Smith 
795989c49a2SBarry Smith    if `PetscCIEnabledPortableErrorOutput` is set, which means the code is running in the PETSc test harness (make test),
796989c49a2SBarry Smith    and `comm` is `MPI_COMM_WORLD` it strives to exit cleanly without calling `MPI_Abort()` and instead calling `MPI_Finalize()`.
797660278c0SBarry Smith 
798989c49a2SBarry Smith    This is currently only used when an error propagates up to the C `main()` program and is detected by a `PetscCall()`, `PetscCallMPI()`,
799989c49a2SBarry Smith    or is set in `main()` with `SETERRQ()`. Abort calls such as `SETERRABORT()`,
800989c49a2SBarry Smith    `PetscCheckAbort()`, `PetscCallMPIAbort()`, and `PetscCallAbort()` always call `MPI_Abort()` and do not have any special
801989c49a2SBarry Smith    handling for the test harness.
802989c49a2SBarry Smith 
803989c49a2SBarry Smith    Developer Note:
804989c49a2SBarry Smith    Should the other abort calls also pass through this call instead of calling `MPI_Abort()` directly?
805989c49a2SBarry Smith 
806989c49a2SBarry Smith .seealso: `PetscError()`, `PetscCall()`, `SETERRABORT()`, `PetscCheckAbort()`, `PetscCallMPIAbort()`, `PetscCall()`, `PetscCallMPI()`,
807989c49a2SBarry Smith           `PetscCallAbort()`, `MPI_Abort()`
8087c66cc67SJunchao Zhang M*/
80929f88feaSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
81029f88feaSJacob Faibussowitsch void PETSCABORT(MPI_Comm, PetscErrorCode);
81129f88feaSJacob Faibussowitsch #else
8129371c9d4SSatish Balay   #define PETSCABORT(comm, ...) \
8139371c9d4SSatish Balay     do { \
8143ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_abort_; \
8153ba16761SJacob Faibussowitsch       if (petscwaitonerrorflg) { ierr_petsc_abort_ = PetscSleep(1000); } \
8163ba16761SJacob Faibussowitsch       if (petscindebugger) { \
8173ba16761SJacob Faibussowitsch         abort(); \
8183ba16761SJacob Faibussowitsch       } else { \
8197da6b3ccSLisandro Dalcin         PetscMPIInt size_; \
8203ba16761SJacob Faibussowitsch         ierr_petsc_abort_ = __VA_ARGS__; \
8217da6b3ccSLisandro Dalcin         MPI_Comm_size(comm, &size_); \
8227da6b3ccSLisandro Dalcin         if (PetscCIEnabledPortableErrorOutput && size_ == PetscGlobalSize && ierr_petsc_abort_ != PETSC_ERR_SIG) { \
8239371c9d4SSatish Balay           MPI_Finalize(); \
8249371c9d4SSatish Balay           exit(0); \
825660278c0SBarry Smith         } else if (PetscCIEnabledPortableErrorOutput && PetscGlobalSize == 1) { \
826660278c0SBarry Smith           exit(0); \
827660278c0SBarry Smith         } else { \
828660278c0SBarry Smith           MPI_Abort(comm, (PetscMPIInt)ierr_petsc_abort_); \
829660278c0SBarry Smith         } \
8303fcd9f07SJacob Faibussowitsch       } \
8317c66cc67SJunchao Zhang     } while (0)
83229f88feaSJacob Faibussowitsch #endif
833986eef2eSBarry Smith 
8349566063dSJacob Faibussowitsch #ifdef PETSC_CLANGUAGE_CXX
835986eef2eSBarry Smith   /*MC
8369566063dSJacob Faibussowitsch   PetscCallThrow - Checks error code, if non-zero it calls the C++ error handler which throws
8379566063dSJacob Faibussowitsch   an exception
838986eef2eSBarry Smith 
839986eef2eSBarry Smith   Synopsis:
8409566063dSJacob Faibussowitsch   #include <petscerror.h>
8419566063dSJacob Faibussowitsch   void PetscCallThrow(PetscErrorCode ierr)
842986eef2eSBarry Smith 
843986eef2eSBarry Smith   Not Collective
844986eef2eSBarry Smith 
8459566063dSJacob Faibussowitsch   Input Parameter:
846986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
847986eef2eSBarry Smith 
848667f096bSBarry Smith   Level: beginner
849667f096bSBarry Smith 
850986eef2eSBarry Smith   Notes:
8519566063dSJacob Faibussowitsch   Requires PETSc to be configured with clanguage = c++. Throws a std::runtime_error() on error.
852986eef2eSBarry Smith 
85387497f52SBarry Smith   Once the error handler throws the exception you can use `PetscCallVoid()` which returns without
85487497f52SBarry Smith   an error code (bad idea since the error is ignored) or `PetscCallAbort()` to have `MPI_Abort()`
8559566063dSJacob Faibussowitsch   called immediately.
8569566063dSJacob Faibussowitsch 
857db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscCall()`, `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`,
858db781477SPatrick Sanan           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`
859986eef2eSBarry Smith M*/
8609371c9d4SSatish Balay   #define PetscCallThrow(...) \
8619371c9d4SSatish Balay     do { \
8623ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
8633ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_throw_ = __VA_ARGS__; \
8643ba16761SJacob 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); \
865ffc4695bSBarry Smith     } while (0)
86685614651SBarry Smith 
867cc26af49SMatthew Knepley   /*MC
868cc26af49SMatthew Knepley   CHKERRXX - Checks error code, if non-zero it calls the C++ error handler which throws an exception
869cc26af49SMatthew Knepley 
870cc26af49SMatthew Knepley   Synopsis:
8719566063dSJacob Faibussowitsch   #include <petscerror.h>
8723af045c5SBarry Smith   void CHKERRXX(PetscErrorCode ierr)
873cc26af49SMatthew Knepley 
874eca87e8dSBarry Smith   Not Collective
875cc26af49SMatthew Knepley 
8769566063dSJacob Faibussowitsch   Input Parameter:
8773af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
878cc26af49SMatthew Knepley 
8799566063dSJacob Faibussowitsch   Level: deprecated
880cc26af49SMatthew Knepley 
881667f096bSBarry Smith   Note:
882667f096bSBarry Smith   Deprecated in favor of `PetscCallThrow()`. This routine behaves identically to it.
883667f096bSBarry Smith 
884db781477SPatrick Sanan .seealso: `PetscCallThrow()`
885cc26af49SMatthew Knepley M*/
8869566063dSJacob Faibussowitsch   #define CHKERRXX(...) PetscCallThrow(__VA_ARGS__)
887fd705b32SMatthew Knepley #endif
888fd705b32SMatthew Knepley 
8893ba16761SJacob Faibussowitsch #define PetscCallCXX_Private(__SETERR_FUNC__, __COMM__, ...) \
8903ba16761SJacob Faibussowitsch   do { \
8913ba16761SJacob Faibussowitsch     PetscStackUpdateLine; \
8923ba16761SJacob Faibussowitsch     try { \
8933ba16761SJacob Faibussowitsch       __VA_ARGS__; \
8943ba16761SJacob Faibussowitsch     } catch (const std::exception &e) { \
8953ba16761SJacob Faibussowitsch       __SETERR_FUNC__(__COMM__, PETSC_ERR_LIB, "%s", e.what()); \
8963ba16761SJacob Faibussowitsch     } \
8973ba16761SJacob Faibussowitsch   } while (0)
8983ba16761SJacob Faibussowitsch 
8993f520e80SJunchao Zhang /*MC
9009566063dSJacob Faibussowitsch   PetscCallCXX - Checks C++ function calls and if they throw an exception, catch it and then
9019566063dSJacob Faibussowitsch   return a PETSc error code
9023f520e80SJunchao Zhang 
9033f520e80SJunchao Zhang   Synopsis:
9049566063dSJacob Faibussowitsch   #include <petscerror.h>
9055687f061SJacob Faibussowitsch   void PetscCallCXX(...) noexcept;
9063f520e80SJunchao Zhang 
9073f520e80SJunchao Zhang   Not Collective
9083f520e80SJunchao Zhang 
9099566063dSJacob Faibussowitsch   Input Parameter:
9105687f061SJacob Faibussowitsch . __VA_ARGS__ - An arbitrary expression
9115687f061SJacob Faibussowitsch 
9125687f061SJacob Faibussowitsch   Level: beginner
9133f520e80SJunchao Zhang 
9143f520e80SJunchao Zhang   Notes:
9155687f061SJacob Faibussowitsch   `PetscCallCXX(...)` is a macro replacement for
9169566063dSJacob Faibussowitsch .vb
9179566063dSJacob Faibussowitsch   try {
9185687f061SJacob Faibussowitsch     __VA_ARGS__;
9199566063dSJacob Faibussowitsch   } catch (const std::exception& e) {
9209566063dSJacob Faibussowitsch     return ConvertToPetscErrorCode(e);
9219566063dSJacob Faibussowitsch   }
9229566063dSJacob Faibussowitsch .ve
9239566063dSJacob Faibussowitsch   Due to the fact that it catches any (reasonable) exception, it is essentially noexcept.
9243f520e80SJunchao Zhang 
9255687f061SJacob Faibussowitsch   If you cannot return a `PetscErrorCode` use `PetscCallCXXAbort()` instead.
9265687f061SJacob Faibussowitsch 
9279566063dSJacob Faibussowitsch   Example Usage:
9289566063dSJacob Faibussowitsch .vb
9299566063dSJacob Faibussowitsch   void foo(void) { throw std::runtime_error("error"); }
9303f520e80SJunchao Zhang 
9319566063dSJacob Faibussowitsch   void bar()
9329566063dSJacob Faibussowitsch   {
933e8952933SJacob Faibussowitsch     PetscCallCXX(foo()); // ERROR bar() does not return PetscErrorCode
9349566063dSJacob Faibussowitsch   }
9359566063dSJacob Faibussowitsch 
9369566063dSJacob Faibussowitsch   PetscErrorCode baz()
9379566063dSJacob Faibussowitsch   {
9389566063dSJacob Faibussowitsch     PetscCallCXX(foo()); // OK
9399566063dSJacob Faibussowitsch 
9409566063dSJacob Faibussowitsch     PetscCallCXX(
9419566063dSJacob Faibussowitsch       bar();
942d5b43468SJose E. Roman       foo(); // OK multiple statements allowed
9439566063dSJacob Faibussowitsch     );
9449566063dSJacob Faibussowitsch   }
9459566063dSJacob Faibussowitsch 
9469566063dSJacob Faibussowitsch   struct bop
9479566063dSJacob Faibussowitsch   {
9489566063dSJacob Faibussowitsch     bop()
9499566063dSJacob Faibussowitsch     {
950e8952933SJacob Faibussowitsch       PetscCallCXX(foo()); // ERROR returns PetscErrorCode, cannot be used in constructors
9519566063dSJacob Faibussowitsch     }
9529566063dSJacob Faibussowitsch   };
9539566063dSJacob Faibussowitsch 
954e8952933SJacob Faibussowitsch   // ERROR contains do-while, cannot be used as function-try block
9559566063dSJacob Faibussowitsch   PetscErrorCode qux() PetscCallCXX(
9569566063dSJacob Faibussowitsch     bar();
9579566063dSJacob Faibussowitsch     baz();
9589566063dSJacob Faibussowitsch     foo();
9599566063dSJacob Faibussowitsch     return 0;
9609566063dSJacob Faibussowitsch   )
9619566063dSJacob Faibussowitsch .ve
9629566063dSJacob Faibussowitsch 
9635687f061SJacob Faibussowitsch .seealso: `PetscCallCXXAbort()`, `PetscCallThrow()`, `SETERRQ()`, `PetscCall()`,
9645687f061SJacob Faibussowitsch           `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
9655687f061SJacob Faibussowitsch           `PetscError()`, `CHKMEMQ`
9663f520e80SJunchao Zhang M*/
9673ba16761SJacob Faibussowitsch #define PetscCallCXX(...) PetscCallCXX_Private(SETERRQ, PETSC_COMM_SELF, __VA_ARGS__)
9685687f061SJacob Faibussowitsch 
9695687f061SJacob Faibussowitsch /*MC
9705687f061SJacob Faibussowitsch   PetscCallCXXAbort - Like `PetscCallCXX()` but calls `MPI_Abort()` instead of returning an
9715687f061SJacob Faibussowitsch   error-code
9725687f061SJacob Faibussowitsch 
9735687f061SJacob Faibussowitsch   Synopsis:
9745687f061SJacob Faibussowitsch   #include <petscerror.h>
9755687f061SJacob Faibussowitsch   void PetscCallCXXAbort(MPI_Comm comm, ...) noexcept;
9765687f061SJacob Faibussowitsch 
97720f4b53cSBarry Smith   Collective; No Fortran Support
9785687f061SJacob Faibussowitsch 
9795687f061SJacob Faibussowitsch   Input Parameters:
9805687f061SJacob Faibussowitsch + comm        - The MPI communicator to abort on
9815687f061SJacob Faibussowitsch - __VA_ARGS__ - An arbitrary expression
9825687f061SJacob Faibussowitsch 
9835687f061SJacob Faibussowitsch   Level: beginner
9845687f061SJacob Faibussowitsch 
9855687f061SJacob Faibussowitsch   Notes:
9865687f061SJacob Faibussowitsch   This macro may be used to check C++ expressions for exceptions in cases where you cannot
9875687f061SJacob Faibussowitsch   return an error code. This includes constructors, destructors, copy/move assignment functions
9885687f061SJacob Faibussowitsch   or constructors among others.
9895687f061SJacob Faibussowitsch 
9905687f061SJacob Faibussowitsch   If an exception is caught, the macro calls `SETERRABORT()` on `comm`. The exception must
9915687f061SJacob Faibussowitsch   derive from `std::exception` in order to be caught.
9925687f061SJacob Faibussowitsch 
9935687f061SJacob Faibussowitsch   If the routine _can_ return an error-code it is highly advised to use `PetscCallCXX()`
9945687f061SJacob Faibussowitsch   instead.
9955687f061SJacob Faibussowitsch 
9965687f061SJacob Faibussowitsch   See `PetscCallCXX()` for additional discussion.
9975687f061SJacob Faibussowitsch 
9985687f061SJacob Faibussowitsch   Example Usage:
9995687f061SJacob Faibussowitsch .vb
10005687f061SJacob Faibussowitsch   class Foo
10015687f061SJacob Faibussowitsch   {
10025687f061SJacob Faibussowitsch     std::vector<int> data_;
10035687f061SJacob Faibussowitsch 
10045687f061SJacob Faibussowitsch   public:
10055687f061SJacob Faibussowitsch     // normally std::vector::reserve() may raise an exception, but since we handle it with
10065687f061SJacob Faibussowitsch     // PetscCallCXXAbort() we may mark this routine as noexcept!
10075687f061SJacob Faibussowitsch     Foo() noexcept
10085687f061SJacob Faibussowitsch     {
10095687f061SJacob Faibussowitsch       PetscCallCXXAbort(PETSC_COMM_SELF, data_.reserve(10));
10105687f061SJacob Faibussowitsch     }
10115687f061SJacob Faibussowitsch   };
10125687f061SJacob Faibussowitsch 
10135687f061SJacob Faibussowitsch   std::vector<int> bar()
10145687f061SJacob Faibussowitsch   {
10155687f061SJacob Faibussowitsch     std::vector<int> v;
10165687f061SJacob Faibussowitsch 
10175687f061SJacob Faibussowitsch     PetscFunctionBegin;
10185687f061SJacob Faibussowitsch     // OK!
10195687f061SJacob Faibussowitsch     PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1));
10205687f061SJacob Faibussowitsch     PetscFunctionReturn(v);
10215687f061SJacob Faibussowitsch   }
10225687f061SJacob Faibussowitsch 
10235687f061SJacob Faibussowitsch   PetscErrorCode baz()
10245687f061SJacob Faibussowitsch   {
10255687f061SJacob Faibussowitsch     std::vector<int> v;
10265687f061SJacob Faibussowitsch 
10275687f061SJacob Faibussowitsch     PetscFunctionBegin;
10285687f061SJacob Faibussowitsch     // WRONG! baz() returns a PetscErrorCode, prefer PetscCallCXX() instead
10295687f061SJacob Faibussowitsch     PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1));
10303ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
10315687f061SJacob Faibussowitsch   }
10325687f061SJacob Faibussowitsch .ve
10335687f061SJacob Faibussowitsch 
10345687f061SJacob Faibussowitsch .seealso: `PetscCallCXX()`, `SETERRABORT()`, `PetscCallAbort()`
10355687f061SJacob Faibussowitsch M*/
10363ba16761SJacob Faibussowitsch #define PetscCallCXXAbort(comm, ...) PetscCallCXX_Private(SETERRABORT, comm, __VA_ARGS__)
10373f520e80SJunchao Zhang 
103830de9b25SBarry Smith /*MC
10399566063dSJacob Faibussowitsch   CHKERRCXX - Checks C++ function calls and if they throw an exception, catch it and then
10409566063dSJacob Faibussowitsch   return a PETSc error code
10419566063dSJacob Faibussowitsch 
10429566063dSJacob Faibussowitsch   Synopsis:
10439566063dSJacob Faibussowitsch   #include <petscerror.h>
10449566063dSJacob Faibussowitsch   void CHKERRCXX(func) noexcept;
10459566063dSJacob Faibussowitsch 
10469566063dSJacob Faibussowitsch   Not Collective
10479566063dSJacob Faibussowitsch 
10489566063dSJacob Faibussowitsch   Input Parameter:
10499566063dSJacob Faibussowitsch . func - C++ function calls
10509566063dSJacob Faibussowitsch 
10519566063dSJacob Faibussowitsch   Level: deprecated
10529566063dSJacob Faibussowitsch 
1053667f096bSBarry Smith   Note:
1054667f096bSBarry Smith   Deprecated in favor of `PetscCallCXX()`. This routine behaves identically to it.
1055667f096bSBarry Smith 
1056db781477SPatrick Sanan .seealso: `PetscCallCXX()`
10579566063dSJacob Faibussowitsch M*/
10589566063dSJacob Faibussowitsch #define CHKERRCXX(...) PetscCallCXX(__VA_ARGS__)
10599566063dSJacob Faibussowitsch 
10609566063dSJacob Faibussowitsch /*MC
106130de9b25SBarry Smith    CHKMEMQ - Checks the memory for corruption, calls error handler if any is detected
106230de9b25SBarry Smith 
106330de9b25SBarry Smith    Synopsis:
1064aaa7dc30SBarry Smith    #include <petscsys.h>
106591d3bdf4SKris Buschelman    CHKMEMQ;
106630de9b25SBarry Smith 
1067eca87e8dSBarry Smith    Not Collective
1068eca87e8dSBarry Smith 
106930de9b25SBarry Smith   Level: beginner
107030de9b25SBarry Smith 
107130de9b25SBarry Smith    Notes:
1072a17b96a8SKyle Gerard Felker     We highly recommend using Valgrind https://petsc.org/release/faq/#valgrind or for NVIDIA CUDA systems
10735ed36255SBarry Smith     https://docs.nvidia.com/cuda/cuda-memcheck/index.html for finding memory problems. The ``CHKMEMQ`` macro is useful on systems that
10745ed36255SBarry Smith     do not have valgrind, but is not as good as valgrind or cuda-memcheck.
10751957e957SBarry Smith 
1076667f096bSBarry Smith     Must run with the option `-malloc_debug` (`-malloc_test` in debug mode; or if `PetscMallocSetDebug()` called) to enable this option
107730de9b25SBarry Smith 
107830de9b25SBarry Smith     Once the error handler is called the calling function is then returned from with the given error code.
107930de9b25SBarry Smith 
108030de9b25SBarry Smith     By defaults prints location where memory that is corrupted was allocated.
108130de9b25SBarry Smith 
108287497f52SBarry Smith     Use `CHKMEMA` for functions that return void
1083f621e05eSBarry Smith 
1084db781477SPatrick Sanan .seealso: `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `SETERRQ()`, `PetscMallocValidate()`
108530de9b25SBarry Smith M*/
10866d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
1087064a246eSJacob Faibussowitsch   #define CHKMEMQ
1088064a246eSJacob Faibussowitsch   #define CHKMEMA
10896d210af2SJacob Faibussowitsch #else
10909371c9d4SSatish Balay   #define CHKMEMQ \
10919371c9d4SSatish Balay     do { \
10923ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_memq_ = PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__); \
10933ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_memq_ != PETSC_SUCCESS)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_memq_, PETSC_ERROR_REPEAT, " "); \
109486d09637SLisandro Dalcin     } while (0)
10956d210af2SJacob Faibussowitsch   #define CHKMEMA PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__)
1096064a246eSJacob Faibussowitsch #endif
10979566063dSJacob Faibussowitsch 
1098668f157eSBarry Smith /*E
1099668f157eSBarry Smith   PetscErrorType - passed to the PETSc error handling routines indicating if this is the first or a later call to the error handlers
1100668f157eSBarry Smith 
1101668f157eSBarry Smith   Level: advanced
1102668f157eSBarry Smith 
1103667f096bSBarry Smith   Note:
110487497f52SBarry Smith   `PETSC_ERROR_IN_CXX` indicates the error was detected in C++ and an exception should be generated
1105d736bfebSBarry Smith 
1106667f096bSBarry Smith   Developer Note:
1107667f096bSBarry Smith     This is currently used to decide when to print the detailed information about the run in `PetscTraceBackErrorHandler()`
1108668f157eSBarry Smith 
110987497f52SBarry Smith .seealso: `PetscError()`, `SETERRQ()`
1110668f157eSBarry Smith E*/
11119371c9d4SSatish Balay typedef enum {
11129371c9d4SSatish Balay   PETSC_ERROR_INITIAL = 0,
11139371c9d4SSatish Balay   PETSC_ERROR_REPEAT  = 1,
11149371c9d4SSatish Balay   PETSC_ERROR_IN_CXX  = 2
11159371c9d4SSatish Balay } PetscErrorType;
11164b209cf6SBarry Smith 
1117eb9e708aSLisandro Dalcin #if defined(__clang_analyzer__)
1118eb9e708aSLisandro Dalcin __attribute__((analyzer_noreturn))
1119eb9e708aSLisandro Dalcin #endif
11209371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode
11219371c9d4SSatish Balay PetscError(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, ...) PETSC_ATTRIBUTE_COLD PETSC_ATTRIBUTE_FORMAT(7, 8);
1122eb9e708aSLisandro Dalcin 
1123014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscErrorPrintfInitialize(void);
11243ba16761SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscErrorMessage(PetscErrorCode, const char *[], char **);
1125d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscTraceBackErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1126d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscIgnoreErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1127d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscEmacsClientErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1128d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscMPIAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1129d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1130d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAttachDebuggerErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1131d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscReturnErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1132efca3c55SSatish Balay PETSC_EXTERN PetscErrorCode PetscPushErrorHandler(PetscErrorCode (*handler)(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *), void *);
1133014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopErrorHandler(void);
11348d359177SBarry Smith PETSC_EXTERN PetscErrorCode PetscSignalHandlerDefault(int, void *);
1135014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPushSignalHandler(PetscErrorCode (*)(int, void *), void *);
1136014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopSignalHandler(void);
113728559dc8SJed Brown PETSC_EXTERN PetscErrorCode PetscCheckPointerSetIntensity(PetscInt);
1138c2a741eeSJunchao Zhang PETSC_EXTERN void           PetscSignalSegvCheckPointerOrMpi(void);
1139edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 13, 0, "PetscSignalSegvCheckPointerOrMpi()", ) static inline void PetscSignalSegvCheckPointer(void)
1140d71ae5a4SJacob Faibussowitsch {
11419371c9d4SSatish Balay   PetscSignalSegvCheckPointerOrMpi();
11429371c9d4SSatish Balay }
1143329f5518SBarry Smith 
1144639ff905SBarry Smith /*MC
1145639ff905SBarry Smith     PetscErrorPrintf - Prints error messages.
1146639ff905SBarry Smith 
1147639ff905SBarry Smith    Synopsis:
1148aaa7dc30SBarry Smith     #include <petscsys.h>
1149639ff905SBarry Smith      PetscErrorCode (*PetscErrorPrintf)(const char format[],...);
1150639ff905SBarry Smith 
11517cdbe19fSJose E. Roman     Not Collective; No Fortran Support
11527cdbe19fSJose E. Roman 
1153f899ff85SJose E. Roman     Input Parameter:
1154cd05f99aSJacob Faibussowitsch .   format - the usual `printf()` format string
1155639ff905SBarry Smith 
1156639ff905SBarry Smith    Options Database Keys:
11571957e957SBarry Smith +    -error_output_stdout - cause error messages to be printed to stdout instead of the (default) stderr
1158e1bc860dSBarry Smith -    -error_output_none - to turn off all printing of error messages (does not change the way the error is handled.)
1159639ff905SBarry Smith 
1160cf53795eSBarry Smith    Level: developer
1161cf53795eSBarry Smith 
116295452b02SPatrick Sanan    Notes:
116395452b02SPatrick Sanan     Use
1164667f096bSBarry Smith .vb
1165667f096bSBarry Smith      PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the
1166667f096bSBarry Smith                         error is handled.) and
1167667f096bSBarry Smith      PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on or you can use your own function
1168667f096bSBarry Smith .ve
1169639ff905SBarry Smith      Use
1170667f096bSBarry Smith .vb
117187497f52SBarry Smith      `PETSC_STDERR` = FILE* obtained from a file open etc. to have stderr printed to the file.
117287497f52SBarry Smith      `PETSC_STDOUT` = FILE* obtained from a file open etc. to have stdout printed to the file.
1173667f096bSBarry Smith .ve
1174639ff905SBarry Smith 
1175639ff905SBarry Smith        Use
117687497f52SBarry Smith       `PetscPushErrorHandler()` to provide your own error handler that determines what kind of messages to print
1177639ff905SBarry Smith 
1178db781477SPatrick Sanan .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscHelpPrintf()`, `PetscPrintf()`, `PetscPushErrorHandler()`, `PetscVFPrintf()`, `PetscHelpPrintf()`
1179639ff905SBarry Smith M*/
11803ca90d2dSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode (*PetscErrorPrintf)(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1181639ff905SBarry Smith 
1182cf0818bdSBarry Smith /*E
1183cf0818bdSBarry Smith      PetscFPTrap - types of floating point exceptions that may be trapped
1184cf0818bdSBarry Smith 
118587497f52SBarry Smith      Currently only `PETSC_FP_TRAP_OFF` and `PETSC_FP_TRAP_ON` are handled. All others are treated as `PETSC_FP_TRAP_ON`.
1186cf0818bdSBarry Smith 
1187cf0818bdSBarry Smith      Level: intermediate
1188cf0818bdSBarry Smith 
11897de69702SBarry Smith .seealso: `PetscSetFPTrap()`, `PetscFPTrapPush()`
1190cf0818bdSBarry Smith  E*/
11919371c9d4SSatish Balay typedef enum {
11929371c9d4SSatish Balay   PETSC_FP_TRAP_OFF      = 0,
11939371c9d4SSatish Balay   PETSC_FP_TRAP_INDIV    = 1,
11949371c9d4SSatish Balay   PETSC_FP_TRAP_FLTOPERR = 2,
11959371c9d4SSatish Balay   PETSC_FP_TRAP_FLTOVF   = 4,
11969371c9d4SSatish Balay   PETSC_FP_TRAP_FLTUND   = 8,
11979371c9d4SSatish Balay   PETSC_FP_TRAP_FLTDIV   = 16,
11989371c9d4SSatish Balay   PETSC_FP_TRAP_FLTINEX  = 32
11999371c9d4SSatish Balay } PetscFPTrap;
1200bd2b07b1SBarry 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)
1201014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSetFPTrap(PetscFPTrap);
1202014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPush(PetscFPTrap);
1203014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPop(void);
1204aba4c478SBarry Smith PETSC_EXTERN PetscErrorCode PetscDetermineInitialFPTrap(void);
120554a8ef01SBarry Smith 
12063a40ed3dSBarry Smith /*
12073a40ed3dSBarry Smith       Allows the code to build a stack frame as it runs
12083a40ed3dSBarry Smith */
12093a40ed3dSBarry Smith 
121099cd645aSJed Brown #define PETSCSTACKSIZE 64
12113a40ed3dSBarry Smith typedef struct {
12120e33f6ddSBarry Smith   const char *function[PETSCSTACKSIZE];
12130e33f6ddSBarry Smith   const char *file[PETSCSTACKSIZE];
1214184914b5SBarry Smith   int         line[PETSCSTACKSIZE];
1215362febeeSStefano Zampini   int         petscroutine[PETSCSTACKSIZE]; /* 0 external called from petsc, 1 petsc functions, 2 petsc user functions */
1216184914b5SBarry Smith   int         currentsize;
1217a2f94806SJed Brown   int         hotdepth;
12184be741a6SBarry Smith   PetscBool   check; /* option to check for correct Push/Pop semantics, true for default petscstack but not other stacks */
12193a40ed3dSBarry Smith } PetscStack;
1220dfb7d7afSStefano Zampini #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
122127104ee2SJacob Faibussowitsch PETSC_EXTERN PetscStack petscstack;
122227104ee2SJacob Faibussowitsch #endif
12233a40ed3dSBarry Smith 
12245d12eec7SSatish Balay #if defined(PETSC_SERIALIZE_FUNCTIONS)
12255d12eec7SSatish Balay   #include <petsc/private/petscfptimpl.h>
12265d12eec7SSatish Balay   /*
12275d12eec7SSatish Balay    Registers the current function into the global function pointer to function name table
12285d12eec7SSatish Balay 
12295d12eec7SSatish Balay    Have to fix this to handle errors but cannot return error since used in PETSC_VIEWER_DRAW_() etc
12305d12eec7SSatish Balay */
12319371c9d4SSatish Balay   #define PetscRegister__FUNCT__() \
12329371c9d4SSatish Balay     do { \
12335d12eec7SSatish Balay       static PetscBool __chked = PETSC_FALSE; \
12345d12eec7SSatish Balay       if (!__chked) { \
12359371c9d4SSatish Balay         void *ptr; \
12363ba16761SJacob Faibussowitsch         PetscCallAbort(PETSC_COMM_SELF, PetscDLSym(NULL, PETSC_FUNCTION_NAME, &ptr)); \
12375d12eec7SSatish Balay         __chked = PETSC_TRUE; \
12389371c9d4SSatish Balay       } \
12399371c9d4SSatish Balay     } while (0)
12405d12eec7SSatish Balay #else
12415d12eec7SSatish Balay   #define PetscRegister__FUNCT__()
12425d12eec7SSatish Balay #endif
12435d12eec7SSatish Balay 
1244eae3dc7dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) || defined(__clang_analyzer__)
12456d210af2SJacob Faibussowitsch   #define PetscStackPushNoCheck(funct, petsc_routine, hot)
124637154d10SBarry Smith   #define PetscStackUpdateLine
1247792fecdfSBarry Smith   #define PetscStackPushExternal(funct)
12486d210af2SJacob Faibussowitsch   #define PetscStackPopNoCheck
12496d210af2SJacob Faibussowitsch   #define PetscStackClearTop
12506d210af2SJacob Faibussowitsch   #define PetscFunctionBegin
12516d210af2SJacob Faibussowitsch   #define PetscFunctionBeginUser
12526d210af2SJacob Faibussowitsch   #define PetscFunctionBeginHot
12530a57284eSJacob Faibussowitsch   #define PetscFunctionReturn(...)  return __VA_ARGS__
12546d210af2SJacob Faibussowitsch   #define PetscFunctionReturnVoid() return
12556d210af2SJacob Faibussowitsch   #define PetscStackPop
12566d210af2SJacob Faibussowitsch   #define PetscStackPush(f)
1257dfb7d7afSStefano Zampini #elif defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1258660278c0SBarry Smith 
12599371c9d4SSatish Balay   #define PetscStackPush_Private(stack__, file__, func__, line__, petsc_routine__, hot__) \
12609371c9d4SSatish Balay     do { \
12615a96b57dSJacob Faibussowitsch       if (stack__.currentsize < PETSCSTACKSIZE) { \
12625a96b57dSJacob Faibussowitsch         stack__.function[stack__.currentsize] = func__; \
1263ef1023bdSBarry Smith         if (petsc_routine__) { \
1264ef1023bdSBarry Smith           stack__.file[stack__.currentsize] = file__; \
12655a96b57dSJacob Faibussowitsch           stack__.line[stack__.currentsize] = line__; \
1266ef1023bdSBarry Smith         } else { \
1267648bc8c4SBarry Smith           stack__.file[stack__.currentsize] = PETSC_NULLPTR; \
1268ef1023bdSBarry Smith           stack__.line[stack__.currentsize] = 0; \
1269ef1023bdSBarry Smith         } \
12705a96b57dSJacob Faibussowitsch         stack__.petscroutine[stack__.currentsize] = petsc_routine__; \
12715a96b57dSJacob Faibussowitsch       } \
12725a96b57dSJacob Faibussowitsch       ++stack__.currentsize; \
12735a96b57dSJacob Faibussowitsch       stack__.hotdepth += (hot__ || stack__.hotdepth); \
12745a96b57dSJacob Faibussowitsch     } while (0)
12755a96b57dSJacob Faibussowitsch 
12764be741a6SBarry Smith   /* uses PetscCheckAbort() because may be used in a function that does not return an error code */
12779371c9d4SSatish Balay   #define PetscStackPop_Private(stack__, func__) \
12789371c9d4SSatish Balay     do { \
12794be741a6SBarry 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__); \
12805a96b57dSJacob Faibussowitsch       if (--stack__.currentsize < PETSCSTACKSIZE) { \
12819371c9d4SSatish 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", \
12829371c9d4SSatish Balay                         stack__.function[stack__.currentsize], stack__.file[stack__.currentsize], stack__.line[stack__.currentsize], func__, __FILE__, __LINE__); \
12835a96b57dSJacob Faibussowitsch         stack__.function[stack__.currentsize]     = PETSC_NULLPTR; \
12845a96b57dSJacob Faibussowitsch         stack__.file[stack__.currentsize]         = PETSC_NULLPTR; \
12855a96b57dSJacob Faibussowitsch         stack__.line[stack__.currentsize]         = 0; \
12865a96b57dSJacob Faibussowitsch         stack__.petscroutine[stack__.currentsize] = 0; \
12875a96b57dSJacob Faibussowitsch       } \
12885a96b57dSJacob Faibussowitsch       stack__.hotdepth = PetscMax(stack__.hotdepth - 1, 0); \
12895a96b57dSJacob Faibussowitsch     } while (0)
12905a96b57dSJacob Faibussowitsch 
1291586f9135SBarry Smith   /*MC
1292586f9135SBarry Smith    PetscStackPushNoCheck - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is
1293586f9135SBarry Smith    currently in the source code.
1294586f9135SBarry Smith 
1295586f9135SBarry Smith    Synopsis:
1296586f9135SBarry Smith    #include <petscsys.h>
1297586f9135SBarry Smith    void PetscStackPushNoCheck(char *funct,int petsc_routine,PetscBool hot);
1298586f9135SBarry Smith 
12997cdbe19fSJose E. Roman    Not Collective
13007cdbe19fSJose E. Roman 
1301586f9135SBarry Smith    Input Parameters:
1302586f9135SBarry Smith +  funct - the function name
1303586f9135SBarry Smith .  petsc_routine - 2 user function, 1 PETSc function, 0 some other function
1304586f9135SBarry Smith -  hot - indicates that the function may be called often so expensive error checking should be turned off inside the function
1305586f9135SBarry Smith 
1306586f9135SBarry Smith    Level: developer
1307586f9135SBarry Smith 
1308586f9135SBarry Smith    Notes:
1309586f9135SBarry 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
131087497f52SBarry 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
1311586f9135SBarry Smith    help debug the problem.
1312586f9135SBarry Smith 
1313ef1023bdSBarry Smith    This version does not check the memory corruption (an expensive operation), use `PetscStackPush()` to check the memory.
1314ef1023bdSBarry Smith 
1315792fecdfSBarry Smith    Use `PetscStackPushExternal()` for a function call that is about to be made to a non-PETSc or user function (such as BLAS etc).
1316ef1023bdSBarry Smith 
1317586f9135SBarry Smith    The default stack is a global variable called `petscstack`.
1318586f9135SBarry Smith 
1319586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1320ef1023bdSBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPush()`, `PetscStackPop`,
1321792fecdfSBarry Smith           `PetscStackPushExternal()`
1322586f9135SBarry Smith M*/
13239371c9d4SSatish Balay   #define PetscStackPushNoCheck(funct, petsc_routine, hot) \
13249371c9d4SSatish Balay     do { \
1325e04113cfSBarry Smith       PetscStackSAWsTakeAccess(); \
13265a96b57dSJacob Faibussowitsch       PetscStackPush_Private(petscstack, __FILE__, funct, __LINE__, petsc_routine, hot); \
1327e04113cfSBarry Smith       PetscStackSAWsGrantAccess(); \
1328441dd030SJed Brown     } while (0)
1329441dd030SJed Brown 
1330586f9135SBarry Smith   /*MC
133187497f52SBarry Smith    PetscStackUpdateLine - in a function that has a `PetscFunctionBegin` or `PetscFunctionBeginUser` updates the stack line number to the
133237154d10SBarry Smith    current line number.
133337154d10SBarry Smith 
133437154d10SBarry Smith    Synopsis:
133537154d10SBarry Smith    #include <petscsys.h>
133637154d10SBarry Smith    void PetscStackUpdateLine
133737154d10SBarry Smith 
13387cdbe19fSJose E. Roman    Not Collective
13397cdbe19fSJose E. Roman 
134037154d10SBarry Smith    Level: developer
134137154d10SBarry Smith 
134237154d10SBarry Smith    Notes:
134387497f52SBarry Smith    Using `PetscCall()` and friends automatically handles this process
134487497f52SBarry Smith 
134537154d10SBarry 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
134637154d10SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
134737154d10SBarry Smith    help debug the problem.
134837154d10SBarry Smith 
134937154d10SBarry Smith    The default stack is a global variable called petscstack.
135037154d10SBarry Smith 
135137154d10SBarry Smith    This is used by `PetscCall()` and is otherwise not like to be needed
135237154d10SBarry Smith 
135337154d10SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`, `PetscCall()`
135437154d10SBarry Smith M*/
135537154d10SBarry Smith   #define PetscStackUpdateLine \
13563ba16761SJacob Faibussowitsch     do { \
13573ba16761SJacob Faibussowitsch       if (petscstack.currentsize > 0 && petscstack.function[petscstack.currentsize - 1] == PETSC_FUNCTION_NAME) { petscstack.line[petscstack.currentsize - 1] = __LINE__; } \
13583ba16761SJacob Faibussowitsch     } while (0)
135937154d10SBarry Smith 
136037154d10SBarry Smith   /*MC
1361792fecdfSBarry Smith    PetscStackPushExternal - Pushes a new function name onto the PETSc default stack that tracks where the running program is
1362ef1023bdSBarry Smith    currently in the source code. Does not include the filename or line number since this is called by the calling routine
1363ef1023bdSBarry Smith    for non-PETSc or user functions.
1364ef1023bdSBarry Smith 
1365ef1023bdSBarry Smith    Synopsis:
1366ef1023bdSBarry Smith    #include <petscsys.h>
1367792fecdfSBarry Smith    void PetscStackPushExternal(char *funct);
1368ef1023bdSBarry Smith 
13697cdbe19fSJose E. Roman    Not Collective
13707cdbe19fSJose E. Roman 
13712fe279fdSBarry Smith    Input Parameter:
1372ef1023bdSBarry Smith .  funct - the function name
1373ef1023bdSBarry Smith 
1374ef1023bdSBarry Smith    Level: developer
1375ef1023bdSBarry Smith 
1376ef1023bdSBarry Smith    Notes:
137787497f52SBarry Smith    Using `PetscCallExternal()` and friends automatically handles this process
137887497f52SBarry Smith 
1379ef1023bdSBarry 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
1380ef1023bdSBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1381ef1023bdSBarry Smith    help debug the problem.
1382ef1023bdSBarry Smith 
1383ef1023bdSBarry Smith    The default stack is a global variable called `petscstack`.
1384ef1023bdSBarry Smith 
1385ef1023bdSBarry Smith    This is to be used when calling an external package function such as a BLAS function.
1386ef1023bdSBarry Smith 
1387ef1023bdSBarry Smith    This also updates the stack line number for the current stack function.
1388ef1023bdSBarry Smith 
1389ef1023bdSBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1390ef1023bdSBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1391ef1023bdSBarry Smith M*/
13929371c9d4SSatish Balay   #define PetscStackPushExternal(funct) \
13939371c9d4SSatish Balay     do { \
13949371c9d4SSatish Balay       PetscStackUpdateLine; \
13959371c9d4SSatish Balay       PetscStackPushNoCheck(funct, 0, PETSC_TRUE); \
1396a8f51744SPierre Jolivet     } while (0)
1397ef1023bdSBarry Smith 
1398ef1023bdSBarry Smith   /*MC
1399586f9135SBarry Smith    PetscStackPopNoCheck - Pops a function name from the PETSc default stack that tracks where the running program is
1400586f9135SBarry Smith    currently in the source code.
1401586f9135SBarry Smith 
1402586f9135SBarry Smith    Synopsis:
1403586f9135SBarry Smith    #include <petscsys.h>
1404586f9135SBarry Smith    void PetscStackPopNoCheck(char *funct);
1405586f9135SBarry Smith 
14067cdbe19fSJose E. Roman    Not Collective
14077cdbe19fSJose E. Roman 
1408586f9135SBarry Smith    Input Parameter:
1409586f9135SBarry Smith .   funct - the function name
1410586f9135SBarry Smith 
1411586f9135SBarry Smith    Level: developer
1412586f9135SBarry Smith 
1413586f9135SBarry Smith    Notes:
141487497f52SBarry Smith    Using `PetscCall()`, `PetscCallExternal()`, `PetscCallBack()` and friends negates the need to call this
141587497f52SBarry Smith 
1416586f9135SBarry 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
1417586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1418586f9135SBarry Smith    help debug the problem.
1419586f9135SBarry Smith 
1420586f9135SBarry Smith    The default stack is a global variable called petscstack.
1421586f9135SBarry Smith 
1422586f9135SBarry Smith    Developer Note:
1423586f9135SBarry Smith    `PetscStackPopNoCheck()` takes a function argument while  `PetscStackPop` does not, this difference is likely just historical.
1424586f9135SBarry Smith 
1425586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1426586f9135SBarry Smith M*/
14279371c9d4SSatish Balay   #define PetscStackPopNoCheck(funct) \
14289371c9d4SSatish Balay     do { \
1429362febeeSStefano Zampini       PetscStackSAWsTakeAccess(); \
14305a96b57dSJacob Faibussowitsch       PetscStackPop_Private(petscstack, funct); \
1431362febeeSStefano Zampini       PetscStackSAWsGrantAccess(); \
1432362febeeSStefano Zampini     } while (0)
1433362febeeSStefano Zampini 
14349371c9d4SSatish Balay   #define PetscStackClearTop \
14359371c9d4SSatish Balay     do { \
1436e04113cfSBarry Smith       PetscStackSAWsTakeAccess(); \
14379371c9d4SSatish Balay       if (petscstack.currentsize > 0 && --petscstack.currentsize < PETSCSTACKSIZE) { \
143827104ee2SJacob Faibussowitsch         petscstack.function[petscstack.currentsize]     = PETSC_NULLPTR; \
143927104ee2SJacob Faibussowitsch         petscstack.file[petscstack.currentsize]         = PETSC_NULLPTR; \
144027104ee2SJacob Faibussowitsch         petscstack.line[petscstack.currentsize]         = 0; \
144127104ee2SJacob Faibussowitsch         petscstack.petscroutine[petscstack.currentsize] = 0; \
1442441dd030SJed Brown       } \
144327104ee2SJacob Faibussowitsch       petscstack.hotdepth = PetscMax(petscstack.hotdepth - 1, 0); \
1444e04113cfSBarry Smith       PetscStackSAWsGrantAccess(); \
1445441dd030SJed Brown     } while (0)
1446441dd030SJed Brown 
144730de9b25SBarry Smith   /*MC
14481957e957SBarry Smith    PetscFunctionBegin - First executable line of each PETSc function,  used for error handling. Final
144987497f52SBarry Smith       line of PETSc functions should be `PetscFunctionReturn`(0);
145030de9b25SBarry Smith 
145130de9b25SBarry Smith    Synopsis:
1452aaa7dc30SBarry Smith    #include <petscsys.h>
145330de9b25SBarry Smith    void PetscFunctionBegin;
145430de9b25SBarry Smith 
145520f4b53cSBarry Smith    Not Collective; No Fortran Support
1456eca87e8dSBarry Smith 
145730de9b25SBarry Smith    Usage:
145830de9b25SBarry Smith .vb
145930de9b25SBarry Smith      int something;
146030de9b25SBarry Smith 
146130de9b25SBarry Smith      PetscFunctionBegin;
146230de9b25SBarry Smith .ve
146330de9b25SBarry Smith 
146430de9b25SBarry Smith    Level: developer
146530de9b25SBarry Smith 
146620f4b53cSBarry Smith    Note:
146720f4b53cSBarry Smith      Use `PetscFunctionBeginUser` for application codes.
146820f4b53cSBarry Smith 
1469586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`
147030de9b25SBarry Smith 
147130de9b25SBarry Smith M*/
14729371c9d4SSatish Balay   #define PetscFunctionBegin \
14739371c9d4SSatish Balay     do { \
1474362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_FALSE); \
1475a2f94806SJed Brown       PetscRegister__FUNCT__(); \
1476a2f94806SJed Brown     } while (0)
1477a2f94806SJed Brown 
1478a2f94806SJed Brown   /*MC
147987497f52SBarry Smith    PetscFunctionBeginHot - Substitute for `PetscFunctionBegin` to be used in functions that are called in
1480a2f94806SJed Brown    performance-critical circumstances.  Use of this function allows for lighter profiling by default.
1481a2f94806SJed Brown 
1482a2f94806SJed Brown    Synopsis:
1483aaa7dc30SBarry Smith    #include <petscsys.h>
1484a2f94806SJed Brown    void PetscFunctionBeginHot;
1485a2f94806SJed Brown 
148620f4b53cSBarry Smith    Not Collective; No Fortran Support
1487a2f94806SJed Brown 
1488a2f94806SJed Brown    Usage:
1489a2f94806SJed Brown .vb
1490a2f94806SJed Brown      int something;
1491a2f94806SJed Brown 
1492a2f94806SJed Brown      PetscFunctionBeginHot;
1493a2f94806SJed Brown .ve
1494a2f94806SJed Brown 
1495a2f94806SJed Brown    Level: developer
1496a2f94806SJed Brown 
1497586f9135SBarry Smith .seealso: `PetscFunctionBegin`, `PetscFunctionReturn()`, `PetscStackPushNoCheck()`
1498a2f94806SJed Brown 
1499a2f94806SJed Brown M*/
15009371c9d4SSatish Balay   #define PetscFunctionBeginHot \
15019371c9d4SSatish Balay     do { \
1502362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_TRUE); \
15032d53ad75SBarry Smith       PetscRegister__FUNCT__(); \
150453c77d0aSJed Brown     } while (0)
150553c77d0aSJed Brown 
1506a8d2bbe5SBarry Smith   /*MC
1507530d85c1SBarry Smith    PetscFunctionBeginUser - First executable line of user provided routines
1508a8d2bbe5SBarry Smith 
1509a8d2bbe5SBarry Smith    Synopsis:
1510aaa7dc30SBarry Smith    #include <petscsys.h>
1511a8d2bbe5SBarry Smith    void PetscFunctionBeginUser;
1512a8d2bbe5SBarry Smith 
151320f4b53cSBarry Smith    Not Collective; No Fortran Support
1514a8d2bbe5SBarry Smith 
1515a8d2bbe5SBarry Smith    Usage:
1516a8d2bbe5SBarry Smith .vb
1517a8d2bbe5SBarry Smith      int something;
1518a8d2bbe5SBarry Smith 
1519ac285190SBarry Smith      PetscFunctionBeginUser;
1520a8d2bbe5SBarry Smith .ve
1521a8d2bbe5SBarry Smith 
152220f4b53cSBarry Smith    Level: intermediate
152320f4b53cSBarry Smith 
1524a8d2bbe5SBarry Smith    Notes:
1525530d85c1SBarry Smith       Functions that incorporate this must call `PetscFunctionReturn()` instead of return except for main().
1526530d85c1SBarry Smith 
1527530d85c1SBarry Smith       May be used before `PetscInitialize()`
15281957e957SBarry Smith 
1529530d85c1SBarry Smith       This is identical to `PetscFunctionBegin` except it labels the routine as a user
1530ac285190SBarry Smith       routine instead of as a PETSc library routine.
1531ac285190SBarry Smith 
1532586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, `PetscFunctionBeginHot`, `PetscStackPushNoCheck()`
1533a8d2bbe5SBarry Smith 
1534a8d2bbe5SBarry Smith M*/
15359371c9d4SSatish Balay   #define PetscFunctionBeginUser \
15369371c9d4SSatish Balay     do { \
1537362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 2, PETSC_FALSE); \
1538a8d2bbe5SBarry Smith       PetscRegister__FUNCT__(); \
1539a8d2bbe5SBarry Smith     } while (0)
1540a8d2bbe5SBarry Smith 
1541586f9135SBarry Smith   /*MC
1542586f9135SBarry Smith    PetscStackPush - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is
1543586f9135SBarry Smith    currently in the source code and verifies the memory is not corrupted.
1544586f9135SBarry Smith 
1545586f9135SBarry Smith    Synopsis:
1546586f9135SBarry Smith    #include <petscsys.h>
1547586f9135SBarry Smith    void PetscStackPush(char *funct)
1548586f9135SBarry Smith 
15497cdbe19fSJose E. Roman    Not Collective
15507cdbe19fSJose E. Roman 
1551586f9135SBarry Smith    Input Parameter:
1552586f9135SBarry Smith .  funct - the function name
1553586f9135SBarry Smith 
1554586f9135SBarry Smith    Level: developer
1555586f9135SBarry Smith 
1556586f9135SBarry Smith    Notes:
1557586f9135SBarry 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
1558586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1559586f9135SBarry Smith    help debug the problem.
1560586f9135SBarry Smith 
1561586f9135SBarry Smith    The default stack is a global variable called petscstack.
1562586f9135SBarry Smith 
1563586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1564586f9135SBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1565586f9135SBarry Smith M*/
15669371c9d4SSatish Balay   #define PetscStackPush(n) \
15679371c9d4SSatish Balay     do { \
1568362febeeSStefano Zampini       PetscStackPushNoCheck(n, 0, PETSC_FALSE); \
156915681b3cSBarry Smith       CHKMEMQ; \
157015681b3cSBarry Smith     } while (0)
15713a40ed3dSBarry Smith 
1572586f9135SBarry Smith   /*MC
1573586f9135SBarry Smith    PetscStackPop - Pops a function name from the PETSc default stack that tracks where the running program is
1574586f9135SBarry Smith    currently in the source code and verifies the memory is not corrupted.
1575586f9135SBarry Smith 
1576586f9135SBarry Smith    Synopsis:
1577586f9135SBarry Smith    #include <petscsys.h>
1578586f9135SBarry Smith    void PetscStackPop
1579586f9135SBarry Smith 
15807cdbe19fSJose E. Roman    Not Collective
15817cdbe19fSJose E. Roman 
1582586f9135SBarry Smith    Level: developer
1583586f9135SBarry Smith 
1584586f9135SBarry Smith    Notes:
1585586f9135SBarry 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
1586586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1587586f9135SBarry Smith    help debug the problem.
1588586f9135SBarry Smith 
1589586f9135SBarry Smith    The default stack is a global variable called petscstack.
1590586f9135SBarry Smith 
1591586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPopNoCheck()`, `PetscStackPush()`
1592586f9135SBarry Smith M*/
15939371c9d4SSatish Balay   #define PetscStackPop \
15949371c9d4SSatish Balay     do { \
1595441dd030SJed Brown       CHKMEMQ; \
1596362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
159715681b3cSBarry Smith     } while (0)
1598d64ed03dSBarry Smith 
159930de9b25SBarry Smith   /*MC
16000a57284eSJacob Faibussowitsch    PetscFunctionReturn - Last executable line of each PETSc function used for error
16010a57284eSJacob Faibussowitsch    handling. Replaces `return()`.
160230de9b25SBarry Smith 
160330de9b25SBarry Smith    Synopsis:
16040a57284eSJacob Faibussowitsch    #include <petscerror.h>
16050a57284eSJacob Faibussowitsch    void PetscFunctionReturn(...)
160630de9b25SBarry Smith 
1607667f096bSBarry Smith    Not Collective; No Fortran Support
1608eca87e8dSBarry Smith 
16095687f061SJacob Faibussowitsch    Level: beginner
16105687f061SJacob Faibussowitsch 
16110a57284eSJacob Faibussowitsch    Notes:
16120a57284eSJacob Faibussowitsch    This routine is a macro, so while it does not "return" anything itself, it does return from
16130a57284eSJacob Faibussowitsch    the function in the literal sense.
16140a57284eSJacob Faibussowitsch 
16150a57284eSJacob Faibussowitsch    Usually the return value is the integer literal `0` (for example in any function returning
16160a57284eSJacob Faibussowitsch    `PetscErrorCode`), however it is possible to return any arbitrary type. The arguments of
16170a57284eSJacob Faibussowitsch    this macro are placed before the `return` statement as-is.
16180a57284eSJacob Faibussowitsch 
16190a57284eSJacob Faibussowitsch    Any routine which returns via `PetscFunctionReturn()` must begin with a corresponding
16200a57284eSJacob Faibussowitsch    `PetscFunctionBegin`.
16210a57284eSJacob Faibussowitsch 
16220a57284eSJacob Faibussowitsch    For routines which return `void` use `PetscFunctionReturnVoid()` instead.
16230a57284eSJacob Faibussowitsch 
16240a57284eSJacob Faibussowitsch    Example Usage:
162530de9b25SBarry Smith .vb
16260a57284eSJacob Faibussowitsch    PetscErrorCode foo(int *x)
16270a57284eSJacob Faibussowitsch    {
16280a57284eSJacob Faibussowitsch      PetscFunctionBegin; // don't forget the begin!
16290a57284eSJacob Faibussowitsch      *x = 10;
16303ba16761SJacob Faibussowitsch      PetscFunctionReturn(PETSC_SUCCESS);
163130de9b25SBarry Smith    }
163230de9b25SBarry Smith .ve
163330de9b25SBarry Smith 
16340a57284eSJacob Faibussowitsch    May return any arbitrary type\:
16350a57284eSJacob Faibussowitsch .vb
16360a57284eSJacob Faibussowitsch   struct Foo
16370a57284eSJacob Faibussowitsch   {
16380a57284eSJacob Faibussowitsch     int x;
16390a57284eSJacob Faibussowitsch   };
16400a57284eSJacob Faibussowitsch 
16410a57284eSJacob Faibussowitsch   struct Foo make_foo(int value)
16420a57284eSJacob Faibussowitsch   {
16430a57284eSJacob Faibussowitsch     struct Foo f;
16440a57284eSJacob Faibussowitsch 
16450a57284eSJacob Faibussowitsch     PetscFunctionBegin;
16460a57284eSJacob Faibussowitsch     f.x = value;
16470a57284eSJacob Faibussowitsch     PetscFunctionReturn(f);
16480a57284eSJacob Faibussowitsch   }
16490a57284eSJacob Faibussowitsch .ve
16500a57284eSJacob Faibussowitsch 
16510a57284eSJacob Faibussowitsch .seealso: `PetscFunctionBegin`, `PetscFunctionBeginUser`, `PetscFunctionReturnVoid()`,
16520a57284eSJacob Faibussowitsch           `PetscStackPopNoCheck()`
165330de9b25SBarry Smith M*/
16545687f061SJacob Faibussowitsch   #define PetscFunctionReturn(...) \
16559371c9d4SSatish Balay     do { \
1656362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
16575687f061SJacob Faibussowitsch       return __VA_ARGS__; \
165827104ee2SJacob Faibussowitsch     } while (0)
1659d64ed03dSBarry Smith 
16600a57284eSJacob Faibussowitsch   /*MC
16610a57284eSJacob Faibussowitsch   PetscFunctionReturnVoid - Like `PetscFunctionReturn()` but returns `void`
16620a57284eSJacob Faibussowitsch 
16630a57284eSJacob Faibussowitsch   Synopsis:
16640a57284eSJacob Faibussowitsch   #include <petscerror.h>
16650a57284eSJacob Faibussowitsch   void PetscFunctionReturnVoid()
16660a57284eSJacob Faibussowitsch 
16670a57284eSJacob Faibussowitsch   Not Collective
16680a57284eSJacob Faibussowitsch 
16695687f061SJacob Faibussowitsch   Level: beginner
16705687f061SJacob Faibussowitsch 
16715687f061SJacob Faibussowitsch   Note:
16720a57284eSJacob Faibussowitsch   Behaves identically to `PetscFunctionReturn()` except that it returns `void`. That is, this
16735687f061SJacob Faibussowitsch   macro culminates with `return`.
16740a57284eSJacob Faibussowitsch 
16750a57284eSJacob Faibussowitsch   Example Usage:
16760a57284eSJacob Faibussowitsch .vb
16770a57284eSJacob Faibussowitsch   void foo()
16780a57284eSJacob Faibussowitsch   {
16790a57284eSJacob Faibussowitsch     PetscFunctionBegin; // must start with PetscFunctionBegin!
16800a57284eSJacob Faibussowitsch     bar();
16810a57284eSJacob Faibussowitsch     baz();
16820a57284eSJacob Faibussowitsch     PetscFunctionReturnVoid();
16830a57284eSJacob Faibussowitsch   }
16840a57284eSJacob Faibussowitsch .ve
16850a57284eSJacob Faibussowitsch 
16860a57284eSJacob Faibussowitsch .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, PetscFunctionBeginUser`
16870a57284eSJacob Faibussowitsch M*/
16889371c9d4SSatish Balay   #define PetscFunctionReturnVoid() \
16899371c9d4SSatish Balay     do { \
1690362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
169127104ee2SJacob Faibussowitsch       return; \
169227104ee2SJacob Faibussowitsch     } while (0)
169327104ee2SJacob Faibussowitsch #else /* PETSC_USE_DEBUG */
169427104ee2SJacob Faibussowitsch   #define PetscStackPushNoCheck(funct, petsc_routine, hot)
169537154d10SBarry Smith   #define PetscStackUpdateLine
1696792fecdfSBarry Smith   #define PetscStackPushExternal(funct)
16973ba16761SJacob Faibussowitsch   #define PetscStackPopNoCheck(...)
169827104ee2SJacob Faibussowitsch   #define PetscStackClearTop
16993a40ed3dSBarry Smith   #define PetscFunctionBegin
17000bdf7c52SPeter Brune   #define PetscFunctionBeginUser
1701a2f94806SJed Brown   #define PetscFunctionBeginHot
17020a57284eSJacob Faibussowitsch   #define PetscFunctionReturn(...)  return __VA_ARGS__
17035665465eSBarry Smith   #define PetscFunctionReturnVoid() return
1704812af9f3SBarry Smith   #define PetscStackPop             CHKMEMQ
1705812af9f3SBarry Smith   #define PetscStackPush(f)         CHKMEMQ
170627104ee2SJacob Faibussowitsch #endif /* PETSC_USE_DEBUG */
17073a40ed3dSBarry Smith 
17086d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
17093ba16761SJacob Faibussowitsch   #define PetscStackCallExternalVoid(...)
17103ba16761SJacob Faibussowitsch template <typename F, typename... Args>
17113ba16761SJacob Faibussowitsch void PetscCallExternal(F, Args...);
17126d210af2SJacob Faibussowitsch #else
1713586f9135SBarry Smith   /*MC
1714e77caa6dSBarry Smith     PetscStackCallExternalVoid - Calls an external library routine or user function after pushing the name of the routine on the stack.
1715eb6b5d47SBarry Smith 
1716eb6b5d47SBarry Smith    Input Parameters:
1717eb6b5d47SBarry Smith +   name - string that gives the name of the function being called
1718586f9135SBarry Smith -   routine - actual call to the routine, for example, functionname(a,b)
1719fd3f9acdSBarry Smith 
1720586f9135SBarry Smith    Level: developer
1721eb6b5d47SBarry Smith 
1722586f9135SBarry Smith    Note:
1723792fecdfSBarry Smith    Often one should use `PetscCallExternal()` instead. This routine is intended for external library routines that DO NOT return error codes
1724eb6b5d47SBarry Smith 
1725586f9135SBarry Smith    In debug mode this also checks the memory for corruption at the end of the function call.
1726586f9135SBarry Smith 
1727586f9135SBarry Smith    Certain external packages, such as BLAS/LAPACK may have their own macros for managing the call, error checking, etc.
1728586f9135SBarry Smith 
1729586f9135SBarry Smith    Developer Note:
1730586f9135SBarry 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.
1731586f9135SBarry Smith 
1732792fecdfSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscCallExternal()`, `PetscCallBLAS()`
1733586f9135SBarry Smith @*/
17343ba16761SJacob Faibussowitsch   #define PetscStackCallExternalVoid(name, ...) \
17359371c9d4SSatish Balay     do { \
173628770333SStefano Zampini       PetscStackPushExternal(name); \
17373ba16761SJacob Faibussowitsch       __VA_ARGS__; \
17389371c9d4SSatish Balay       PetscStackPop; \
17399371c9d4SSatish Balay     } while (0)
1740eb6b5d47SBarry Smith 
1741586f9135SBarry Smith   /*MC
1742792fecdfSBarry Smith     PetscCallExternal - Calls an external library routine that returns an error code after pushing the name of the routine on the stack.
1743fd3f9acdSBarry Smith 
1744fd3f9acdSBarry Smith    Input Parameters:
1745fd3f9acdSBarry Smith +   func-  name of the routine
1746586f9135SBarry Smith -   args - arguments to the routine
1747586f9135SBarry Smith 
1748586f9135SBarry Smith    Level: developer
1749fd3f9acdSBarry Smith 
175095452b02SPatrick Sanan    Notes:
1751e77caa6dSBarry Smith    This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not.
1752dbf62e16SBarry Smith 
1753586f9135SBarry Smith    In debug mode this also checks the memory for corruption at the end of the function call.
1754fd3f9acdSBarry Smith 
175587497f52SBarry Smith    Assumes the error return code of the function is an integer and that a value of 0 indicates success
175687497f52SBarry Smith 
1757586f9135SBarry Smith    Developer Note:
1758d5b43468SJose E. Roman    This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc.
1759586f9135SBarry Smith 
1760e77caa6dSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()`
1761586f9135SBarry Smith M*/
17629371c9d4SSatish Balay   #define PetscCallExternal(func, ...) \
17639371c9d4SSatish Balay     do { \
1764a74df02fSJacob Faibussowitsch       PetscStackPush(PetscStringize(func)); \
17653ba16761SJacob Faibussowitsch       int ierr_petsc_call_external_ = func(__VA_ARGS__); \
17661d4906efSStefano Zampini       PetscStackPop; \
17673ba16761SJacob 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_); \
1768fd3f9acdSBarry Smith     } while (0)
17696d210af2SJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */
1770