xref: /petsc/include/petscerror.h (revision edd03b47ed9fa9b7ef9dd67df42d5194593a501d)
154a8ef01SBarry Smith /*
2f621e05eSBarry Smith     Contains all error handling interfaces for PETSc.
354a8ef01SBarry Smith */
46524c165SJacob Faibussowitsch #ifndef PETSCERROR_H
526bd1501SBarry Smith #define PETSCERROR_H
66c7e564aSBarry Smith 
7e1bf4ed2SJacob Faibussowitsch #include <petscmacros.h>
8e1bf4ed2SJacob Faibussowitsch #include <petscsystypes.h>
9e1bf4ed2SJacob Faibussowitsch 
10af75f258SJacob Faibussowitsch #if defined(__cplusplus)
11af75f258SJacob Faibussowitsch   #include <exception> // std::exception
12af75f258SJacob Faibussowitsch #endif
13af75f258SJacob Faibussowitsch 
14ac09b921SBarry Smith /* SUBMANSEC = Sys */
15ac09b921SBarry Smith 
16*edd03b47SJacob Faibussowitsch #define SETERRQ1(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
17*edd03b47SJacob Faibussowitsch #define SETERRQ2(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
18*edd03b47SJacob Faibussowitsch #define SETERRQ3(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
19*edd03b47SJacob Faibussowitsch #define SETERRQ4(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
20*edd03b47SJacob Faibussowitsch #define SETERRQ5(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
21*edd03b47SJacob Faibussowitsch #define SETERRQ6(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
22*edd03b47SJacob Faibussowitsch #define SETERRQ7(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
23*edd03b47SJacob Faibussowitsch #define SETERRQ8(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
24*edd03b47SJacob Faibussowitsch #define SETERRQ9(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
2598921bdaSJacob Faibussowitsch 
2630de9b25SBarry Smith /*MC
271957e957SBarry Smith    SETERRQ - Macro to be called when an error has been detected,
2830de9b25SBarry Smith 
2930de9b25SBarry Smith    Synopsis:
30aaa7dc30SBarry Smith    #include <petscsys.h>
3198921bdaSJacob Faibussowitsch    PetscErrorCode SETERRQ(MPI_Comm comm,PetscErrorCode ierr,char *message,...)
3230de9b25SBarry Smith 
33d083f849SBarry Smith    Collective
3430de9b25SBarry Smith 
3530de9b25SBarry Smith    Input Parameters:
3687497f52SBarry Smith +  comm - A communicator, use `PETSC_COMM_SELF` unless you know all ranks of another communicator will detect the error
373af045c5SBarry Smith .  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
3830de9b25SBarry Smith -  message - error message
3930de9b25SBarry Smith 
4030de9b25SBarry Smith   Level: beginner
4130de9b25SBarry Smith 
4230de9b25SBarry Smith    Notes:
4387497f52SBarry Smith     This is rarely needed, one should use `PetscCheck()` and `PetscCall()` and friends to automatically handle error conditions.
4430de9b25SBarry Smith     Once the error handler is called the calling function is then returned from with the given error code.
4530de9b25SBarry Smith 
4687497f52SBarry Smith     Experienced users can set the error handler with `PetscPushErrorHandler()`.
4730de9b25SBarry Smith 
4816a05f60SBarry Smith    Fortran Note:
49989c49a2SBarry Smith    `SETERRQ()` may be called from Fortran subroutines but `SETERRA()` must be called from the
503b1008b8SBarry Smith    Fortran main program.
513b1008b8SBarry Smith 
52db781477SPatrick Sanan .seealso: `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
53db781477SPatrick Sanan           `PetscError()`, `PetscCall()`, `CHKMEMQ`, `CHKERRA()`, `PetscCallMPI()`
5430de9b25SBarry Smith M*/
55e809461dSJacob Faibussowitsch #define SETERRQ(comm, ierr, ...) \
56e809461dSJacob Faibussowitsch   do { \
57e809461dSJacob Faibussowitsch     PetscErrorCode ierr_seterrq_petsc_ = PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \
58e809461dSJacob Faibussowitsch     return ierr_seterrq_petsc_ ? ierr_seterrq_petsc_ : PETSC_ERR_RETURN; \
59e809461dSJacob Faibussowitsch   } while (0)
60986eef2eSBarry Smith 
61ffc4695bSBarry Smith /*
62ffc4695bSBarry Smith     Returned from PETSc functions that are called from MPI, such as related to attributes
63ffc4695bSBarry Smith       Do not confuse PETSC_MPI_ERROR_CODE and PETSC_ERR_MPI, the first is registered with MPI and returned to MPI as
64888a1fb3SBarry 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.
65ffc4695bSBarry Smith */
66ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CLASS;
67ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CODE;
68ffc4695bSBarry Smith 
69986eef2eSBarry Smith /*MC
70986eef2eSBarry Smith    SETERRMPI - Macro to be called when an error has been detected within an MPI callback function
71986eef2eSBarry Smith 
72989c49a2SBarry Smith    No Fortran Support
73989c49a2SBarry Smith 
74986eef2eSBarry Smith    Synopsis:
75986eef2eSBarry Smith    #include <petscsys.h>
7698921bdaSJacob Faibussowitsch    PetscErrorCode SETERRMPI(MPI_Comm comm,PetscErrorCode ierr,char *message,...)
77986eef2eSBarry Smith 
78d083f849SBarry Smith    Collective
79986eef2eSBarry Smith 
80986eef2eSBarry Smith    Input Parameters:
8187497f52SBarry Smith +  comm - A communicator, use `PETSC_COMM_SELF` unless you know all ranks of another communicator will detect the error
82986eef2eSBarry Smith .  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
83986eef2eSBarry Smith -  message - error message
84986eef2eSBarry Smith 
85986eef2eSBarry Smith   Level: developer
86986eef2eSBarry Smith 
87986eef2eSBarry Smith    Notes:
8887497f52SBarry 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`
8987497f52SBarry Smith     which is registered with `MPI_Add_error_code()` when PETSc is initialized.
90986eef2eSBarry Smith 
91db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscCall()`, `PetscCallMPI()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`
92986eef2eSBarry Smith M*/
933ba16761SJacob Faibussowitsch #define SETERRMPI(comm, ierr, ...) return ((void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__), PETSC_MPI_ERROR_CODE)
94ee8199e6SMatthew G. Knepley 
95ee8199e6SMatthew G. Knepley /*MC
96f388eb8bSPatrick Sanan    SETERRA - Fortran-only macro that can be called when an error has been detected from the main program
97f388eb8bSPatrick Sanan 
98f388eb8bSPatrick Sanan    Synopsis:
99f388eb8bSPatrick Sanan    #include <petscsys.h>
100f388eb8bSPatrick Sanan    PetscErrorCode SETERRA(MPI_Comm comm,PetscErrorCode ierr,char *message)
101f388eb8bSPatrick Sanan 
102f388eb8bSPatrick Sanan    Collective
103f388eb8bSPatrick Sanan 
104f388eb8bSPatrick Sanan    Input Parameters:
105f388eb8bSPatrick Sanan +  comm - A communicator, so that the error can be collective
106f388eb8bSPatrick Sanan .  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
107f388eb8bSPatrick Sanan -  message - error message in the printf format
108f388eb8bSPatrick Sanan 
109f388eb8bSPatrick Sanan   Level: beginner
110f388eb8bSPatrick Sanan 
111f388eb8bSPatrick Sanan    Notes:
11287497f52SBarry Smith    This should only be used with Fortran. With C/C++, use `SETERRQ()`.
113f388eb8bSPatrick Sanan 
11487497f52SBarry Smith    `SETERRQ()` may be called from Fortran subroutines but `SETERRA()` must be called from the
115f388eb8bSPatrick Sanan     Fortran main program.
116f388eb8bSPatrick Sanan 
117db781477SPatrick Sanan .seealso: `SETERRQ()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()`
118f388eb8bSPatrick Sanan M*/
119f388eb8bSPatrick Sanan 
120f388eb8bSPatrick Sanan /*MC
121fa190f98SMatthew G. Knepley    SETERRABORT - Macro that can be called when an error has been detected,
122fa190f98SMatthew G. Knepley 
123fa190f98SMatthew G. Knepley    Synopsis:
124fa190f98SMatthew G. Knepley    #include <petscsys.h>
12598921bdaSJacob Faibussowitsch    PetscErrorCode SETERRABORT(MPI_Comm comm,PetscErrorCode ierr,char *message,...)
126fa190f98SMatthew G. Knepley 
127d083f849SBarry Smith    Collective
128fa190f98SMatthew G. Knepley 
129fa190f98SMatthew G. Knepley    Input Parameters:
130fa190f98SMatthew G. Knepley +  comm - A communicator, so that the error can be collective
1313af045c5SBarry Smith .  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
132fa190f98SMatthew G. Knepley -  message - error message in the printf format
133fa190f98SMatthew G. Knepley 
134fa190f98SMatthew G. Knepley   Level: beginner
135fa190f98SMatthew G. Knepley 
136fa190f98SMatthew G. Knepley    Notes:
13787497f52SBarry Smith    This function just calls `MPI_Abort()`.
13887497f52SBarry Smith 
13987497f52SBarry Smith    This should only be called in routines that cannot return an error code, such as in C++ constructors.
140fa190f98SMatthew G. Knepley 
141989c49a2SBarry Smith    Fortran Note:
142989c49a2SBarry Smith    Use `SETERRA()` in Fortran main program and `SETERRQ()` in Fortran subroutines
143989c49a2SBarry Smith 
144989c49a2SBarry Smith    Developer Note:
145989c49a2SBarry Smith    In Fortran `SETERRA()` could be called `SETERRABORT()` since they serve the same purpose
146989c49a2SBarry Smith 
147db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `PetscCall()`, `CHKMEMQ`
148fa190f98SMatthew G. Knepley M*/
1499371c9d4SSatish Balay #define SETERRABORT(comm, ierr, ...) \
1509371c9d4SSatish Balay   do { \
1513ba16761SJacob Faibussowitsch     (void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \
15298921bdaSJacob Faibussowitsch     MPI_Abort(comm, ierr); \
15398921bdaSJacob Faibussowitsch   } while (0)
1549a00fa46SSatish Balay 
15530de9b25SBarry Smith /*MC
1562c71b3e2SJacob Faibussowitsch   PetscCheck - Check that a particular condition is true
1572c71b3e2SJacob Faibussowitsch 
1582c71b3e2SJacob Faibussowitsch   Synopsis:
1592c71b3e2SJacob Faibussowitsch   #include <petscerror.h>
1602c71b3e2SJacob Faibussowitsch   void PetscCheck(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
1612c71b3e2SJacob Faibussowitsch 
1627cdbe19fSJose E. Roman   Collective; No Fortran Support
1632c71b3e2SJacob Faibussowitsch 
1642c71b3e2SJacob Faibussowitsch   Input Parameters:
1652c71b3e2SJacob Faibussowitsch + cond    - The boolean condition
1662c71b3e2SJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
1672c71b3e2SJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
1682c71b3e2SJacob Faibussowitsch - message - Error message in printf format
1692c71b3e2SJacob Faibussowitsch 
170667f096bSBarry Smith   Level: beginner
171667f096bSBarry Smith 
1722c71b3e2SJacob Faibussowitsch   Notes:
1732c71b3e2SJacob Faibussowitsch   Enabled in both optimized and debug builds.
1742c71b3e2SJacob Faibussowitsch 
17587497f52SBarry Smith   Calls `SETERRQ()` if the assertion fails, so can only be called from functions returning a
17687497f52SBarry Smith   `PetscErrorCode` (or equivalent type after conversion).
1772c71b3e2SJacob Faibussowitsch 
1784be741a6SBarry Smith  .seealso: `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheckAbort()`
1799566063dSJacob Faibussowitsch M*/
1809371c9d4SSatish Balay #define PetscCheck(cond, comm, ierr, ...) \
1819371c9d4SSatish Balay   do { \
1829371c9d4SSatish Balay     if (PetscUnlikely(!(cond))) SETERRQ(comm, ierr, __VA_ARGS__); \
1839371c9d4SSatish Balay   } while (0)
1842c71b3e2SJacob Faibussowitsch 
1852c71b3e2SJacob Faibussowitsch /*MC
1864be741a6SBarry Smith   PetscCheckAbort - Check that a particular condition is true, otherwise prints error and aborts
1874be741a6SBarry Smith 
1884be741a6SBarry Smith   Synopsis:
1894be741a6SBarry Smith   #include <petscerror.h>
1904be741a6SBarry Smith   void PetscCheckAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
1914be741a6SBarry Smith 
1927cdbe19fSJose E. Roman   Collective; No Fortran Support
1934be741a6SBarry Smith 
1944be741a6SBarry Smith   Input Parameters:
1954be741a6SBarry Smith + cond    - The boolean condition
1964be741a6SBarry Smith . comm    - The communicator on which the check can be collective on
1974be741a6SBarry Smith . ierr    - A nonzero error code, see include/petscerror.h for the complete list
1984be741a6SBarry Smith - message - Error message in printf format
1994be741a6SBarry Smith 
200667f096bSBarry Smith   Level: developer
201667f096bSBarry Smith 
2024be741a6SBarry Smith   Notes:
2034be741a6SBarry Smith   Enabled in both optimized and debug builds.
2044be741a6SBarry Smith 
20587497f52SBarry Smith   Calls `SETERRABORT()` if the assertion fails, can be called from a function that does not return an
20687497f52SBarry Smith   error code, such as a C++ constructor. usually `PetscCheck()` should be used.
2074be741a6SBarry Smith 
208989c49a2SBarry Smith .seealso: `PetscAssertAbort()`, `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheck()`, `SETERRABORT()`
2094be741a6SBarry Smith M*/
2109371c9d4SSatish Balay #define PetscCheckAbort(cond, comm, ierr, ...) \
2110e6b6b59SJacob Faibussowitsch   do { \
2120e6b6b59SJacob Faibussowitsch     if (PetscUnlikely(!(cond))) SETERRABORT(comm, ierr, __VA_ARGS__); \
213c7481402SJacob Faibussowitsch   } while (0)
2144be741a6SBarry Smith 
2154be741a6SBarry Smith /*MC
2169ace16cdSJacob Faibussowitsch   PetscAssert - Assert that a particular condition is true
2179ace16cdSJacob Faibussowitsch 
2189ace16cdSJacob Faibussowitsch   Synopsis:
2199ace16cdSJacob Faibussowitsch   #include <petscerror.h>
2209ace16cdSJacob Faibussowitsch   void PetscAssert(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
2219ace16cdSJacob Faibussowitsch 
2227cdbe19fSJose E. Roman   Collective; No Fortran Support
2239ace16cdSJacob Faibussowitsch 
2249ace16cdSJacob Faibussowitsch   Input Parameters:
2259ace16cdSJacob Faibussowitsch + cond    - The boolean condition
2269ace16cdSJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
2279ace16cdSJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
2289ace16cdSJacob Faibussowitsch - message - Error message in printf format
2299ace16cdSJacob Faibussowitsch 
230667f096bSBarry Smith   Level: beginner
231667f096bSBarry Smith 
2329ace16cdSJacob Faibussowitsch   Notes:
233c7481402SJacob Faibussowitsch   Equivalent to `PetscCheck()` if debugging is enabled, and `PetscAssume(cond)` otherwise.
2342c71b3e2SJacob Faibussowitsch 
23587497f52SBarry Smith   See `PetscCheck()` for usage and behaviour.
23687497f52SBarry Smith 
23787497f52SBarry Smith   This is needed instead of simply using `assert()` because this correctly handles the collective nature of errors under MPI
2389ace16cdSJacob Faibussowitsch 
2390e6b6b59SJacob Faibussowitsch .seealso: `PetscCheck()`, `SETERRQ()`, `PetscError()`, `PetscAssertAbort()`
2409566063dSJacob Faibussowitsch M*/
241c7481402SJacob Faibussowitsch #if PetscDefined(USE_DEBUG)
242c7481402SJacob Faibussowitsch   #define PetscAssert(cond, comm, ierr, ...) PetscCheck(cond, comm, ierr, __VA_ARGS__)
243c7481402SJacob Faibussowitsch #else
244c7481402SJacob Faibussowitsch   #define PetscAssert(cond, ...) PetscAssume(cond)
245c7481402SJacob Faibussowitsch #endif
2469ace16cdSJacob Faibussowitsch 
2479ace16cdSJacob Faibussowitsch /*MC
2480e6b6b59SJacob Faibussowitsch   PetscAssertAbort - Assert that a particular condition is true, otherwise prints error and aborts
2490e6b6b59SJacob Faibussowitsch 
2500e6b6b59SJacob Faibussowitsch   Synopsis:
2510e6b6b59SJacob Faibussowitsch   #include <petscerror.h>
2520e6b6b59SJacob Faibussowitsch   void PetscAssertAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
2530e6b6b59SJacob Faibussowitsch 
2547cdbe19fSJose E. Roman   Collective; No Fortran Support
2550e6b6b59SJacob Faibussowitsch 
2560e6b6b59SJacob Faibussowitsch   Input Parameters:
2570e6b6b59SJacob Faibussowitsch + cond    - The boolean condition
2580e6b6b59SJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
2590e6b6b59SJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
2600e6b6b59SJacob Faibussowitsch - message - Error message in printf format
2610e6b6b59SJacob Faibussowitsch 
262667f096bSBarry Smith   Level: beginner
263667f096bSBarry Smith 
2640e6b6b59SJacob Faibussowitsch   Notes:
2650e6b6b59SJacob Faibussowitsch   Enabled only in debug builds. See `PetscCheckAbort()` for usage.
2660e6b6b59SJacob Faibussowitsch 
2670e6b6b59SJacob Faibussowitsch .seealso: `PetscCheckAbort()`, `PetscAssert()`, `PetscCheck()`, `SETERRABORT()`, `PetscError()`
2680e6b6b59SJacob Faibussowitsch M*/
2693ba16761SJacob Faibussowitsch #if PetscDefined(USE_DEBUG)
2703ba16761SJacob Faibussowitsch   #define PetscAssertAbort(cond, comm, ierr, ...) PetscCheckAbort(cond, comm, ierr, __VA_ARGS__)
2713ba16761SJacob Faibussowitsch #else
2723ba16761SJacob Faibussowitsch   #define PetscAssertAbort(cond, comm, ierr, ...) PetscAssume(cond)
2733ba16761SJacob Faibussowitsch #endif
2740e6b6b59SJacob Faibussowitsch 
2750e6b6b59SJacob Faibussowitsch /*MC
2763ba16761SJacob Faibussowitsch   PetscCall - Calls a PETSc function and then checks the resulting error code, if it is
2773ba16761SJacob Faibussowitsch   non-zero it calls the error handler and returns from the current function with the error
2783ba16761SJacob Faibussowitsch   code.
2799566063dSJacob Faibussowitsch 
2809566063dSJacob Faibussowitsch   Synopsis:
2819566063dSJacob Faibussowitsch   #include <petscerror.h>
28249c86fc7SBarry Smith   void PetscCall(PetscFunction(args))
2839566063dSJacob Faibussowitsch 
2849566063dSJacob Faibussowitsch   Not Collective
2859566063dSJacob Faibussowitsch 
2869566063dSJacob Faibussowitsch   Input Parameter:
28749c86fc7SBarry Smith . PetscFunction - any PETSc function that returns an error code
2889566063dSJacob Faibussowitsch 
289667f096bSBarry Smith   Level: beginner
290667f096bSBarry Smith 
2919566063dSJacob Faibussowitsch   Notes:
2929566063dSJacob Faibussowitsch   Once the error handler is called the calling function is then returned from with the given
29387497f52SBarry Smith   error code. Experienced users can set the error handler with `PetscPushErrorHandler()`.
2949566063dSJacob Faibussowitsch 
29587497f52SBarry Smith   `PetscCall()` cannot be used in functions returning a datatype not convertible to
29687497f52SBarry Smith   `PetscErrorCode`. For example, `PetscCall()` may not be used in functions returning void, use
2973ba16761SJacob Faibussowitsch   `PetscCallAbort()` or `PetscCallVoid()` in this case.
2989566063dSJacob Faibussowitsch 
2999566063dSJacob Faibussowitsch   Example Usage:
3009566063dSJacob Faibussowitsch .vb
3019566063dSJacob Faibussowitsch   PetscCall(PetscInitiailize(...)); // OK to call even when PETSc is not yet initialized!
3029566063dSJacob Faibussowitsch 
3039566063dSJacob Faibussowitsch   struct my_struct
3049566063dSJacob Faibussowitsch   {
3059566063dSJacob Faibussowitsch     void *data;
3069566063dSJacob Faibussowitsch   } my_complex_type;
3079566063dSJacob Faibussowitsch 
3089566063dSJacob Faibussowitsch   struct my_struct bar(void)
3099566063dSJacob Faibussowitsch   {
3106aad120cSJose E. Roman     PetscCall(foo(15)); // ERROR PetscErrorCode not convertible to struct my_struct!
3119566063dSJacob Faibussowitsch   }
3129566063dSJacob Faibussowitsch 
3136aad120cSJose E. Roman   PetscCall(bar()) // ERROR input not convertible to PetscErrorCode
3149566063dSJacob Faibussowitsch .ve
3159566063dSJacob Faibussowitsch 
31687497f52SBarry Smith   It is also possible to call this directly on a `PetscErrorCode` variable
31749c86fc7SBarry Smith .vb
31849c86fc7SBarry Smith   PetscCall(ierr);  // check if ierr is nonzero
31949c86fc7SBarry Smith .ve
32049c86fc7SBarry Smith 
321792fecdfSBarry Smith   Should not be used to call callback functions provided by users, `PetscCallBack()` should be used in that situation.
322ef1023bdSBarry Smith 
3236a8be23eSBarry Smith   `PetscUseTypeMethod()` or `PetscTryTypeMethod()` should be used when calling functions pointers contained in a PETSc object's `ops` array
3246a8be23eSBarry Smith 
32549c86fc7SBarry Smith   Fortran Notes:
32649c86fc7SBarry Smith     The Fortran function from which this is used must declare a variable PetscErrorCode ierr and ierr must be
32787497f52SBarry Smith     the final argument to the PETSc function being called.
32849c86fc7SBarry Smith 
32949c86fc7SBarry Smith     In the main program and in Fortran subroutines that do not have ierr as the final return parameter one
33087497f52SBarry Smith     should use `PetscCallA()`
33149c86fc7SBarry Smith 
33249c86fc7SBarry Smith   Example Fortran Usage:
33349c86fc7SBarry Smith .vb
33449c86fc7SBarry Smith   PetscErrorCode ierr
33549c86fc7SBarry Smith   Vec v
33649c86fc7SBarry Smith 
33749c86fc7SBarry Smith   ...
33849c86fc7SBarry Smith   PetscCall(VecShift(v,1.0,ierr))
33949c86fc7SBarry Smith   PetscCallA(VecShift(v,1.0,ierr))
34049c86fc7SBarry Smith .ve
34149c86fc7SBarry Smith 
342667f096bSBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`,
343667f096bSBarry Smith           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`,
3443ba16761SJacob Faibussowitsch           `CHKERRMPI()`, `PetscCallBack()`, `PetscCallAbort()`, `PetscCallVoid()`
3459566063dSJacob Faibussowitsch M*/
346ef1023bdSBarry Smith 
347ef1023bdSBarry Smith /*MC
348989c49a2SBarry Smith    PetscCallA - Fortran-only macro that should be used in the main program to call PETSc functions instead of using
349989c49a2SBarry Smith    PetscCall() which should be used in other Fortran subroutines
350989c49a2SBarry Smith 
351989c49a2SBarry Smith    Synopsis:
352989c49a2SBarry Smith    #include <petscsys.h>
353989c49a2SBarry Smith    PetscErrorCode PetscCallA(PetscFunction(arguments,ierr))
354989c49a2SBarry Smith 
355989c49a2SBarry Smith    Collective
356989c49a2SBarry Smith 
357989c49a2SBarry Smith    Input Parameter:
358989c49a2SBarry Smith .  PetscFunction(arguments,ierr) - the call to the function
359989c49a2SBarry Smith 
360989c49a2SBarry Smith   Level: beginner
361989c49a2SBarry Smith 
362989c49a2SBarry Smith    Notes:
363989c49a2SBarry Smith    This should only be used with Fortran. With C/C++, use `PetscCall()` always.
364989c49a2SBarry Smith 
365989c49a2SBarry Smith    Use `SETERRA()` to set an error in a Fortran main program and `SETERRQ()` in Fortran subroutines
366989c49a2SBarry Smith 
367989c49a2SBarry Smith .seealso: `SETERRQ()`, `SETERRA()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()`
368989c49a2SBarry Smith M*/
369989c49a2SBarry Smith 
370989c49a2SBarry Smith /*MC
371792fecdfSBarry 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
372ef1023bdSBarry Smith   handler and returns from the current function with the error code.
373ef1023bdSBarry Smith 
374ef1023bdSBarry Smith   Synopsis:
375ef1023bdSBarry Smith   #include <petscerror.h>
376792fecdfSBarry Smith   void PetscCallBack(const char *functionname,PetscFunction(args))
377ef1023bdSBarry Smith 
3787cdbe19fSJose E. Roman   Not Collective; No Fortran Support
379ef1023bdSBarry Smith 
380ef1023bdSBarry Smith   Input Parameters:
381ef1023bdSBarry Smith + functionname - the name of the function being called, this can be a string with spaces that describes the meaning of the callback
38287497f52SBarry Smith - PetscFunction - user provided callback function that returns an error code
383ef1023bdSBarry Smith 
384ef1023bdSBarry Smith   Example Usage:
385ef1023bdSBarry Smith .vb
386792fecdfSBarry Smith   PetscCallBack("XXX callback to do something",a->callback(...));
387ef1023bdSBarry Smith .ve
388ef1023bdSBarry Smith 
389ef1023bdSBarry Smith   Level: developer
390ef1023bdSBarry Smith 
391667f096bSBarry Smith   Notes:
392667f096bSBarry Smith   Once the error handler is called the calling function is then returned from with the given
393667f096bSBarry Smith   error code. Experienced users can set the error handler with `PetscPushErrorHandler()`.
394667f096bSBarry Smith 
395667f096bSBarry Smith   `PetscCallBack()` should only be called in PETSc when a call is being made to a user provided call-back routine.
396667f096bSBarry Smith 
39787497f52SBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscCall()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`
398ef1023bdSBarry Smith           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, `CHKERRMPI()`, `PetscCall()`
399ef1023bdSBarry Smith M*/
400ef1023bdSBarry Smith 
4013ba16761SJacob Faibussowitsch /*MC
4023ba16761SJacob Faibussowitsch   PetscCallVoid - Like `PetscCall()` but for functions returning `void`
4033ba16761SJacob Faibussowitsch 
4043ba16761SJacob Faibussowitsch   Synopsis:
4053ba16761SJacob Faibussowitsch   #include <petscerror.h>
4063ba16761SJacob Faibussowitsch   void PetscCall(PetscFunction(args))
4073ba16761SJacob Faibussowitsch 
4087cdbe19fSJose E. Roman   Not Collective; No Fortran Support
4093ba16761SJacob Faibussowitsch 
4103ba16761SJacob Faibussowitsch   Input Parameter:
4113ba16761SJacob Faibussowitsch . PetscFunction - any PETSc function that returns an error code
4123ba16761SJacob Faibussowitsch 
4133ba16761SJacob Faibussowitsch   Example Usage:
4143ba16761SJacob Faibussowitsch .vb
4153ba16761SJacob Faibussowitsch   void foo()
4163ba16761SJacob Faibussowitsch   {
4173ba16761SJacob Faibussowitsch     KSP ksp;
4183ba16761SJacob Faibussowitsch 
4193ba16761SJacob Faibussowitsch     PetscFunctionBeginUser;
4203ba16761SJacob Faibussowitsch     // OK, properly handles PETSc error codes
4213ba16761SJacob Faibussowitsch     PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp));
4223ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
4233ba16761SJacob Faibussowitsch   }
4243ba16761SJacob Faibussowitsch 
4253ba16761SJacob Faibussowitsch   PetscErrorCode bar()
4263ba16761SJacob Faibussowitsch   {
4273ba16761SJacob Faibussowitsch     KSP ksp;
4283ba16761SJacob Faibussowitsch 
4293ba16761SJacob Faibussowitsch     PetscFunctionBeginUser;
4303ba16761SJacob Faibussowitsch     // ERROR, Non-void function 'bar' should return a value
4313ba16761SJacob Faibussowitsch     PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp));
4323ba16761SJacob Faibussowitsch     // OK, returning PetscErrorCode
4333ba16761SJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
4343ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
4353ba16761SJacob Faibussowitsch   }
436667f096bSBarry Smith .ve
4373ba16761SJacob Faibussowitsch 
4383ba16761SJacob Faibussowitsch   Level: beginner
4393ba16761SJacob Faibussowitsch 
440667f096bSBarry Smith   Notes:
441667f096bSBarry Smith   Has identical usage to `PetscCall()`, except that it returns `void` on error instead of a
442667f096bSBarry Smith   `PetscErrorCode`. See `PetscCall()` for more detailed discussion.
443667f096bSBarry Smith 
444667f096bSBarry Smith   Note that users should prefer `PetscCallAbort()` to this routine. While this routine does
445667f096bSBarry Smith   "handle" errors by returning from the enclosing function, it effectively gobbles the
446667f096bSBarry Smith   error. Since the enclosing function itself returns `void`, its callers have no way of knowing
447667f096bSBarry Smith   that the routine returned early due to an error. `PetscCallAbort()` at least ensures that the
448667f096bSBarry Smith   program crashes gracefully.
449667f096bSBarry Smith 
4503ba16761SJacob Faibussowitsch .seealso: `PetscCall()`, `PetscErrorCode`
4513ba16761SJacob Faibussowitsch M*/
4529566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
4539566063dSJacob Faibussowitsch void PetscCall(PetscErrorCode);
454792fecdfSBarry Smith void PetscCallBack(const char *, PetscErrorCode);
4559566063dSJacob Faibussowitsch void PetscCallVoid(PetscErrorCode);
4569566063dSJacob Faibussowitsch #else
4579371c9d4SSatish Balay   #define PetscCall(...) \
4589371c9d4SSatish Balay     do { \
4593ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_q_; \
46037154d10SBarry Smith       PetscStackUpdateLine; \
4613ba16761SJacob Faibussowitsch       ierr_petsc_call_q_ = __VA_ARGS__; \
4623ba16761SJacob 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, " "); \
4639566063dSJacob Faibussowitsch     } while (0)
4649371c9d4SSatish Balay   #define PetscCallBack(function, ...) \
4659371c9d4SSatish Balay     do { \
4663ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_q_; \
467e33ced7fSLisandro Dalcin       PetscStackUpdateLine; \
468e33ced7fSLisandro Dalcin       PetscStackPushExternal(function); \
4693ba16761SJacob Faibussowitsch       ierr_petsc_call_q_ = __VA_ARGS__; \
470ef1023bdSBarry Smith       PetscStackPop; \
4713ba16761SJacob 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, " "); \
472ef1023bdSBarry Smith     } while (0)
4739371c9d4SSatish Balay   #define PetscCallVoid(...) \
4749371c9d4SSatish Balay     do { \
4753ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_void_; \
476e33ced7fSLisandro Dalcin       PetscStackUpdateLine; \
4773ba16761SJacob Faibussowitsch       ierr_petsc_call_void_ = __VA_ARGS__; \
4783ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_void_ != PETSC_SUCCESS)) { \
4793ba16761SJacob Faibussowitsch         ierr_petsc_call_void_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_void_, PETSC_ERROR_REPEAT, " "); \
4803ba16761SJacob Faibussowitsch         (void)ierr_petsc_call_void_; \
4819371c9d4SSatish Balay         return; \
4829371c9d4SSatish Balay       } \
4839566063dSJacob Faibussowitsch     } while (0)
4849566063dSJacob Faibussowitsch #endif
4859566063dSJacob Faibussowitsch 
4869566063dSJacob Faibussowitsch /*MC
4879566063dSJacob Faibussowitsch   CHKERRQ - Checks error code returned from PETSc function
48830de9b25SBarry Smith 
48930de9b25SBarry Smith   Synopsis:
490aaa7dc30SBarry Smith   #include <petscsys.h>
4919566063dSJacob Faibussowitsch   void CHKERRQ(PetscErrorCode ierr)
4929566063dSJacob Faibussowitsch 
4939566063dSJacob Faibussowitsch   Not Collective
4949566063dSJacob Faibussowitsch 
4952fe279fdSBarry Smith   Input Parameter:
4969566063dSJacob Faibussowitsch . ierr - nonzero error code
4979566063dSJacob Faibussowitsch 
4989566063dSJacob Faibussowitsch   Level: deprecated
4999566063dSJacob Faibussowitsch 
500667f096bSBarry Smith   Note:
501667f096bSBarry Smith   Deprecated in favor of `PetscCall()`. This routine behaves identically to it.
502667f096bSBarry Smith 
503db781477SPatrick Sanan .seealso: `PetscCall()`
5049566063dSJacob Faibussowitsch M*/
5059566063dSJacob Faibussowitsch #define CHKERRQ(...) PetscCall(__VA_ARGS__)
5069566063dSJacob Faibussowitsch #define CHKERRV(...) PetscCallVoid(__VA_ARGS__)
5079566063dSJacob Faibussowitsch 
508db9cea48SBarry Smith PETSC_EXTERN void PetscMPIErrorString(PetscMPIInt, char *);
509db9cea48SBarry Smith 
5109566063dSJacob Faibussowitsch /*MC
5119566063dSJacob Faibussowitsch   PetscCallMPI - Checks error code returned from MPI calls, if non-zero it calls the error
5129566063dSJacob Faibussowitsch   handler and then returns
5139566063dSJacob Faibussowitsch 
5149566063dSJacob Faibussowitsch   Synopsis:
5159566063dSJacob Faibussowitsch   #include <petscerror.h>
51649c86fc7SBarry Smith   void PetscCallMPI(MPI_Function(args))
51730de9b25SBarry Smith 
518eca87e8dSBarry Smith   Not Collective
51930de9b25SBarry Smith 
5202fe279fdSBarry Smith   Input Parameter:
52149c86fc7SBarry Smith . MPI_Function - an MPI function that returns an MPI error code
52230de9b25SBarry Smith 
523667f096bSBarry Smith   Level: beginner
524667f096bSBarry Smith 
5259566063dSJacob Faibussowitsch   Notes:
52687497f52SBarry Smith   Always returns the error code `PETSC_ERR_MPI`; the MPI error code and string are embedded in
5279566063dSJacob Faibussowitsch   the string error message. Do not use this to call any other routines (for example PETSc
5283ba16761SJacob Faibussowitsch   routines), it should only be used for direct MPI calls. The user may configure PETSc with the
5293ba16761SJacob Faibussowitsch   `--with-strict-petscerrorcode` option to check this at compile-time, otherwise they must
5309566063dSJacob Faibussowitsch   check this themselves.
5319566063dSJacob Faibussowitsch 
532aaa8cc7dSPierre Jolivet   This routine can only be used in functions returning `PetscErrorCode` themselves. If the
5333ba16761SJacob Faibussowitsch   calling function returns a different type, use `PetscCallMPIAbort()` instead.
5343ba16761SJacob Faibussowitsch 
5359566063dSJacob Faibussowitsch   Example Usage:
5369566063dSJacob Faibussowitsch .vb
5379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(...)); // OK, calling MPI function
5389566063dSJacob Faibussowitsch 
5399566063dSJacob Faibussowitsch   PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead!
5409566063dSJacob Faibussowitsch .ve
5419566063dSJacob Faibussowitsch 
54249c86fc7SBarry Smith   Fortran Notes:
54387497f52SBarry Smith     The Fortran function from which this is used must declare a variable `PetscErrorCode` ierr and ierr must be
54449c86fc7SBarry Smith     the final argument to the MPI function being called.
54549c86fc7SBarry Smith 
54649c86fc7SBarry Smith     In the main program and in Fortran subroutines that do not have ierr as the final return parameter one
54787497f52SBarry Smith     should use `PetscCallMPIA()`
54849c86fc7SBarry Smith 
54949c86fc7SBarry Smith   Fortran Usage:
55049c86fc7SBarry Smith .vb
55149c86fc7SBarry Smith   PetscErrorCode ierr or integer ierr
55249c86fc7SBarry Smith   ...
55349c86fc7SBarry Smith   PetscCallMPI(MPI_Comm_size(...,ierr))
55449c86fc7SBarry Smith   PetscCallMPIA(MPI_Comm_size(...,ierr)) ! Will abort after calling error handler
55549c86fc7SBarry Smith 
55649c86fc7SBarry Smith   PetscCallMPI(MPI_Comm_size(...,eflag)) ! ERROR, final argument must be ierr
55749c86fc7SBarry Smith .ve
55849c86fc7SBarry Smith 
559db781477SPatrick Sanan .seealso: `SETERRMPI()`, `PetscCall()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`,
5603ba16761SJacob Faibussowitsch           `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
5613ba16761SJacob Faibussowitsch           `PetscError()`, `CHKMEMQ`
5623ba16761SJacob Faibussowitsch M*/
5633ba16761SJacob Faibussowitsch 
5643ba16761SJacob Faibussowitsch /*MC
5653ba16761SJacob Faibussowitsch   PetscCallMPIAbort - Like `PetscCallMPI()` but calls `MPI_Abort()` on error
5663ba16761SJacob Faibussowitsch 
5673ba16761SJacob Faibussowitsch   Synopsis:
5683ba16761SJacob Faibussowitsch   #include <petscerror.h>
5693ba16761SJacob Faibussowitsch   void PetscCallMPIAbort(MPI_Comm comm, MPI_Function(args))
5703ba16761SJacob Faibussowitsch 
5713ba16761SJacob Faibussowitsch   Not Collective
5723ba16761SJacob Faibussowitsch 
5733ba16761SJacob Faibussowitsch   Input Parameters:
5743ba16761SJacob Faibussowitsch + comm         - the MPI communicator to abort on
5753ba16761SJacob Faibussowitsch - MPI_Function - an MPI function that returns an MPI error code
5763ba16761SJacob Faibussowitsch 
577667f096bSBarry Smith   Level: beginner
578667f096bSBarry Smith 
5793ba16761SJacob Faibussowitsch   Notes:
5803ba16761SJacob Faibussowitsch   Usage is identical to `PetscCallMPI()`. See `PetscCallMPI()` for detailed discussion.
5813ba16761SJacob Faibussowitsch 
5823ba16761SJacob Faibussowitsch   This routine may be used in functions returning `void` or other non-`PetscErrorCode` types.
5833ba16761SJacob Faibussowitsch 
584989c49a2SBarry Smith   Fortran Note:
585989c49a2SBarry Smith   In Fortran this is called `PetscCallMPIA()` and is intended to be used in the main program while `PetscCallMPI()` is
586989c49a2SBarry Smith   used in Fortran subroutines.
587989c49a2SBarry Smith 
588989c49a2SBarry Smith   Developer Note:
589989c49a2SBarry Smith   This should have the same name in Fortran.
590989c49a2SBarry Smith 
5913ba16761SJacob Faibussowitsch .seealso: `PetscCallMPI()`, `PetscCallAbort()`, `SETERRABORT()`
59230de9b25SBarry Smith M*/
5933fcd9f07SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
59463a3b9bcSJacob Faibussowitsch void PetscCallMPI(PetscMPIInt);
5953ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm, PetscMPIInt);
596064a246eSJacob Faibussowitsch #else
5973ba16761SJacob Faibussowitsch   #define PetscCallMPI_Private(__PETSC_STACK_POP_FUNC__, __SETERR_FUNC__, __COMM__, ...) \
5989371c9d4SSatish Balay     do { \
5993ba16761SJacob Faibussowitsch       PetscMPIInt ierr_petsc_call_mpi_; \
600ef1023bdSBarry Smith       PetscStackUpdateLine; \
601792fecdfSBarry Smith       PetscStackPushExternal("MPI function"); \
602d71ae5a4SJacob Faibussowitsch       { \
6033ba16761SJacob Faibussowitsch         ierr_petsc_call_mpi_ = __VA_ARGS__; \
604d71ae5a4SJacob Faibussowitsch       } \
6053ba16761SJacob Faibussowitsch       __PETSC_STACK_POP_FUNC__; \
6063ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_mpi_ != MPI_SUCCESS)) { \
6073ba16761SJacob Faibussowitsch         char petsc_mpi_7_errorstring[2 * MPI_MAX_ERROR_STRING]; \
6083ba16761SJacob Faibussowitsch         PetscMPIErrorString(ierr_petsc_call_mpi_, (char *)petsc_mpi_7_errorstring); \
6093ba16761SJacob Faibussowitsch         __SETERR_FUNC__(__COMM__, PETSC_ERR_MPI, "MPI error %d %s", (int)ierr_petsc_call_mpi_, petsc_mpi_7_errorstring); \
6103fcd9f07SJacob Faibussowitsch       } \
6113fcd9f07SJacob Faibussowitsch     } while (0)
6123ba16761SJacob Faibussowitsch 
6133ba16761SJacob Faibussowitsch   #define PetscCallMPI(...)            PetscCallMPI_Private(PetscStackPop, SETERRQ, PETSC_COMM_SELF, __VA_ARGS__)
6143ba16761SJacob Faibussowitsch   #define PetscCallMPIAbort(comm, ...) PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRABORT, comm, __VA_ARGS__)
6159566063dSJacob Faibussowitsch #endif
616064a246eSJacob Faibussowitsch 
6177037b0edSPatrick Sanan /*MC
6189566063dSJacob Faibussowitsch   CHKERRMPI - Checks error code returned from MPI calls, if non-zero it calls the error
6199566063dSJacob Faibussowitsch   handler and then returns
6209566063dSJacob Faibussowitsch 
6219566063dSJacob Faibussowitsch   Synopsis:
6229566063dSJacob Faibussowitsch   #include <petscerror.h>
6239566063dSJacob Faibussowitsch   void CHKERRMPI(PetscErrorCode ierr)
6249566063dSJacob Faibussowitsch 
6259566063dSJacob Faibussowitsch   Not Collective
6269566063dSJacob Faibussowitsch 
6279566063dSJacob Faibussowitsch   Input Parameter:
6289566063dSJacob Faibussowitsch . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
6299566063dSJacob Faibussowitsch 
6309566063dSJacob Faibussowitsch   Level: deprecated
6319566063dSJacob Faibussowitsch 
632667f096bSBarry Smith   Note:
633667f096bSBarry Smith   Deprecated in favor of `PetscCallMPI()`. This routine behaves identically to it.
634667f096bSBarry Smith 
635db781477SPatrick Sanan .seealso: `PetscCallMPI()`
6369566063dSJacob Faibussowitsch M*/
6379566063dSJacob Faibussowitsch #define CHKERRMPI(...) PetscCallMPI(__VA_ARGS__)
6389566063dSJacob Faibussowitsch 
6399566063dSJacob Faibussowitsch /*MC
640989c49a2SBarry Smith   PetscCallAbort - Checks error code returned from PETSc function, if non-zero it aborts immediately by calling `MPI_Abort()`
6419566063dSJacob Faibussowitsch 
6429566063dSJacob Faibussowitsch   Synopsis:
6439566063dSJacob Faibussowitsch   #include <petscerror.h>
6449566063dSJacob Faibussowitsch   void PetscCallAbort(MPI_Comm comm, PetscErrorCode ierr)
6459566063dSJacob Faibussowitsch 
646c3339decSBarry Smith   Collective
6479566063dSJacob Faibussowitsch 
6489566063dSJacob Faibussowitsch   Input Parameters:
6499566063dSJacob Faibussowitsch + comm - the MPI communicator on which to abort
6509566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
6519566063dSJacob Faibussowitsch 
652667f096bSBarry Smith   Level: intermediate
653667f096bSBarry Smith 
6549566063dSJacob Faibussowitsch   Notes:
65587497f52SBarry Smith   This macro has identical type and usage semantics to `PetscCall()` with the important caveat
6569566063dSJacob Faibussowitsch   that this macro does not return. Instead, if ierr is nonzero it calls the PETSc error handler
65787497f52SBarry Smith   and then immediately calls `MPI_Abort()`. It can therefore be used anywhere.
6589566063dSJacob Faibussowitsch 
65987497f52SBarry Smith   As per `MPI_Abort()` semantics the communicator passed must be valid, although there is currently
66087497f52SBarry Smith   no attempt made at handling any potential errors from `MPI_Abort()`. Note that while
66187497f52SBarry Smith   `MPI_Abort()` is required to terminate only those processes which reside on comm, it is often
66287497f52SBarry Smith   the case that `MPI_Abort()` terminates *all* processes.
6639566063dSJacob Faibussowitsch 
6649566063dSJacob Faibussowitsch   Example Usage:
6659566063dSJacob Faibussowitsch .vb
6669566063dSJacob Faibussowitsch   PetscErrorCode boom(void) { return PETSC_ERR_MEM; }
6679566063dSJacob Faibussowitsch 
6689566063dSJacob Faibussowitsch   void foo(void)
6699566063dSJacob Faibussowitsch   {
6709566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type
6719566063dSJacob Faibussowitsch   }
6729566063dSJacob Faibussowitsch 
6739566063dSJacob Faibussowitsch   double bar(void)
6749566063dSJacob Faibussowitsch   {
6759566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type
6769566063dSJacob Faibussowitsch   }
6779566063dSJacob Faibussowitsch 
6789566063dSJacob Faibussowitsch   PetscCallAbort(MPI_COMM_NULL,boom()); // ERROR, communicator should be valid
6799566063dSJacob Faibussowitsch 
6809566063dSJacob Faibussowitsch   struct baz
6819566063dSJacob Faibussowitsch   {
6829566063dSJacob Faibussowitsch     baz()
6839566063dSJacob Faibussowitsch     {
6849566063dSJacob Faibussowitsch       PetscCallAbort(PETSC_COMM_SELF,boom()); // OK
6859566063dSJacob Faibussowitsch     }
6869566063dSJacob Faibussowitsch 
6879566063dSJacob Faibussowitsch     ~baz()
6889566063dSJacob Faibussowitsch     {
6899566063dSJacob Faibussowitsch       PetscCallAbort(PETSC_COMM_SELF,boom()); // OK (in fact the only way to handle PETSc errors)
6909566063dSJacob Faibussowitsch     }
6919566063dSJacob Faibussowitsch   };
6929566063dSJacob Faibussowitsch .ve
6939566063dSJacob Faibussowitsch 
694989c49a2SBarry Smith   Fortran Note:
695989c49a2SBarry Smith   Use `PetscCallA()`.
696989c49a2SBarry Smith 
697989c49a2SBarry Smith   Developer Note:
698989c49a2SBarry Smith   This should have the same name in Fortran as in C.
699989c49a2SBarry Smith 
700db781477SPatrick Sanan .seealso: `SETERRABORT()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`,
7015687f061SJacob Faibussowitsch           `SETERRQ()`, `CHKMEMQ`, `PetscCallMPI()`, `PetscCallCXXAbort()`
7029566063dSJacob Faibussowitsch M*/
7039566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
7049566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm, PetscErrorCode);
7059566063dSJacob Faibussowitsch void PetscCallContinue(PetscErrorCode);
7069566063dSJacob Faibussowitsch #else
7079371c9d4SSatish Balay   #define PetscCallAbort(comm, ...) \
7089371c9d4SSatish Balay     do { \
7097da6b3ccSLisandro Dalcin       PetscErrorCode ierr_petsc_call_abort_; \
7103ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
7117da6b3ccSLisandro Dalcin       ierr_petsc_call_abort_ = __VA_ARGS__; \
7123ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_abort_ != PETSC_SUCCESS)) { \
7133ba16761SJacob Faibussowitsch         ierr_petsc_call_abort_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_abort_, PETSC_ERROR_REPEAT, " "); \
7143ba16761SJacob Faibussowitsch         (void)MPI_Abort(comm, (PetscMPIInt)ierr_petsc_call_abort_); \
7159566063dSJacob Faibussowitsch       } \
7169566063dSJacob Faibussowitsch     } while (0)
7179371c9d4SSatish Balay   #define PetscCallContinue(...) \
7189371c9d4SSatish Balay     do { \
7197da6b3ccSLisandro Dalcin       PetscErrorCode ierr_petsc_call_continue_; \
7203ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
7217da6b3ccSLisandro Dalcin       ierr_petsc_call_continue_ = __VA_ARGS__; \
7223ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_continue_ != PETSC_SUCCESS)) { \
7233ba16761SJacob Faibussowitsch         ierr_petsc_call_continue_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_continue_, PETSC_ERROR_REPEAT, " "); \
7243ba16761SJacob Faibussowitsch         (void)ierr_petsc_call_continue_; \
7253ba16761SJacob Faibussowitsch       } \
7269566063dSJacob Faibussowitsch     } while (0)
7279566063dSJacob Faibussowitsch #endif
7289566063dSJacob Faibussowitsch 
7299566063dSJacob Faibussowitsch /*MC
7309566063dSJacob Faibussowitsch   CHKERRABORT - Checks error code returned from PETSc function. If non-zero it aborts immediately.
7319566063dSJacob Faibussowitsch 
7329566063dSJacob Faibussowitsch   Synopsis:
7339566063dSJacob Faibussowitsch   #include <petscerror.h>
7349566063dSJacob Faibussowitsch   void CHKERRABORT(MPI_Comm comm, PetscErrorCode ierr)
7359566063dSJacob Faibussowitsch 
7369566063dSJacob Faibussowitsch   Not Collective
7379566063dSJacob Faibussowitsch 
7389566063dSJacob Faibussowitsch   Input Parameters:
7399566063dSJacob Faibussowitsch + comm - the MPI communicator
7409566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
7419566063dSJacob Faibussowitsch 
7429566063dSJacob Faibussowitsch   Level: deprecated
7439566063dSJacob Faibussowitsch 
744667f096bSBarry Smith   Note:
745667f096bSBarry Smith   Deprecated in favor of `PetscCallAbort()`. This routine behaves identically to it.
746667f096bSBarry Smith 
747db781477SPatrick Sanan .seealso: `PetscCallAbort()`
7489566063dSJacob Faibussowitsch M*/
7499566063dSJacob Faibussowitsch #define CHKERRABORT(comm, ...) PetscCallAbort(comm, __VA_ARGS__)
7509566063dSJacob Faibussowitsch #define CHKERRCONTINUE(...)    PetscCallContinue(__VA_ARGS__)
7519566063dSJacob Faibussowitsch 
7529566063dSJacob Faibussowitsch /*MC
75387497f52SBarry Smith    CHKERRA - Fortran-only replacement for use of `CHKERRQ()` in the main program, which aborts immediately
754f388eb8bSPatrick Sanan 
755f388eb8bSPatrick Sanan    Synopsis:
756f388eb8bSPatrick Sanan    #include <petscsys.h>
757f388eb8bSPatrick Sanan    PetscErrorCode CHKERRA(PetscErrorCode ierr)
758f388eb8bSPatrick Sanan 
759f388eb8bSPatrick Sanan    Not Collective
760f388eb8bSPatrick Sanan 
7612fe279fdSBarry Smith    Input Parameter:
762f388eb8bSPatrick Sanan .  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
763f388eb8bSPatrick Sanan 
76487497f52SBarry Smith   Level: deprecated
765f388eb8bSPatrick Sanan 
76687497f52SBarry Smith    Note:
76787497f52SBarry Smith    This macro is rarely needed, normal usage is `PetscCallA()` in the main Fortran program.
768f388eb8bSPatrick Sanan 
769989c49a2SBarry Smith    Developer Note:
770989c49a2SBarry Smith    Why isn't this named `CHKERRABORT()` in Fortran?
771989c49a2SBarry Smith 
77287497f52SBarry Smith .seealso: `PetscCall()`, `PetscCallA()`, `PetscCallAbort()`, `CHKERRQ()`, `SETERRA()`, `SETERRQ()`, `SETERRABORT()`
773f388eb8bSPatrick Sanan M*/
774f388eb8bSPatrick Sanan 
7751e4f893aSSatish Balay PETSC_EXTERN PetscBool petscwaitonerrorflg;
7761e4f893aSSatish Balay PETSC_EXTERN PetscBool petscindebugger;
7777c66cc67SJunchao Zhang 
7787c66cc67SJunchao Zhang /*MC
779667f096bSBarry Smith    PETSCABORT - Call `MPI_Abort()` with an informative error code
7807c66cc67SJunchao Zhang 
7817c66cc67SJunchao Zhang    Synopsis:
7827c66cc67SJunchao Zhang    #include <petscsys.h>
7837c66cc67SJunchao Zhang    PETSCABORT(MPI_Comm comm, PetscErrorCode ierr)
7847c66cc67SJunchao Zhang 
7857cdbe19fSJose E. Roman    Collective; No Fortran Support
7867c66cc67SJunchao Zhang 
7877c66cc67SJunchao Zhang    Input Parameters:
7887c66cc67SJunchao Zhang +  comm - A communicator, so that the error can be collective
7897c66cc67SJunchao Zhang -  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
7907c66cc67SJunchao Zhang 
791bf4d2887SBarry Smith    Level: advanced
7927c66cc67SJunchao Zhang 
793bf4d2887SBarry Smith    Notes:
794989c49a2SBarry Smith    If the option `-start_in_debugger` was used then this calls `abort()` to stop the program in the debugger.
795bf4d2887SBarry Smith 
796989c49a2SBarry Smith    if `PetscCIEnabledPortableErrorOutput` is set, which means the code is running in the PETSc test harness (make test),
797989c49a2SBarry Smith    and `comm` is `MPI_COMM_WORLD` it strives to exit cleanly without calling `MPI_Abort()` and instead calling `MPI_Finalize()`.
798660278c0SBarry Smith 
799989c49a2SBarry Smith    This is currently only used when an error propagates up to the C `main()` program and is detected by a `PetscCall()`, `PetscCallMPI()`,
800989c49a2SBarry Smith    or is set in `main()` with `SETERRQ()`. Abort calls such as `SETERRABORT()`,
801989c49a2SBarry Smith    `PetscCheckAbort()`, `PetscCallMPIAbort()`, and `PetscCallAbort()` always call `MPI_Abort()` and do not have any special
802989c49a2SBarry Smith    handling for the test harness.
803989c49a2SBarry Smith 
804989c49a2SBarry Smith    Developer Note:
805989c49a2SBarry Smith    Should the other abort calls also pass through this call instead of calling `MPI_Abort()` directly?
806989c49a2SBarry Smith 
807989c49a2SBarry Smith .seealso: `PetscError()`, `PetscCall()`, `SETERRABORT()`, `PetscCheckAbort()`, `PetscCallMPIAbort()`, `PetscCall()`, `PetscCallMPI()`,
808989c49a2SBarry Smith           `PetscCallAbort()`, `MPI_Abort()`
8097c66cc67SJunchao Zhang M*/
81029f88feaSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
81129f88feaSJacob Faibussowitsch void PETSCABORT(MPI_Comm, PetscErrorCode);
81229f88feaSJacob Faibussowitsch #else
8139371c9d4SSatish Balay   #define PETSCABORT(comm, ...) \
8149371c9d4SSatish Balay     do { \
8153ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_abort_; \
8163ba16761SJacob Faibussowitsch       if (petscwaitonerrorflg) { ierr_petsc_abort_ = PetscSleep(1000); } \
8173ba16761SJacob Faibussowitsch       if (petscindebugger) { \
8183ba16761SJacob Faibussowitsch         abort(); \
8193ba16761SJacob Faibussowitsch       } else { \
8207da6b3ccSLisandro Dalcin         PetscMPIInt size_; \
8213ba16761SJacob Faibussowitsch         ierr_petsc_abort_ = __VA_ARGS__; \
8227da6b3ccSLisandro Dalcin         MPI_Comm_size(comm, &size_); \
8237da6b3ccSLisandro Dalcin         if (PetscCIEnabledPortableErrorOutput && size_ == PetscGlobalSize && ierr_petsc_abort_ != PETSC_ERR_SIG) { \
8249371c9d4SSatish Balay           MPI_Finalize(); \
8259371c9d4SSatish Balay           exit(0); \
826660278c0SBarry Smith         } else if (PetscCIEnabledPortableErrorOutput && PetscGlobalSize == 1) { \
827660278c0SBarry Smith           exit(0); \
828660278c0SBarry Smith         } else { \
829660278c0SBarry Smith           MPI_Abort(comm, (PetscMPIInt)ierr_petsc_abort_); \
830660278c0SBarry Smith         } \
8313fcd9f07SJacob Faibussowitsch       } \
8327c66cc67SJunchao Zhang     } while (0)
83329f88feaSJacob Faibussowitsch #endif
834986eef2eSBarry Smith 
8359566063dSJacob Faibussowitsch #ifdef PETSC_CLANGUAGE_CXX
836986eef2eSBarry Smith   /*MC
8379566063dSJacob Faibussowitsch   PetscCallThrow - Checks error code, if non-zero it calls the C++ error handler which throws
8389566063dSJacob Faibussowitsch   an exception
839986eef2eSBarry Smith 
840986eef2eSBarry Smith   Synopsis:
8419566063dSJacob Faibussowitsch   #include <petscerror.h>
8429566063dSJacob Faibussowitsch   void PetscCallThrow(PetscErrorCode ierr)
843986eef2eSBarry Smith 
844986eef2eSBarry Smith   Not Collective
845986eef2eSBarry Smith 
8469566063dSJacob Faibussowitsch   Input Parameter:
847986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
848986eef2eSBarry Smith 
849667f096bSBarry Smith   Level: beginner
850667f096bSBarry Smith 
851986eef2eSBarry Smith   Notes:
8529566063dSJacob Faibussowitsch   Requires PETSc to be configured with clanguage = c++. Throws a std::runtime_error() on error.
853986eef2eSBarry Smith 
85487497f52SBarry Smith   Once the error handler throws the exception you can use `PetscCallVoid()` which returns without
85587497f52SBarry Smith   an error code (bad idea since the error is ignored) or `PetscCallAbort()` to have `MPI_Abort()`
8569566063dSJacob Faibussowitsch   called immediately.
8579566063dSJacob Faibussowitsch 
858db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscCall()`, `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`,
859db781477SPatrick Sanan           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`
860986eef2eSBarry Smith M*/
8619371c9d4SSatish Balay   #define PetscCallThrow(...) \
8629371c9d4SSatish Balay     do { \
8633ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
8643ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_throw_ = __VA_ARGS__; \
8653ba16761SJacob 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); \
866ffc4695bSBarry Smith     } while (0)
86785614651SBarry Smith 
868cc26af49SMatthew Knepley   /*MC
869cc26af49SMatthew Knepley   CHKERRXX - Checks error code, if non-zero it calls the C++ error handler which throws an exception
870cc26af49SMatthew Knepley 
871cc26af49SMatthew Knepley   Synopsis:
8729566063dSJacob Faibussowitsch   #include <petscerror.h>
8733af045c5SBarry Smith   void CHKERRXX(PetscErrorCode ierr)
874cc26af49SMatthew Knepley 
875eca87e8dSBarry Smith   Not Collective
876cc26af49SMatthew Knepley 
8779566063dSJacob Faibussowitsch   Input Parameter:
8783af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
879cc26af49SMatthew Knepley 
8809566063dSJacob Faibussowitsch   Level: deprecated
881cc26af49SMatthew Knepley 
882667f096bSBarry Smith   Note:
883667f096bSBarry Smith   Deprecated in favor of `PetscCallThrow()`. This routine behaves identically to it.
884667f096bSBarry Smith 
885db781477SPatrick Sanan .seealso: `PetscCallThrow()`
886cc26af49SMatthew Knepley M*/
8879566063dSJacob Faibussowitsch   #define CHKERRXX(...) PetscCallThrow(__VA_ARGS__)
888fd705b32SMatthew Knepley #endif
889fd705b32SMatthew Knepley 
8903ba16761SJacob Faibussowitsch #define PetscCallCXX_Private(__SETERR_FUNC__, __COMM__, ...) \
8913ba16761SJacob Faibussowitsch   do { \
8923ba16761SJacob Faibussowitsch     PetscStackUpdateLine; \
8933ba16761SJacob Faibussowitsch     try { \
8943ba16761SJacob Faibussowitsch       __VA_ARGS__; \
8953ba16761SJacob Faibussowitsch     } catch (const std::exception &e) { \
8963ba16761SJacob Faibussowitsch       __SETERR_FUNC__(__COMM__, PETSC_ERR_LIB, "%s", e.what()); \
8973ba16761SJacob Faibussowitsch     } \
8983ba16761SJacob Faibussowitsch   } while (0)
8993ba16761SJacob Faibussowitsch 
9003f520e80SJunchao Zhang /*MC
9019566063dSJacob Faibussowitsch   PetscCallCXX - Checks C++ function calls and if they throw an exception, catch it and then
9029566063dSJacob Faibussowitsch   return a PETSc error code
9033f520e80SJunchao Zhang 
9043f520e80SJunchao Zhang   Synopsis:
9059566063dSJacob Faibussowitsch   #include <petscerror.h>
9065687f061SJacob Faibussowitsch   void PetscCallCXX(...) noexcept;
9073f520e80SJunchao Zhang 
9083f520e80SJunchao Zhang   Not Collective
9093f520e80SJunchao Zhang 
9109566063dSJacob Faibussowitsch   Input Parameter:
9115687f061SJacob Faibussowitsch . __VA_ARGS__ - An arbitrary expression
9125687f061SJacob Faibussowitsch 
9135687f061SJacob Faibussowitsch   Level: beginner
9143f520e80SJunchao Zhang 
9153f520e80SJunchao Zhang   Notes:
9165687f061SJacob Faibussowitsch   `PetscCallCXX(...)` is a macro replacement for
9179566063dSJacob Faibussowitsch .vb
9189566063dSJacob Faibussowitsch   try {
9195687f061SJacob Faibussowitsch     __VA_ARGS__;
9209566063dSJacob Faibussowitsch   } catch (const std::exception& e) {
9219566063dSJacob Faibussowitsch     return ConvertToPetscErrorCode(e);
9229566063dSJacob Faibussowitsch   }
9239566063dSJacob Faibussowitsch .ve
9249566063dSJacob Faibussowitsch   Due to the fact that it catches any (reasonable) exception, it is essentially noexcept.
9253f520e80SJunchao Zhang 
9265687f061SJacob Faibussowitsch   If you cannot return a `PetscErrorCode` use `PetscCallCXXAbort()` instead.
9275687f061SJacob Faibussowitsch 
9289566063dSJacob Faibussowitsch   Example Usage:
9299566063dSJacob Faibussowitsch .vb
9309566063dSJacob Faibussowitsch   void foo(void) { throw std::runtime_error("error"); }
9313f520e80SJunchao Zhang 
9329566063dSJacob Faibussowitsch   void bar()
9339566063dSJacob Faibussowitsch   {
934e8952933SJacob Faibussowitsch     PetscCallCXX(foo()); // ERROR bar() does not return PetscErrorCode
9359566063dSJacob Faibussowitsch   }
9369566063dSJacob Faibussowitsch 
9379566063dSJacob Faibussowitsch   PetscErrorCode baz()
9389566063dSJacob Faibussowitsch   {
9399566063dSJacob Faibussowitsch     PetscCallCXX(foo()); // OK
9409566063dSJacob Faibussowitsch 
9419566063dSJacob Faibussowitsch     PetscCallCXX(
9429566063dSJacob Faibussowitsch       bar();
943d5b43468SJose E. Roman       foo(); // OK multiple statements allowed
9449566063dSJacob Faibussowitsch     );
9459566063dSJacob Faibussowitsch   }
9469566063dSJacob Faibussowitsch 
9479566063dSJacob Faibussowitsch   struct bop
9489566063dSJacob Faibussowitsch   {
9499566063dSJacob Faibussowitsch     bop()
9509566063dSJacob Faibussowitsch     {
951e8952933SJacob Faibussowitsch       PetscCallCXX(foo()); // ERROR returns PetscErrorCode, cannot be used in constructors
9529566063dSJacob Faibussowitsch     }
9539566063dSJacob Faibussowitsch   };
9549566063dSJacob Faibussowitsch 
955e8952933SJacob Faibussowitsch   // ERROR contains do-while, cannot be used as function-try block
9569566063dSJacob Faibussowitsch   PetscErrorCode qux() PetscCallCXX(
9579566063dSJacob Faibussowitsch     bar();
9589566063dSJacob Faibussowitsch     baz();
9599566063dSJacob Faibussowitsch     foo();
9609566063dSJacob Faibussowitsch     return 0;
9619566063dSJacob Faibussowitsch   )
9629566063dSJacob Faibussowitsch .ve
9639566063dSJacob Faibussowitsch 
9645687f061SJacob Faibussowitsch .seealso: `PetscCallCXXAbort()`, `PetscCallThrow()`, `SETERRQ()`, `PetscCall()`,
9655687f061SJacob Faibussowitsch           `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
9665687f061SJacob Faibussowitsch           `PetscError()`, `CHKMEMQ`
9673f520e80SJunchao Zhang M*/
9683ba16761SJacob Faibussowitsch #define PetscCallCXX(...) PetscCallCXX_Private(SETERRQ, PETSC_COMM_SELF, __VA_ARGS__)
9695687f061SJacob Faibussowitsch 
9705687f061SJacob Faibussowitsch /*MC
9715687f061SJacob Faibussowitsch   PetscCallCXXAbort - Like `PetscCallCXX()` but calls `MPI_Abort()` instead of returning an
9725687f061SJacob Faibussowitsch   error-code
9735687f061SJacob Faibussowitsch 
9745687f061SJacob Faibussowitsch   Synopsis:
9755687f061SJacob Faibussowitsch   #include <petscerror.h>
9765687f061SJacob Faibussowitsch   void PetscCallCXXAbort(MPI_Comm comm, ...) noexcept;
9775687f061SJacob Faibussowitsch 
97820f4b53cSBarry Smith   Collective; No Fortran Support
9795687f061SJacob Faibussowitsch 
9805687f061SJacob Faibussowitsch   Input Parameters:
9815687f061SJacob Faibussowitsch + comm        - The MPI communicator to abort on
9825687f061SJacob Faibussowitsch - __VA_ARGS__ - An arbitrary expression
9835687f061SJacob Faibussowitsch 
9845687f061SJacob Faibussowitsch   Level: beginner
9855687f061SJacob Faibussowitsch 
9865687f061SJacob Faibussowitsch   Notes:
9875687f061SJacob Faibussowitsch   This macro may be used to check C++ expressions for exceptions in cases where you cannot
9885687f061SJacob Faibussowitsch   return an error code. This includes constructors, destructors, copy/move assignment functions
9895687f061SJacob Faibussowitsch   or constructors among others.
9905687f061SJacob Faibussowitsch 
9915687f061SJacob Faibussowitsch   If an exception is caught, the macro calls `SETERRABORT()` on `comm`. The exception must
9925687f061SJacob Faibussowitsch   derive from `std::exception` in order to be caught.
9935687f061SJacob Faibussowitsch 
9945687f061SJacob Faibussowitsch   If the routine _can_ return an error-code it is highly advised to use `PetscCallCXX()`
9955687f061SJacob Faibussowitsch   instead.
9965687f061SJacob Faibussowitsch 
9975687f061SJacob Faibussowitsch   See `PetscCallCXX()` for additional discussion.
9985687f061SJacob Faibussowitsch 
9995687f061SJacob Faibussowitsch   Example Usage:
10005687f061SJacob Faibussowitsch .vb
10015687f061SJacob Faibussowitsch   class Foo
10025687f061SJacob Faibussowitsch   {
10035687f061SJacob Faibussowitsch     std::vector<int> data_;
10045687f061SJacob Faibussowitsch 
10055687f061SJacob Faibussowitsch   public:
10065687f061SJacob Faibussowitsch     // normally std::vector::reserve() may raise an exception, but since we handle it with
10075687f061SJacob Faibussowitsch     // PetscCallCXXAbort() we may mark this routine as noexcept!
10085687f061SJacob Faibussowitsch     Foo() noexcept
10095687f061SJacob Faibussowitsch     {
10105687f061SJacob Faibussowitsch       PetscCallCXXAbort(PETSC_COMM_SELF, data_.reserve(10));
10115687f061SJacob Faibussowitsch     }
10125687f061SJacob Faibussowitsch   };
10135687f061SJacob Faibussowitsch 
10145687f061SJacob Faibussowitsch   std::vector<int> bar()
10155687f061SJacob Faibussowitsch   {
10165687f061SJacob Faibussowitsch     std::vector<int> v;
10175687f061SJacob Faibussowitsch 
10185687f061SJacob Faibussowitsch     PetscFunctionBegin;
10195687f061SJacob Faibussowitsch     // OK!
10205687f061SJacob Faibussowitsch     PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1));
10215687f061SJacob Faibussowitsch     PetscFunctionReturn(v);
10225687f061SJacob Faibussowitsch   }
10235687f061SJacob Faibussowitsch 
10245687f061SJacob Faibussowitsch   PetscErrorCode baz()
10255687f061SJacob Faibussowitsch   {
10265687f061SJacob Faibussowitsch     std::vector<int> v;
10275687f061SJacob Faibussowitsch 
10285687f061SJacob Faibussowitsch     PetscFunctionBegin;
10295687f061SJacob Faibussowitsch     // WRONG! baz() returns a PetscErrorCode, prefer PetscCallCXX() instead
10305687f061SJacob Faibussowitsch     PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1));
10313ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
10325687f061SJacob Faibussowitsch   }
10335687f061SJacob Faibussowitsch .ve
10345687f061SJacob Faibussowitsch 
10355687f061SJacob Faibussowitsch .seealso: `PetscCallCXX()`, `SETERRABORT()`, `PetscCallAbort()`
10365687f061SJacob Faibussowitsch M*/
10373ba16761SJacob Faibussowitsch #define PetscCallCXXAbort(comm, ...) PetscCallCXX_Private(SETERRABORT, comm, __VA_ARGS__)
10383f520e80SJunchao Zhang 
103930de9b25SBarry Smith /*MC
10409566063dSJacob Faibussowitsch   CHKERRCXX - Checks C++ function calls and if they throw an exception, catch it and then
10419566063dSJacob Faibussowitsch   return a PETSc error code
10429566063dSJacob Faibussowitsch 
10439566063dSJacob Faibussowitsch   Synopsis:
10449566063dSJacob Faibussowitsch   #include <petscerror.h>
10459566063dSJacob Faibussowitsch   void CHKERRCXX(func) noexcept;
10469566063dSJacob Faibussowitsch 
10479566063dSJacob Faibussowitsch   Not Collective
10489566063dSJacob Faibussowitsch 
10499566063dSJacob Faibussowitsch   Input Parameter:
10509566063dSJacob Faibussowitsch . func - C++ function calls
10519566063dSJacob Faibussowitsch 
10529566063dSJacob Faibussowitsch   Level: deprecated
10539566063dSJacob Faibussowitsch 
1054667f096bSBarry Smith   Note:
1055667f096bSBarry Smith   Deprecated in favor of `PetscCallCXX()`. This routine behaves identically to it.
1056667f096bSBarry Smith 
1057db781477SPatrick Sanan .seealso: `PetscCallCXX()`
10589566063dSJacob Faibussowitsch M*/
10599566063dSJacob Faibussowitsch #define CHKERRCXX(...) PetscCallCXX(__VA_ARGS__)
10609566063dSJacob Faibussowitsch 
10619566063dSJacob Faibussowitsch /*MC
106230de9b25SBarry Smith    CHKMEMQ - Checks the memory for corruption, calls error handler if any is detected
106330de9b25SBarry Smith 
106430de9b25SBarry Smith    Synopsis:
1065aaa7dc30SBarry Smith    #include <petscsys.h>
106691d3bdf4SKris Buschelman    CHKMEMQ;
106730de9b25SBarry Smith 
1068eca87e8dSBarry Smith    Not Collective
1069eca87e8dSBarry Smith 
107030de9b25SBarry Smith   Level: beginner
107130de9b25SBarry Smith 
107230de9b25SBarry Smith    Notes:
1073a17b96a8SKyle Gerard Felker     We highly recommend using Valgrind https://petsc.org/release/faq/#valgrind or for NVIDIA CUDA systems
10745ed36255SBarry Smith     https://docs.nvidia.com/cuda/cuda-memcheck/index.html for finding memory problems. The ``CHKMEMQ`` macro is useful on systems that
10755ed36255SBarry Smith     do not have valgrind, but is not as good as valgrind or cuda-memcheck.
10761957e957SBarry Smith 
1077667f096bSBarry Smith     Must run with the option `-malloc_debug` (`-malloc_test` in debug mode; or if `PetscMallocSetDebug()` called) to enable this option
107830de9b25SBarry Smith 
107930de9b25SBarry Smith     Once the error handler is called the calling function is then returned from with the given error code.
108030de9b25SBarry Smith 
108130de9b25SBarry Smith     By defaults prints location where memory that is corrupted was allocated.
108230de9b25SBarry Smith 
108387497f52SBarry Smith     Use `CHKMEMA` for functions that return void
1084f621e05eSBarry Smith 
1085db781477SPatrick Sanan .seealso: `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `SETERRQ()`, `PetscMallocValidate()`
108630de9b25SBarry Smith M*/
10876d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
1088064a246eSJacob Faibussowitsch   #define CHKMEMQ
1089064a246eSJacob Faibussowitsch   #define CHKMEMA
10906d210af2SJacob Faibussowitsch #else
10919371c9d4SSatish Balay   #define CHKMEMQ \
10929371c9d4SSatish Balay     do { \
10933ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_memq_ = PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__); \
10943ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_memq_ != PETSC_SUCCESS)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_memq_, PETSC_ERROR_REPEAT, " "); \
109586d09637SLisandro Dalcin     } while (0)
10966d210af2SJacob Faibussowitsch   #define CHKMEMA PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__)
1097064a246eSJacob Faibussowitsch #endif
10989566063dSJacob Faibussowitsch 
1099668f157eSBarry Smith /*E
1100668f157eSBarry Smith   PetscErrorType - passed to the PETSc error handling routines indicating if this is the first or a later call to the error handlers
1101668f157eSBarry Smith 
1102668f157eSBarry Smith   Level: advanced
1103668f157eSBarry Smith 
1104667f096bSBarry Smith   Note:
110587497f52SBarry Smith   `PETSC_ERROR_IN_CXX` indicates the error was detected in C++ and an exception should be generated
1106d736bfebSBarry Smith 
1107667f096bSBarry Smith   Developer Note:
1108667f096bSBarry Smith     This is currently used to decide when to print the detailed information about the run in `PetscTraceBackErrorHandler()`
1109668f157eSBarry Smith 
111087497f52SBarry Smith .seealso: `PetscError()`, `SETERRQ()`
1111668f157eSBarry Smith E*/
11129371c9d4SSatish Balay typedef enum {
11139371c9d4SSatish Balay   PETSC_ERROR_INITIAL = 0,
11149371c9d4SSatish Balay   PETSC_ERROR_REPEAT  = 1,
11159371c9d4SSatish Balay   PETSC_ERROR_IN_CXX  = 2
11169371c9d4SSatish Balay } PetscErrorType;
11174b209cf6SBarry Smith 
1118eb9e708aSLisandro Dalcin #if defined(__clang_analyzer__)
1119eb9e708aSLisandro Dalcin __attribute__((analyzer_noreturn))
1120eb9e708aSLisandro Dalcin #endif
11219371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode
11229371c9d4SSatish Balay PetscError(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, ...) PETSC_ATTRIBUTE_COLD PETSC_ATTRIBUTE_FORMAT(7, 8);
1123eb9e708aSLisandro Dalcin 
1124014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscErrorPrintfInitialize(void);
11253ba16761SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscErrorMessage(PetscErrorCode, const char *[], char **);
1126d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscTraceBackErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1127d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscIgnoreErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1128d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscEmacsClientErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1129d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscMPIAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1130d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1131d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAttachDebuggerErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1132d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscReturnErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1133efca3c55SSatish Balay PETSC_EXTERN PetscErrorCode PetscPushErrorHandler(PetscErrorCode (*handler)(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *), void *);
1134014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopErrorHandler(void);
11358d359177SBarry Smith PETSC_EXTERN PetscErrorCode PetscSignalHandlerDefault(int, void *);
1136014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPushSignalHandler(PetscErrorCode (*)(int, void *), void *);
1137014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopSignalHandler(void);
113828559dc8SJed Brown PETSC_EXTERN PetscErrorCode PetscCheckPointerSetIntensity(PetscInt);
1139c2a741eeSJunchao Zhang PETSC_EXTERN void           PetscSignalSegvCheckPointerOrMpi(void);
1140*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 13, 0, "PetscSignalSegvCheckPointerOrMpi()", ) static inline void PetscSignalSegvCheckPointer(void)
1141d71ae5a4SJacob Faibussowitsch {
11429371c9d4SSatish Balay   PetscSignalSegvCheckPointerOrMpi();
11439371c9d4SSatish Balay }
1144329f5518SBarry Smith 
1145639ff905SBarry Smith /*MC
1146639ff905SBarry Smith     PetscErrorPrintf - Prints error messages.
1147639ff905SBarry Smith 
1148639ff905SBarry Smith    Synopsis:
1149aaa7dc30SBarry Smith     #include <petscsys.h>
1150639ff905SBarry Smith      PetscErrorCode (*PetscErrorPrintf)(const char format[],...);
1151639ff905SBarry Smith 
11527cdbe19fSJose E. Roman     Not Collective; No Fortran Support
11537cdbe19fSJose E. Roman 
1154f899ff85SJose E. Roman     Input Parameter:
1155cd05f99aSJacob Faibussowitsch .   format - the usual `printf()` format string
1156639ff905SBarry Smith 
1157639ff905SBarry Smith    Options Database Keys:
11581957e957SBarry Smith +    -error_output_stdout - cause error messages to be printed to stdout instead of the (default) stderr
1159e1bc860dSBarry Smith -    -error_output_none - to turn off all printing of error messages (does not change the way the error is handled.)
1160639ff905SBarry Smith 
1161cf53795eSBarry Smith    Level: developer
1162cf53795eSBarry Smith 
116395452b02SPatrick Sanan    Notes:
116495452b02SPatrick Sanan     Use
1165667f096bSBarry Smith .vb
1166667f096bSBarry Smith      PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the
1167667f096bSBarry Smith                         error is handled.) and
1168667f096bSBarry Smith      PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on or you can use your own function
1169667f096bSBarry Smith .ve
1170639ff905SBarry Smith      Use
1171667f096bSBarry Smith .vb
117287497f52SBarry Smith      `PETSC_STDERR` = FILE* obtained from a file open etc. to have stderr printed to the file.
117387497f52SBarry Smith      `PETSC_STDOUT` = FILE* obtained from a file open etc. to have stdout printed to the file.
1174667f096bSBarry Smith .ve
1175639ff905SBarry Smith 
1176639ff905SBarry Smith        Use
117787497f52SBarry Smith       `PetscPushErrorHandler()` to provide your own error handler that determines what kind of messages to print
1178639ff905SBarry Smith 
1179db781477SPatrick Sanan .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscHelpPrintf()`, `PetscPrintf()`, `PetscPushErrorHandler()`, `PetscVFPrintf()`, `PetscHelpPrintf()`
1180639ff905SBarry Smith M*/
11813ca90d2dSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode (*PetscErrorPrintf)(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1182639ff905SBarry Smith 
1183cf0818bdSBarry Smith /*E
1184cf0818bdSBarry Smith      PetscFPTrap - types of floating point exceptions that may be trapped
1185cf0818bdSBarry Smith 
118687497f52SBarry Smith      Currently only `PETSC_FP_TRAP_OFF` and `PETSC_FP_TRAP_ON` are handled. All others are treated as `PETSC_FP_TRAP_ON`.
1187cf0818bdSBarry Smith 
1188cf0818bdSBarry Smith      Level: intermediate
1189cf0818bdSBarry Smith 
11907de69702SBarry Smith .seealso: `PetscSetFPTrap()`, `PetscFPTrapPush()`
1191cf0818bdSBarry Smith  E*/
11929371c9d4SSatish Balay typedef enum {
11939371c9d4SSatish Balay   PETSC_FP_TRAP_OFF      = 0,
11949371c9d4SSatish Balay   PETSC_FP_TRAP_INDIV    = 1,
11959371c9d4SSatish Balay   PETSC_FP_TRAP_FLTOPERR = 2,
11969371c9d4SSatish Balay   PETSC_FP_TRAP_FLTOVF   = 4,
11979371c9d4SSatish Balay   PETSC_FP_TRAP_FLTUND   = 8,
11989371c9d4SSatish Balay   PETSC_FP_TRAP_FLTDIV   = 16,
11999371c9d4SSatish Balay   PETSC_FP_TRAP_FLTINEX  = 32
12009371c9d4SSatish Balay } PetscFPTrap;
1201bd2b07b1SBarry 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)
1202014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSetFPTrap(PetscFPTrap);
1203014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPush(PetscFPTrap);
1204014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPop(void);
1205aba4c478SBarry Smith PETSC_EXTERN PetscErrorCode PetscDetermineInitialFPTrap(void);
120654a8ef01SBarry Smith 
12073a40ed3dSBarry Smith /*
12083a40ed3dSBarry Smith       Allows the code to build a stack frame as it runs
12093a40ed3dSBarry Smith */
12103a40ed3dSBarry Smith 
121199cd645aSJed Brown #define PETSCSTACKSIZE 64
12123a40ed3dSBarry Smith typedef struct {
12130e33f6ddSBarry Smith   const char *function[PETSCSTACKSIZE];
12140e33f6ddSBarry Smith   const char *file[PETSCSTACKSIZE];
1215184914b5SBarry Smith   int         line[PETSCSTACKSIZE];
1216362febeeSStefano Zampini   int         petscroutine[PETSCSTACKSIZE]; /* 0 external called from petsc, 1 petsc functions, 2 petsc user functions */
1217184914b5SBarry Smith   int         currentsize;
1218a2f94806SJed Brown   int         hotdepth;
12194be741a6SBarry Smith   PetscBool   check; /* option to check for correct Push/Pop semantics, true for default petscstack but not other stacks */
12203a40ed3dSBarry Smith } PetscStack;
1221dfb7d7afSStefano Zampini #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
122227104ee2SJacob Faibussowitsch PETSC_EXTERN PetscStack petscstack;
122327104ee2SJacob Faibussowitsch #endif
12243a40ed3dSBarry Smith 
12255d12eec7SSatish Balay #if defined(PETSC_SERIALIZE_FUNCTIONS)
12265d12eec7SSatish Balay   #include <petsc/private/petscfptimpl.h>
12275d12eec7SSatish Balay   /*
12285d12eec7SSatish Balay    Registers the current function into the global function pointer to function name table
12295d12eec7SSatish Balay 
12305d12eec7SSatish Balay    Have to fix this to handle errors but cannot return error since used in PETSC_VIEWER_DRAW_() etc
12315d12eec7SSatish Balay */
12329371c9d4SSatish Balay   #define PetscRegister__FUNCT__() \
12339371c9d4SSatish Balay     do { \
12345d12eec7SSatish Balay       static PetscBool __chked = PETSC_FALSE; \
12355d12eec7SSatish Balay       if (!__chked) { \
12369371c9d4SSatish Balay         void *ptr; \
12373ba16761SJacob Faibussowitsch         PetscCallAbort(PETSC_COMM_SELF, PetscDLSym(NULL, PETSC_FUNCTION_NAME, &ptr)); \
12385d12eec7SSatish Balay         __chked = PETSC_TRUE; \
12399371c9d4SSatish Balay       } \
12409371c9d4SSatish Balay     } while (0)
12415d12eec7SSatish Balay #else
12425d12eec7SSatish Balay   #define PetscRegister__FUNCT__()
12435d12eec7SSatish Balay #endif
12445d12eec7SSatish Balay 
1245eae3dc7dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) || defined(__clang_analyzer__)
12466d210af2SJacob Faibussowitsch   #define PetscStackPushNoCheck(funct, petsc_routine, hot)
124737154d10SBarry Smith   #define PetscStackUpdateLine
1248792fecdfSBarry Smith   #define PetscStackPushExternal(funct)
12496d210af2SJacob Faibussowitsch   #define PetscStackPopNoCheck
12506d210af2SJacob Faibussowitsch   #define PetscStackClearTop
12516d210af2SJacob Faibussowitsch   #define PetscFunctionBegin
12526d210af2SJacob Faibussowitsch   #define PetscFunctionBeginUser
12536d210af2SJacob Faibussowitsch   #define PetscFunctionBeginHot
12540a57284eSJacob Faibussowitsch   #define PetscFunctionReturn(...)  return __VA_ARGS__
12556d210af2SJacob Faibussowitsch   #define PetscFunctionReturnVoid() return
12566d210af2SJacob Faibussowitsch   #define PetscStackPop
12576d210af2SJacob Faibussowitsch   #define PetscStackPush(f)
1258dfb7d7afSStefano Zampini #elif defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1259660278c0SBarry Smith 
12609371c9d4SSatish Balay   #define PetscStackPush_Private(stack__, file__, func__, line__, petsc_routine__, hot__) \
12619371c9d4SSatish Balay     do { \
12625a96b57dSJacob Faibussowitsch       if (stack__.currentsize < PETSCSTACKSIZE) { \
12635a96b57dSJacob Faibussowitsch         stack__.function[stack__.currentsize] = func__; \
1264ef1023bdSBarry Smith         if (petsc_routine__) { \
1265ef1023bdSBarry Smith           stack__.file[stack__.currentsize] = file__; \
12665a96b57dSJacob Faibussowitsch           stack__.line[stack__.currentsize] = line__; \
1267ef1023bdSBarry Smith         } else { \
1268648bc8c4SBarry Smith           stack__.file[stack__.currentsize] = PETSC_NULLPTR; \
1269ef1023bdSBarry Smith           stack__.line[stack__.currentsize] = 0; \
1270ef1023bdSBarry Smith         } \
12715a96b57dSJacob Faibussowitsch         stack__.petscroutine[stack__.currentsize] = petsc_routine__; \
12725a96b57dSJacob Faibussowitsch       } \
12735a96b57dSJacob Faibussowitsch       ++stack__.currentsize; \
12745a96b57dSJacob Faibussowitsch       stack__.hotdepth += (hot__ || stack__.hotdepth); \
12755a96b57dSJacob Faibussowitsch     } while (0)
12765a96b57dSJacob Faibussowitsch 
12774be741a6SBarry Smith   /* uses PetscCheckAbort() because may be used in a function that does not return an error code */
12789371c9d4SSatish Balay   #define PetscStackPop_Private(stack__, func__) \
12799371c9d4SSatish Balay     do { \
12804be741a6SBarry 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__); \
12815a96b57dSJacob Faibussowitsch       if (--stack__.currentsize < PETSCSTACKSIZE) { \
12829371c9d4SSatish 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", \
12839371c9d4SSatish Balay                         stack__.function[stack__.currentsize], stack__.file[stack__.currentsize], stack__.line[stack__.currentsize], func__, __FILE__, __LINE__); \
12845a96b57dSJacob Faibussowitsch         stack__.function[stack__.currentsize]     = PETSC_NULLPTR; \
12855a96b57dSJacob Faibussowitsch         stack__.file[stack__.currentsize]         = PETSC_NULLPTR; \
12865a96b57dSJacob Faibussowitsch         stack__.line[stack__.currentsize]         = 0; \
12875a96b57dSJacob Faibussowitsch         stack__.petscroutine[stack__.currentsize] = 0; \
12885a96b57dSJacob Faibussowitsch       } \
12895a96b57dSJacob Faibussowitsch       stack__.hotdepth = PetscMax(stack__.hotdepth - 1, 0); \
12905a96b57dSJacob Faibussowitsch     } while (0)
12915a96b57dSJacob Faibussowitsch 
1292586f9135SBarry Smith   /*MC
1293586f9135SBarry Smith    PetscStackPushNoCheck - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is
1294586f9135SBarry Smith    currently in the source code.
1295586f9135SBarry Smith 
1296586f9135SBarry Smith    Synopsis:
1297586f9135SBarry Smith    #include <petscsys.h>
1298586f9135SBarry Smith    void PetscStackPushNoCheck(char *funct,int petsc_routine,PetscBool hot);
1299586f9135SBarry Smith 
13007cdbe19fSJose E. Roman    Not Collective
13017cdbe19fSJose E. Roman 
1302586f9135SBarry Smith    Input Parameters:
1303586f9135SBarry Smith +  funct - the function name
1304586f9135SBarry Smith .  petsc_routine - 2 user function, 1 PETSc function, 0 some other function
1305586f9135SBarry Smith -  hot - indicates that the function may be called often so expensive error checking should be turned off inside the function
1306586f9135SBarry Smith 
1307586f9135SBarry Smith    Level: developer
1308586f9135SBarry Smith 
1309586f9135SBarry Smith    Notes:
1310586f9135SBarry 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
131187497f52SBarry 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
1312586f9135SBarry Smith    help debug the problem.
1313586f9135SBarry Smith 
1314ef1023bdSBarry Smith    This version does not check the memory corruption (an expensive operation), use `PetscStackPush()` to check the memory.
1315ef1023bdSBarry Smith 
1316792fecdfSBarry Smith    Use `PetscStackPushExternal()` for a function call that is about to be made to a non-PETSc or user function (such as BLAS etc).
1317ef1023bdSBarry Smith 
1318586f9135SBarry Smith    The default stack is a global variable called `petscstack`.
1319586f9135SBarry Smith 
1320586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1321ef1023bdSBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPush()`, `PetscStackPop`,
1322792fecdfSBarry Smith           `PetscStackPushExternal()`
1323586f9135SBarry Smith M*/
13249371c9d4SSatish Balay   #define PetscStackPushNoCheck(funct, petsc_routine, hot) \
13259371c9d4SSatish Balay     do { \
1326e04113cfSBarry Smith       PetscStackSAWsTakeAccess(); \
13275a96b57dSJacob Faibussowitsch       PetscStackPush_Private(petscstack, __FILE__, funct, __LINE__, petsc_routine, hot); \
1328e04113cfSBarry Smith       PetscStackSAWsGrantAccess(); \
1329441dd030SJed Brown     } while (0)
1330441dd030SJed Brown 
1331586f9135SBarry Smith   /*MC
133287497f52SBarry Smith    PetscStackUpdateLine - in a function that has a `PetscFunctionBegin` or `PetscFunctionBeginUser` updates the stack line number to the
133337154d10SBarry Smith    current line number.
133437154d10SBarry Smith 
133537154d10SBarry Smith    Synopsis:
133637154d10SBarry Smith    #include <petscsys.h>
133737154d10SBarry Smith    void PetscStackUpdateLine
133837154d10SBarry Smith 
13397cdbe19fSJose E. Roman    Not Collective
13407cdbe19fSJose E. Roman 
134137154d10SBarry Smith    Level: developer
134237154d10SBarry Smith 
134337154d10SBarry Smith    Notes:
134487497f52SBarry Smith    Using `PetscCall()` and friends automatically handles this process
134587497f52SBarry Smith 
134637154d10SBarry 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
134737154d10SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
134837154d10SBarry Smith    help debug the problem.
134937154d10SBarry Smith 
135037154d10SBarry Smith    The default stack is a global variable called petscstack.
135137154d10SBarry Smith 
135237154d10SBarry Smith    This is used by `PetscCall()` and is otherwise not like to be needed
135337154d10SBarry Smith 
135437154d10SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`, `PetscCall()`
135537154d10SBarry Smith M*/
135637154d10SBarry Smith   #define PetscStackUpdateLine \
13573ba16761SJacob Faibussowitsch     do { \
13583ba16761SJacob Faibussowitsch       if (petscstack.currentsize > 0 && petscstack.function[petscstack.currentsize - 1] == PETSC_FUNCTION_NAME) { petscstack.line[petscstack.currentsize - 1] = __LINE__; } \
13593ba16761SJacob Faibussowitsch     } while (0)
136037154d10SBarry Smith 
136137154d10SBarry Smith   /*MC
1362792fecdfSBarry Smith    PetscStackPushExternal - Pushes a new function name onto the PETSc default stack that tracks where the running program is
1363ef1023bdSBarry Smith    currently in the source code. Does not include the filename or line number since this is called by the calling routine
1364ef1023bdSBarry Smith    for non-PETSc or user functions.
1365ef1023bdSBarry Smith 
1366ef1023bdSBarry Smith    Synopsis:
1367ef1023bdSBarry Smith    #include <petscsys.h>
1368792fecdfSBarry Smith    void PetscStackPushExternal(char *funct);
1369ef1023bdSBarry Smith 
13707cdbe19fSJose E. Roman    Not Collective
13717cdbe19fSJose E. Roman 
13722fe279fdSBarry Smith    Input Parameter:
1373ef1023bdSBarry Smith .  funct - the function name
1374ef1023bdSBarry Smith 
1375ef1023bdSBarry Smith    Level: developer
1376ef1023bdSBarry Smith 
1377ef1023bdSBarry Smith    Notes:
137887497f52SBarry Smith    Using `PetscCallExternal()` and friends automatically handles this process
137987497f52SBarry Smith 
1380ef1023bdSBarry 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
1381ef1023bdSBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1382ef1023bdSBarry Smith    help debug the problem.
1383ef1023bdSBarry Smith 
1384ef1023bdSBarry Smith    The default stack is a global variable called `petscstack`.
1385ef1023bdSBarry Smith 
1386ef1023bdSBarry Smith    This is to be used when calling an external package function such as a BLAS function.
1387ef1023bdSBarry Smith 
1388ef1023bdSBarry Smith    This also updates the stack line number for the current stack function.
1389ef1023bdSBarry Smith 
1390ef1023bdSBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1391ef1023bdSBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1392ef1023bdSBarry Smith M*/
13939371c9d4SSatish Balay   #define PetscStackPushExternal(funct) \
13949371c9d4SSatish Balay     do { \
13959371c9d4SSatish Balay       PetscStackUpdateLine; \
13969371c9d4SSatish Balay       PetscStackPushNoCheck(funct, 0, PETSC_TRUE); \
13979371c9d4SSatish Balay     } while (0);
1398ef1023bdSBarry Smith 
1399ef1023bdSBarry Smith   /*MC
1400586f9135SBarry Smith    PetscStackPopNoCheck - Pops a function name from the PETSc default stack that tracks where the running program is
1401586f9135SBarry Smith    currently in the source code.
1402586f9135SBarry Smith 
1403586f9135SBarry Smith    Synopsis:
1404586f9135SBarry Smith    #include <petscsys.h>
1405586f9135SBarry Smith    void PetscStackPopNoCheck(char *funct);
1406586f9135SBarry Smith 
14077cdbe19fSJose E. Roman    Not Collective
14087cdbe19fSJose E. Roman 
1409586f9135SBarry Smith    Input Parameter:
1410586f9135SBarry Smith .   funct - the function name
1411586f9135SBarry Smith 
1412586f9135SBarry Smith    Level: developer
1413586f9135SBarry Smith 
1414586f9135SBarry Smith    Notes:
141587497f52SBarry Smith    Using `PetscCall()`, `PetscCallExternal()`, `PetscCallBack()` and friends negates the need to call this
141687497f52SBarry Smith 
1417586f9135SBarry 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
1418586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1419586f9135SBarry Smith    help debug the problem.
1420586f9135SBarry Smith 
1421586f9135SBarry Smith    The default stack is a global variable called petscstack.
1422586f9135SBarry Smith 
1423586f9135SBarry Smith    Developer Note:
1424586f9135SBarry Smith    `PetscStackPopNoCheck()` takes a function argument while  `PetscStackPop` does not, this difference is likely just historical.
1425586f9135SBarry Smith 
1426586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1427586f9135SBarry Smith M*/
14289371c9d4SSatish Balay   #define PetscStackPopNoCheck(funct) \
14299371c9d4SSatish Balay     do { \
1430362febeeSStefano Zampini       PetscStackSAWsTakeAccess(); \
14315a96b57dSJacob Faibussowitsch       PetscStackPop_Private(petscstack, funct); \
1432362febeeSStefano Zampini       PetscStackSAWsGrantAccess(); \
1433362febeeSStefano Zampini     } while (0)
1434362febeeSStefano Zampini 
14359371c9d4SSatish Balay   #define PetscStackClearTop \
14369371c9d4SSatish Balay     do { \
1437e04113cfSBarry Smith       PetscStackSAWsTakeAccess(); \
14389371c9d4SSatish Balay       if (petscstack.currentsize > 0 && --petscstack.currentsize < PETSCSTACKSIZE) { \
143927104ee2SJacob Faibussowitsch         petscstack.function[petscstack.currentsize]     = PETSC_NULLPTR; \
144027104ee2SJacob Faibussowitsch         petscstack.file[petscstack.currentsize]         = PETSC_NULLPTR; \
144127104ee2SJacob Faibussowitsch         petscstack.line[petscstack.currentsize]         = 0; \
144227104ee2SJacob Faibussowitsch         petscstack.petscroutine[petscstack.currentsize] = 0; \
1443441dd030SJed Brown       } \
144427104ee2SJacob Faibussowitsch       petscstack.hotdepth = PetscMax(petscstack.hotdepth - 1, 0); \
1445e04113cfSBarry Smith       PetscStackSAWsGrantAccess(); \
1446441dd030SJed Brown     } while (0)
1447441dd030SJed Brown 
144830de9b25SBarry Smith   /*MC
14491957e957SBarry Smith    PetscFunctionBegin - First executable line of each PETSc function,  used for error handling. Final
145087497f52SBarry Smith       line of PETSc functions should be `PetscFunctionReturn`(0);
145130de9b25SBarry Smith 
145230de9b25SBarry Smith    Synopsis:
1453aaa7dc30SBarry Smith    #include <petscsys.h>
145430de9b25SBarry Smith    void PetscFunctionBegin;
145530de9b25SBarry Smith 
145620f4b53cSBarry Smith    Not Collective; No Fortran Support
1457eca87e8dSBarry Smith 
145830de9b25SBarry Smith    Usage:
145930de9b25SBarry Smith .vb
146030de9b25SBarry Smith      int something;
146130de9b25SBarry Smith 
146230de9b25SBarry Smith      PetscFunctionBegin;
146330de9b25SBarry Smith .ve
146430de9b25SBarry Smith 
146530de9b25SBarry Smith    Level: developer
146630de9b25SBarry Smith 
146720f4b53cSBarry Smith    Note:
146820f4b53cSBarry Smith      Use `PetscFunctionBeginUser` for application codes.
146920f4b53cSBarry Smith 
1470586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`
147130de9b25SBarry Smith 
147230de9b25SBarry Smith M*/
14739371c9d4SSatish Balay   #define PetscFunctionBegin \
14749371c9d4SSatish Balay     do { \
1475362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_FALSE); \
1476a2f94806SJed Brown       PetscRegister__FUNCT__(); \
1477a2f94806SJed Brown     } while (0)
1478a2f94806SJed Brown 
1479a2f94806SJed Brown   /*MC
148087497f52SBarry Smith    PetscFunctionBeginHot - Substitute for `PetscFunctionBegin` to be used in functions that are called in
1481a2f94806SJed Brown    performance-critical circumstances.  Use of this function allows for lighter profiling by default.
1482a2f94806SJed Brown 
1483a2f94806SJed Brown    Synopsis:
1484aaa7dc30SBarry Smith    #include <petscsys.h>
1485a2f94806SJed Brown    void PetscFunctionBeginHot;
1486a2f94806SJed Brown 
148720f4b53cSBarry Smith    Not Collective; No Fortran Support
1488a2f94806SJed Brown 
1489a2f94806SJed Brown    Usage:
1490a2f94806SJed Brown .vb
1491a2f94806SJed Brown      int something;
1492a2f94806SJed Brown 
1493a2f94806SJed Brown      PetscFunctionBeginHot;
1494a2f94806SJed Brown .ve
1495a2f94806SJed Brown 
1496a2f94806SJed Brown    Level: developer
1497a2f94806SJed Brown 
1498586f9135SBarry Smith .seealso: `PetscFunctionBegin`, `PetscFunctionReturn()`, `PetscStackPushNoCheck()`
1499a2f94806SJed Brown 
1500a2f94806SJed Brown M*/
15019371c9d4SSatish Balay   #define PetscFunctionBeginHot \
15029371c9d4SSatish Balay     do { \
1503362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_TRUE); \
15042d53ad75SBarry Smith       PetscRegister__FUNCT__(); \
150553c77d0aSJed Brown     } while (0)
150653c77d0aSJed Brown 
1507a8d2bbe5SBarry Smith   /*MC
1508530d85c1SBarry Smith    PetscFunctionBeginUser - First executable line of user provided routines
1509a8d2bbe5SBarry Smith 
1510a8d2bbe5SBarry Smith    Synopsis:
1511aaa7dc30SBarry Smith    #include <petscsys.h>
1512a8d2bbe5SBarry Smith    void PetscFunctionBeginUser;
1513a8d2bbe5SBarry Smith 
151420f4b53cSBarry Smith    Not Collective; No Fortran Support
1515a8d2bbe5SBarry Smith 
1516a8d2bbe5SBarry Smith    Usage:
1517a8d2bbe5SBarry Smith .vb
1518a8d2bbe5SBarry Smith      int something;
1519a8d2bbe5SBarry Smith 
1520ac285190SBarry Smith      PetscFunctionBeginUser;
1521a8d2bbe5SBarry Smith .ve
1522a8d2bbe5SBarry Smith 
152320f4b53cSBarry Smith    Level: intermediate
152420f4b53cSBarry Smith 
1525a8d2bbe5SBarry Smith    Notes:
1526530d85c1SBarry Smith       Functions that incorporate this must call `PetscFunctionReturn()` instead of return except for main().
1527530d85c1SBarry Smith 
1528530d85c1SBarry Smith       May be used before `PetscInitialize()`
15291957e957SBarry Smith 
1530530d85c1SBarry Smith       This is identical to `PetscFunctionBegin` except it labels the routine as a user
1531ac285190SBarry Smith       routine instead of as a PETSc library routine.
1532ac285190SBarry Smith 
1533586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, `PetscFunctionBeginHot`, `PetscStackPushNoCheck()`
1534a8d2bbe5SBarry Smith 
1535a8d2bbe5SBarry Smith M*/
15369371c9d4SSatish Balay   #define PetscFunctionBeginUser \
15379371c9d4SSatish Balay     do { \
1538362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 2, PETSC_FALSE); \
1539a8d2bbe5SBarry Smith       PetscRegister__FUNCT__(); \
1540a8d2bbe5SBarry Smith     } while (0)
1541a8d2bbe5SBarry Smith 
1542586f9135SBarry Smith   /*MC
1543586f9135SBarry Smith    PetscStackPush - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is
1544586f9135SBarry Smith    currently in the source code and verifies the memory is not corrupted.
1545586f9135SBarry Smith 
1546586f9135SBarry Smith    Synopsis:
1547586f9135SBarry Smith    #include <petscsys.h>
1548586f9135SBarry Smith    void PetscStackPush(char *funct)
1549586f9135SBarry Smith 
15507cdbe19fSJose E. Roman    Not Collective
15517cdbe19fSJose E. Roman 
1552586f9135SBarry Smith    Input Parameter:
1553586f9135SBarry Smith .  funct - the function name
1554586f9135SBarry Smith 
1555586f9135SBarry Smith    Level: developer
1556586f9135SBarry Smith 
1557586f9135SBarry Smith    Notes:
1558586f9135SBarry 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
1559586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1560586f9135SBarry Smith    help debug the problem.
1561586f9135SBarry Smith 
1562586f9135SBarry Smith    The default stack is a global variable called petscstack.
1563586f9135SBarry Smith 
1564586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1565586f9135SBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1566586f9135SBarry Smith M*/
15679371c9d4SSatish Balay   #define PetscStackPush(n) \
15689371c9d4SSatish Balay     do { \
1569362febeeSStefano Zampini       PetscStackPushNoCheck(n, 0, PETSC_FALSE); \
157015681b3cSBarry Smith       CHKMEMQ; \
157115681b3cSBarry Smith     } while (0)
15723a40ed3dSBarry Smith 
1573586f9135SBarry Smith   /*MC
1574586f9135SBarry Smith    PetscStackPop - Pops a function name from the PETSc default stack that tracks where the running program is
1575586f9135SBarry Smith    currently in the source code and verifies the memory is not corrupted.
1576586f9135SBarry Smith 
1577586f9135SBarry Smith    Synopsis:
1578586f9135SBarry Smith    #include <petscsys.h>
1579586f9135SBarry Smith    void PetscStackPop
1580586f9135SBarry Smith 
15817cdbe19fSJose E. Roman    Not Collective
15827cdbe19fSJose E. Roman 
1583586f9135SBarry Smith    Level: developer
1584586f9135SBarry Smith 
1585586f9135SBarry Smith    Notes:
1586586f9135SBarry 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
1587586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1588586f9135SBarry Smith    help debug the problem.
1589586f9135SBarry Smith 
1590586f9135SBarry Smith    The default stack is a global variable called petscstack.
1591586f9135SBarry Smith 
1592586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPopNoCheck()`, `PetscStackPush()`
1593586f9135SBarry Smith M*/
15949371c9d4SSatish Balay   #define PetscStackPop \
15959371c9d4SSatish Balay     do { \
1596441dd030SJed Brown       CHKMEMQ; \
1597362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
159815681b3cSBarry Smith     } while (0)
1599d64ed03dSBarry Smith 
160030de9b25SBarry Smith   /*MC
16010a57284eSJacob Faibussowitsch    PetscFunctionReturn - Last executable line of each PETSc function used for error
16020a57284eSJacob Faibussowitsch    handling. Replaces `return()`.
160330de9b25SBarry Smith 
160430de9b25SBarry Smith    Synopsis:
16050a57284eSJacob Faibussowitsch    #include <petscerror.h>
16060a57284eSJacob Faibussowitsch    void PetscFunctionReturn(...)
160730de9b25SBarry Smith 
1608667f096bSBarry Smith    Not Collective; No Fortran Support
1609eca87e8dSBarry Smith 
16105687f061SJacob Faibussowitsch    Level: beginner
16115687f061SJacob Faibussowitsch 
16120a57284eSJacob Faibussowitsch    Notes:
16130a57284eSJacob Faibussowitsch    This routine is a macro, so while it does not "return" anything itself, it does return from
16140a57284eSJacob Faibussowitsch    the function in the literal sense.
16150a57284eSJacob Faibussowitsch 
16160a57284eSJacob Faibussowitsch    Usually the return value is the integer literal `0` (for example in any function returning
16170a57284eSJacob Faibussowitsch    `PetscErrorCode`), however it is possible to return any arbitrary type. The arguments of
16180a57284eSJacob Faibussowitsch    this macro are placed before the `return` statement as-is.
16190a57284eSJacob Faibussowitsch 
16200a57284eSJacob Faibussowitsch    Any routine which returns via `PetscFunctionReturn()` must begin with a corresponding
16210a57284eSJacob Faibussowitsch    `PetscFunctionBegin`.
16220a57284eSJacob Faibussowitsch 
16230a57284eSJacob Faibussowitsch    For routines which return `void` use `PetscFunctionReturnVoid()` instead.
16240a57284eSJacob Faibussowitsch 
16250a57284eSJacob Faibussowitsch    Example Usage:
162630de9b25SBarry Smith .vb
16270a57284eSJacob Faibussowitsch    PetscErrorCode foo(int *x)
16280a57284eSJacob Faibussowitsch    {
16290a57284eSJacob Faibussowitsch      PetscFunctionBegin; // don't forget the begin!
16300a57284eSJacob Faibussowitsch      *x = 10;
16313ba16761SJacob Faibussowitsch      PetscFunctionReturn(PETSC_SUCCESS);
163230de9b25SBarry Smith    }
163330de9b25SBarry Smith .ve
163430de9b25SBarry Smith 
16350a57284eSJacob Faibussowitsch    May return any arbitrary type\:
16360a57284eSJacob Faibussowitsch .vb
16370a57284eSJacob Faibussowitsch   struct Foo
16380a57284eSJacob Faibussowitsch   {
16390a57284eSJacob Faibussowitsch     int x;
16400a57284eSJacob Faibussowitsch   };
16410a57284eSJacob Faibussowitsch 
16420a57284eSJacob Faibussowitsch   struct Foo make_foo(int value)
16430a57284eSJacob Faibussowitsch   {
16440a57284eSJacob Faibussowitsch     struct Foo f;
16450a57284eSJacob Faibussowitsch 
16460a57284eSJacob Faibussowitsch     PetscFunctionBegin;
16470a57284eSJacob Faibussowitsch     f.x = value;
16480a57284eSJacob Faibussowitsch     PetscFunctionReturn(f);
16490a57284eSJacob Faibussowitsch   }
16500a57284eSJacob Faibussowitsch .ve
16510a57284eSJacob Faibussowitsch 
16520a57284eSJacob Faibussowitsch .seealso: `PetscFunctionBegin`, `PetscFunctionBeginUser`, `PetscFunctionReturnVoid()`,
16530a57284eSJacob Faibussowitsch           `PetscStackPopNoCheck()`
165430de9b25SBarry Smith M*/
16555687f061SJacob Faibussowitsch   #define PetscFunctionReturn(...) \
16569371c9d4SSatish Balay     do { \
1657362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
16585687f061SJacob Faibussowitsch       return __VA_ARGS__; \
165927104ee2SJacob Faibussowitsch     } while (0)
1660d64ed03dSBarry Smith 
16610a57284eSJacob Faibussowitsch   /*MC
16620a57284eSJacob Faibussowitsch   PetscFunctionReturnVoid - Like `PetscFunctionReturn()` but returns `void`
16630a57284eSJacob Faibussowitsch 
16640a57284eSJacob Faibussowitsch   Synopsis:
16650a57284eSJacob Faibussowitsch   #include <petscerror.h>
16660a57284eSJacob Faibussowitsch   void PetscFunctionReturnVoid()
16670a57284eSJacob Faibussowitsch 
16680a57284eSJacob Faibussowitsch   Not Collective
16690a57284eSJacob Faibussowitsch 
16705687f061SJacob Faibussowitsch   Level: beginner
16715687f061SJacob Faibussowitsch 
16725687f061SJacob Faibussowitsch   Note:
16730a57284eSJacob Faibussowitsch   Behaves identically to `PetscFunctionReturn()` except that it returns `void`. That is, this
16745687f061SJacob Faibussowitsch   macro culminates with `return`.
16750a57284eSJacob Faibussowitsch 
16760a57284eSJacob Faibussowitsch   Example Usage:
16770a57284eSJacob Faibussowitsch .vb
16780a57284eSJacob Faibussowitsch   void foo()
16790a57284eSJacob Faibussowitsch   {
16800a57284eSJacob Faibussowitsch     PetscFunctionBegin; // must start with PetscFunctionBegin!
16810a57284eSJacob Faibussowitsch     bar();
16820a57284eSJacob Faibussowitsch     baz();
16830a57284eSJacob Faibussowitsch     PetscFunctionReturnVoid();
16840a57284eSJacob Faibussowitsch   }
16850a57284eSJacob Faibussowitsch .ve
16860a57284eSJacob Faibussowitsch 
16870a57284eSJacob Faibussowitsch .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, PetscFunctionBeginUser`
16880a57284eSJacob Faibussowitsch M*/
16899371c9d4SSatish Balay   #define PetscFunctionReturnVoid() \
16909371c9d4SSatish Balay     do { \
1691362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
169227104ee2SJacob Faibussowitsch       return; \
169327104ee2SJacob Faibussowitsch     } while (0)
169427104ee2SJacob Faibussowitsch #else /* PETSC_USE_DEBUG */
169527104ee2SJacob Faibussowitsch   #define PetscStackPushNoCheck(funct, petsc_routine, hot)
169637154d10SBarry Smith   #define PetscStackUpdateLine
1697792fecdfSBarry Smith   #define PetscStackPushExternal(funct)
16983ba16761SJacob Faibussowitsch   #define PetscStackPopNoCheck(...)
169927104ee2SJacob Faibussowitsch   #define PetscStackClearTop
17003a40ed3dSBarry Smith   #define PetscFunctionBegin
17010bdf7c52SPeter Brune   #define PetscFunctionBeginUser
1702a2f94806SJed Brown   #define PetscFunctionBeginHot
17030a57284eSJacob Faibussowitsch   #define PetscFunctionReturn(...)  return __VA_ARGS__
17045665465eSBarry Smith   #define PetscFunctionReturnVoid() return
1705812af9f3SBarry Smith   #define PetscStackPop             CHKMEMQ
1706812af9f3SBarry Smith   #define PetscStackPush(f)         CHKMEMQ
170727104ee2SJacob Faibussowitsch #endif /* PETSC_USE_DEBUG */
17083a40ed3dSBarry Smith 
17096d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
17103ba16761SJacob Faibussowitsch   #define PetscStackCallExternalVoid(...)
17113ba16761SJacob Faibussowitsch template <typename F, typename... Args>
17123ba16761SJacob Faibussowitsch void PetscCallExternal(F, Args...);
17136d210af2SJacob Faibussowitsch #else
1714586f9135SBarry Smith   /*MC
1715e77caa6dSBarry Smith     PetscStackCallExternalVoid - Calls an external library routine or user function after pushing the name of the routine on the stack.
1716eb6b5d47SBarry Smith 
1717eb6b5d47SBarry Smith    Input Parameters:
1718eb6b5d47SBarry Smith +   name - string that gives the name of the function being called
1719586f9135SBarry Smith -   routine - actual call to the routine, for example, functionname(a,b)
1720fd3f9acdSBarry Smith 
1721586f9135SBarry Smith    Level: developer
1722eb6b5d47SBarry Smith 
1723586f9135SBarry Smith    Note:
1724792fecdfSBarry Smith    Often one should use `PetscCallExternal()` instead. This routine is intended for external library routines that DO NOT return error codes
1725eb6b5d47SBarry Smith 
1726586f9135SBarry Smith    In debug mode this also checks the memory for corruption at the end of the function call.
1727586f9135SBarry Smith 
1728586f9135SBarry Smith    Certain external packages, such as BLAS/LAPACK may have their own macros for managing the call, error checking, etc.
1729586f9135SBarry Smith 
1730586f9135SBarry Smith    Developer Note:
1731586f9135SBarry 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.
1732586f9135SBarry Smith 
1733792fecdfSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscCallExternal()`, `PetscCallBLAS()`
1734586f9135SBarry Smith @*/
17353ba16761SJacob Faibussowitsch   #define PetscStackCallExternalVoid(name, ...) \
17369371c9d4SSatish Balay     do { \
173728770333SStefano Zampini       PetscStackPushExternal(name); \
17383ba16761SJacob Faibussowitsch       __VA_ARGS__; \
17399371c9d4SSatish Balay       PetscStackPop; \
17409371c9d4SSatish Balay     } while (0)
1741eb6b5d47SBarry Smith 
1742586f9135SBarry Smith   /*MC
1743792fecdfSBarry Smith     PetscCallExternal - Calls an external library routine that returns an error code after pushing the name of the routine on the stack.
1744fd3f9acdSBarry Smith 
1745fd3f9acdSBarry Smith    Input Parameters:
1746fd3f9acdSBarry Smith +   func-  name of the routine
1747586f9135SBarry Smith -   args - arguments to the routine
1748586f9135SBarry Smith 
1749586f9135SBarry Smith    Level: developer
1750fd3f9acdSBarry Smith 
175195452b02SPatrick Sanan    Notes:
1752e77caa6dSBarry Smith    This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not.
1753dbf62e16SBarry Smith 
1754586f9135SBarry Smith    In debug mode this also checks the memory for corruption at the end of the function call.
1755fd3f9acdSBarry Smith 
175687497f52SBarry Smith    Assumes the error return code of the function is an integer and that a value of 0 indicates success
175787497f52SBarry Smith 
1758586f9135SBarry Smith    Developer Note:
1759d5b43468SJose E. Roman    This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc.
1760586f9135SBarry Smith 
1761e77caa6dSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()`
1762586f9135SBarry Smith M*/
17639371c9d4SSatish Balay   #define PetscCallExternal(func, ...) \
17649371c9d4SSatish Balay     do { \
1765a74df02fSJacob Faibussowitsch       PetscStackPush(PetscStringize(func)); \
17663ba16761SJacob Faibussowitsch       int ierr_petsc_call_external_ = func(__VA_ARGS__); \
17671d4906efSStefano Zampini       PetscStackPop; \
17683ba16761SJacob 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_); \
1769fd3f9acdSBarry Smith     } while (0)
17706d210af2SJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */
177106d1fe2cSBarry Smith 
177206d1fe2cSBarry Smith #endif
1773