xref: /petsc/include/petscerror.h (revision aaa8cc7d2a5c3913edcbb923e20f154fe9c4aa65)
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 
1698921bdaSJacob Faibussowitsch #define SETERRQ1(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__)
1798921bdaSJacob Faibussowitsch #define SETERRQ2(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__)
1898921bdaSJacob Faibussowitsch #define SETERRQ3(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__)
1998921bdaSJacob Faibussowitsch #define SETERRQ4(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__)
2098921bdaSJacob Faibussowitsch #define SETERRQ5(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__)
2198921bdaSJacob Faibussowitsch #define SETERRQ6(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__)
2298921bdaSJacob Faibussowitsch #define SETERRQ7(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__)
2398921bdaSJacob Faibussowitsch #define SETERRQ8(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") SETERRQ(__VA_ARGS__)
2498921bdaSJacob Faibussowitsch #define SETERRQ9(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use SETERRQ() (since version 3.17)\"") 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 
483b1008b8SBarry Smith    Fortran Notes:
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 
158989c49a2SBarry Smith   No Fortran Support
159989c49a2SBarry Smith 
1602c71b3e2SJacob Faibussowitsch   Synopsis:
1612c71b3e2SJacob Faibussowitsch   #include <petscerror.h>
1622c71b3e2SJacob Faibussowitsch   void PetscCheck(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
1632c71b3e2SJacob Faibussowitsch 
1642c71b3e2SJacob Faibussowitsch   Collective
1652c71b3e2SJacob Faibussowitsch 
1662c71b3e2SJacob Faibussowitsch   Input Parameters:
1672c71b3e2SJacob Faibussowitsch + cond    - The boolean condition
1682c71b3e2SJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
1692c71b3e2SJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
1702c71b3e2SJacob Faibussowitsch - message - Error message in printf format
1712c71b3e2SJacob Faibussowitsch 
172667f096bSBarry Smith   Level: beginner
173667f096bSBarry Smith 
1742c71b3e2SJacob Faibussowitsch   Notes:
1752c71b3e2SJacob Faibussowitsch   Enabled in both optimized and debug builds.
1762c71b3e2SJacob Faibussowitsch 
17787497f52SBarry Smith   Calls `SETERRQ()` if the assertion fails, so can only be called from functions returning a
17887497f52SBarry Smith   `PetscErrorCode` (or equivalent type after conversion).
1792c71b3e2SJacob Faibussowitsch 
1804be741a6SBarry Smith  .seealso: `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheckAbort()`
1819566063dSJacob Faibussowitsch M*/
1829371c9d4SSatish Balay #define PetscCheck(cond, comm, ierr, ...) \
1839371c9d4SSatish Balay   do { \
1849371c9d4SSatish Balay     if (PetscUnlikely(!(cond))) SETERRQ(comm, ierr, __VA_ARGS__); \
1859371c9d4SSatish Balay   } while (0)
1862c71b3e2SJacob Faibussowitsch 
1872c71b3e2SJacob Faibussowitsch /*MC
1884be741a6SBarry Smith   PetscCheckAbort - Check that a particular condition is true, otherwise prints error and aborts
1894be741a6SBarry Smith 
190989c49a2SBarry Smith   No Fortran Support
191989c49a2SBarry Smith 
1924be741a6SBarry Smith   Synopsis:
1934be741a6SBarry Smith   #include <petscerror.h>
1944be741a6SBarry Smith   void PetscCheckAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
1954be741a6SBarry Smith 
1964be741a6SBarry Smith   Collective
1974be741a6SBarry Smith 
1984be741a6SBarry Smith   Input Parameters:
1994be741a6SBarry Smith + cond    - The boolean condition
2004be741a6SBarry Smith . comm    - The communicator on which the check can be collective on
2014be741a6SBarry Smith . ierr    - A nonzero error code, see include/petscerror.h for the complete list
2024be741a6SBarry Smith - message - Error message in printf format
2034be741a6SBarry Smith 
204667f096bSBarry Smith   Level: developer
205667f096bSBarry Smith 
2064be741a6SBarry Smith   Notes:
2074be741a6SBarry Smith   Enabled in both optimized and debug builds.
2084be741a6SBarry Smith 
20987497f52SBarry Smith   Calls `SETERRABORT()` if the assertion fails, can be called from a function that does not return an
21087497f52SBarry Smith   error code, such as a C++ constructor. usually `PetscCheck()` should be used.
2114be741a6SBarry Smith 
212989c49a2SBarry Smith .seealso: `PetscAssertAbort()`, `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheck()`, `SETERRABORT()`
2134be741a6SBarry Smith M*/
2149371c9d4SSatish Balay #define PetscCheckAbort(cond, comm, ierr, ...) \
2150e6b6b59SJacob Faibussowitsch   do { \
2160e6b6b59SJacob Faibussowitsch     if (PetscUnlikely(!(cond))) SETERRABORT(comm, ierr, __VA_ARGS__); \
217c7481402SJacob Faibussowitsch   } while (0)
2184be741a6SBarry Smith 
2194be741a6SBarry Smith /*MC
2209ace16cdSJacob Faibussowitsch   PetscAssert - Assert that a particular condition is true
2219ace16cdSJacob Faibussowitsch 
222989c49a2SBarry Smith   No Fortran Support
223989c49a2SBarry Smith 
2249ace16cdSJacob Faibussowitsch   Synopsis:
2259ace16cdSJacob Faibussowitsch   #include <petscerror.h>
2269ace16cdSJacob Faibussowitsch   void PetscAssert(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
2279ace16cdSJacob Faibussowitsch 
2289ace16cdSJacob Faibussowitsch   Collective
2299ace16cdSJacob Faibussowitsch 
2309ace16cdSJacob Faibussowitsch   Input Parameters:
2319ace16cdSJacob Faibussowitsch + cond    - The boolean condition
2329ace16cdSJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
2339ace16cdSJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
2349ace16cdSJacob Faibussowitsch - message - Error message in printf format
2359ace16cdSJacob Faibussowitsch 
236667f096bSBarry Smith   Level: beginner
237667f096bSBarry Smith 
2389ace16cdSJacob Faibussowitsch   Notes:
239c7481402SJacob Faibussowitsch   Equivalent to `PetscCheck()` if debugging is enabled, and `PetscAssume(cond)` otherwise.
2402c71b3e2SJacob Faibussowitsch 
24187497f52SBarry Smith   See `PetscCheck()` for usage and behaviour.
24287497f52SBarry Smith 
24387497f52SBarry Smith   This is needed instead of simply using `assert()` because this correctly handles the collective nature of errors under MPI
2449ace16cdSJacob Faibussowitsch 
2450e6b6b59SJacob Faibussowitsch .seealso: `PetscCheck()`, `SETERRQ()`, `PetscError()`, `PetscAssertAbort()`
2469566063dSJacob Faibussowitsch M*/
247c7481402SJacob Faibussowitsch #if PetscDefined(USE_DEBUG)
248c7481402SJacob Faibussowitsch   #define PetscAssert(cond, comm, ierr, ...) PetscCheck(cond, comm, ierr, __VA_ARGS__)
249c7481402SJacob Faibussowitsch #else
250c7481402SJacob Faibussowitsch   #define PetscAssert(cond, ...) PetscAssume(cond)
251c7481402SJacob Faibussowitsch #endif
2529ace16cdSJacob Faibussowitsch 
2539ace16cdSJacob Faibussowitsch /*MC
2540e6b6b59SJacob Faibussowitsch   PetscAssertAbort - Assert that a particular condition is true, otherwise prints error and aborts
2550e6b6b59SJacob Faibussowitsch 
256989c49a2SBarry Smith   No Fortran Support
257989c49a2SBarry Smith 
2580e6b6b59SJacob Faibussowitsch   Synopsis:
2590e6b6b59SJacob Faibussowitsch   #include <petscerror.h>
2600e6b6b59SJacob Faibussowitsch   void PetscAssertAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
2610e6b6b59SJacob Faibussowitsch 
2620e6b6b59SJacob Faibussowitsch   Collective
2630e6b6b59SJacob Faibussowitsch 
2640e6b6b59SJacob Faibussowitsch   Input Parameters:
2650e6b6b59SJacob Faibussowitsch + cond    - The boolean condition
2660e6b6b59SJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
2670e6b6b59SJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
2680e6b6b59SJacob Faibussowitsch - message - Error message in printf format
2690e6b6b59SJacob Faibussowitsch 
270667f096bSBarry Smith   Level: beginner
271667f096bSBarry Smith 
2720e6b6b59SJacob Faibussowitsch   Notes:
2730e6b6b59SJacob Faibussowitsch   Enabled only in debug builds. See `PetscCheckAbort()` for usage.
2740e6b6b59SJacob Faibussowitsch 
2750e6b6b59SJacob Faibussowitsch .seealso: `PetscCheckAbort()`, `PetscAssert()`, `PetscCheck()`, `SETERRABORT()`, `PetscError()`
2760e6b6b59SJacob Faibussowitsch M*/
2773ba16761SJacob Faibussowitsch #if PetscDefined(USE_DEBUG)
2783ba16761SJacob Faibussowitsch   #define PetscAssertAbort(cond, comm, ierr, ...) PetscCheckAbort(cond, comm, ierr, __VA_ARGS__)
2793ba16761SJacob Faibussowitsch #else
2803ba16761SJacob Faibussowitsch   #define PetscAssertAbort(cond, comm, ierr, ...) PetscAssume(cond)
2813ba16761SJacob Faibussowitsch #endif
2820e6b6b59SJacob Faibussowitsch 
2830e6b6b59SJacob Faibussowitsch /*MC
2843ba16761SJacob Faibussowitsch   PetscCall - Calls a PETSc function and then checks the resulting error code, if it is
2853ba16761SJacob Faibussowitsch   non-zero it calls the error handler and returns from the current function with the error
2863ba16761SJacob Faibussowitsch   code.
2879566063dSJacob Faibussowitsch 
2889566063dSJacob Faibussowitsch   Synopsis:
2899566063dSJacob Faibussowitsch   #include <petscerror.h>
29049c86fc7SBarry Smith   void PetscCall(PetscFunction(args))
2919566063dSJacob Faibussowitsch 
2929566063dSJacob Faibussowitsch   Not Collective
2939566063dSJacob Faibussowitsch 
2949566063dSJacob Faibussowitsch   Input Parameter:
29549c86fc7SBarry Smith . PetscFunction - any PETSc function that returns an error code
2969566063dSJacob Faibussowitsch 
297667f096bSBarry Smith   Level: beginner
298667f096bSBarry Smith 
2999566063dSJacob Faibussowitsch   Notes:
3009566063dSJacob Faibussowitsch   Once the error handler is called the calling function is then returned from with the given
30187497f52SBarry Smith   error code. Experienced users can set the error handler with `PetscPushErrorHandler()`.
3029566063dSJacob Faibussowitsch 
30387497f52SBarry Smith   `PetscCall()` cannot be used in functions returning a datatype not convertible to
30487497f52SBarry Smith   `PetscErrorCode`. For example, `PetscCall()` may not be used in functions returning void, use
3053ba16761SJacob Faibussowitsch   `PetscCallAbort()` or `PetscCallVoid()` in this case.
3069566063dSJacob Faibussowitsch 
3079566063dSJacob Faibussowitsch   Example Usage:
3089566063dSJacob Faibussowitsch .vb
3099566063dSJacob Faibussowitsch   PetscCall(PetscInitiailize(...)); // OK to call even when PETSc is not yet initialized!
3109566063dSJacob Faibussowitsch 
3119566063dSJacob Faibussowitsch   struct my_struct
3129566063dSJacob Faibussowitsch   {
3139566063dSJacob Faibussowitsch     void *data;
3149566063dSJacob Faibussowitsch   } my_complex_type;
3159566063dSJacob Faibussowitsch 
3169566063dSJacob Faibussowitsch   struct my_struct bar(void)
3179566063dSJacob Faibussowitsch   {
3186aad120cSJose E. Roman     PetscCall(foo(15)); // ERROR PetscErrorCode not convertible to struct my_struct!
3199566063dSJacob Faibussowitsch   }
3209566063dSJacob Faibussowitsch 
3216aad120cSJose E. Roman   PetscCall(bar()) // ERROR input not convertible to PetscErrorCode
3229566063dSJacob Faibussowitsch .ve
3239566063dSJacob Faibussowitsch 
32487497f52SBarry Smith   It is also possible to call this directly on a `PetscErrorCode` variable
32549c86fc7SBarry Smith .vb
32649c86fc7SBarry Smith   PetscCall(ierr);  // check if ierr is nonzero
32749c86fc7SBarry Smith .ve
32849c86fc7SBarry Smith 
329792fecdfSBarry Smith   Should not be used to call callback functions provided by users, `PetscCallBack()` should be used in that situation.
330ef1023bdSBarry Smith 
3316a8be23eSBarry Smith   `PetscUseTypeMethod()` or `PetscTryTypeMethod()` should be used when calling functions pointers contained in a PETSc object's `ops` array
3326a8be23eSBarry Smith 
33349c86fc7SBarry Smith   Fortran Notes:
33449c86fc7SBarry Smith     The Fortran function from which this is used must declare a variable PetscErrorCode ierr and ierr must be
33587497f52SBarry Smith     the final argument to the PETSc function being called.
33649c86fc7SBarry Smith 
33749c86fc7SBarry Smith     In the main program and in Fortran subroutines that do not have ierr as the final return parameter one
33887497f52SBarry Smith     should use `PetscCallA()`
33949c86fc7SBarry Smith 
34049c86fc7SBarry Smith   Example Fortran Usage:
34149c86fc7SBarry Smith .vb
34249c86fc7SBarry Smith   PetscErrorCode ierr
34349c86fc7SBarry Smith   Vec v
34449c86fc7SBarry Smith 
34549c86fc7SBarry Smith   ...
34649c86fc7SBarry Smith   PetscCall(VecShift(v,1.0,ierr))
34749c86fc7SBarry Smith   PetscCallA(VecShift(v,1.0,ierr))
34849c86fc7SBarry Smith .ve
34949c86fc7SBarry Smith 
350667f096bSBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`,
351667f096bSBarry Smith           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`,
3523ba16761SJacob Faibussowitsch           `CHKERRMPI()`, `PetscCallBack()`, `PetscCallAbort()`, `PetscCallVoid()`
3539566063dSJacob Faibussowitsch M*/
354ef1023bdSBarry Smith 
355ef1023bdSBarry Smith /*MC
356989c49a2SBarry Smith    PetscCallA - Fortran-only macro that should be used in the main program to call PETSc functions instead of using
357989c49a2SBarry Smith    PetscCall() which should be used in other Fortran subroutines
358989c49a2SBarry Smith 
359989c49a2SBarry Smith    Synopsis:
360989c49a2SBarry Smith    #include <petscsys.h>
361989c49a2SBarry Smith    PetscErrorCode PetscCallA(PetscFunction(arguments,ierr))
362989c49a2SBarry Smith 
363989c49a2SBarry Smith    Collective
364989c49a2SBarry Smith 
365989c49a2SBarry Smith    Input Parameter:
366989c49a2SBarry Smith .  PetscFunction(arguments,ierr) - the call to the function
367989c49a2SBarry Smith 
368989c49a2SBarry Smith   Level: beginner
369989c49a2SBarry Smith 
370989c49a2SBarry Smith    Notes:
371989c49a2SBarry Smith    This should only be used with Fortran. With C/C++, use `PetscCall()` always.
372989c49a2SBarry Smith 
373989c49a2SBarry Smith    Use `SETERRA()` to set an error in a Fortran main program and `SETERRQ()` in Fortran subroutines
374989c49a2SBarry Smith 
375989c49a2SBarry Smith .seealso: `SETERRQ()`, `SETERRA()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()`
376989c49a2SBarry Smith M*/
377989c49a2SBarry Smith 
378989c49a2SBarry Smith /*MC
379792fecdfSBarry 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
380ef1023bdSBarry Smith   handler and returns from the current function with the error code.
381ef1023bdSBarry Smith 
382989c49a2SBarry Smith   No Fortran Support
383989c49a2SBarry Smith 
384ef1023bdSBarry Smith   Synopsis:
385ef1023bdSBarry Smith   #include <petscerror.h>
386792fecdfSBarry Smith   void PetscCallBack(const char *functionname,PetscFunction(args))
387ef1023bdSBarry Smith 
388ef1023bdSBarry Smith   Not Collective
389ef1023bdSBarry Smith 
390ef1023bdSBarry Smith   Input Parameters:
391ef1023bdSBarry Smith + functionname - the name of the function being called, this can be a string with spaces that describes the meaning of the callback
39287497f52SBarry Smith - PetscFunction - user provided callback function that returns an error code
393ef1023bdSBarry Smith 
394ef1023bdSBarry Smith   Example Usage:
395ef1023bdSBarry Smith .vb
396792fecdfSBarry Smith   PetscCallBack("XXX callback to do something",a->callback(...));
397ef1023bdSBarry Smith .ve
398ef1023bdSBarry Smith 
399ef1023bdSBarry Smith   Level: developer
400ef1023bdSBarry Smith 
401667f096bSBarry Smith   Notes:
402667f096bSBarry Smith   Once the error handler is called the calling function is then returned from with the given
403667f096bSBarry Smith   error code. Experienced users can set the error handler with `PetscPushErrorHandler()`.
404667f096bSBarry Smith 
405667f096bSBarry Smith   `PetscCallBack()` should only be called in PETSc when a call is being made to a user provided call-back routine.
406667f096bSBarry Smith 
40787497f52SBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscCall()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`
408ef1023bdSBarry Smith           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, `CHKERRMPI()`, `PetscCall()`
409ef1023bdSBarry Smith M*/
410ef1023bdSBarry Smith 
4113ba16761SJacob Faibussowitsch /*MC
4123ba16761SJacob Faibussowitsch   PetscCallVoid - Like `PetscCall()` but for functions returning `void`
4133ba16761SJacob Faibussowitsch 
414989c49a2SBarry Smith   No Fortran Support
415989c49a2SBarry Smith 
4163ba16761SJacob Faibussowitsch   Synopsis:
4173ba16761SJacob Faibussowitsch   #include <petscerror.h>
4183ba16761SJacob Faibussowitsch   void PetscCall(PetscFunction(args))
4193ba16761SJacob Faibussowitsch 
4203ba16761SJacob Faibussowitsch   Not Collective
4213ba16761SJacob Faibussowitsch 
4223ba16761SJacob Faibussowitsch   Input Parameter:
4233ba16761SJacob Faibussowitsch . PetscFunction - any PETSc function that returns an error code
4243ba16761SJacob Faibussowitsch 
4253ba16761SJacob Faibussowitsch   Example Usage:
4263ba16761SJacob Faibussowitsch .vb
4273ba16761SJacob Faibussowitsch   void foo()
4283ba16761SJacob Faibussowitsch   {
4293ba16761SJacob Faibussowitsch     KSP ksp;
4303ba16761SJacob Faibussowitsch 
4313ba16761SJacob Faibussowitsch     PetscFunctionBeginUser;
4323ba16761SJacob Faibussowitsch     // OK, properly handles PETSc error codes
4333ba16761SJacob Faibussowitsch     PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp));
4343ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
4353ba16761SJacob Faibussowitsch   }
4363ba16761SJacob Faibussowitsch 
4373ba16761SJacob Faibussowitsch   PetscErrorCode bar()
4383ba16761SJacob Faibussowitsch   {
4393ba16761SJacob Faibussowitsch     KSP ksp;
4403ba16761SJacob Faibussowitsch 
4413ba16761SJacob Faibussowitsch     PetscFunctionBeginUser;
4423ba16761SJacob Faibussowitsch     // ERROR, Non-void function 'bar' should return a value
4433ba16761SJacob Faibussowitsch     PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp));
4443ba16761SJacob Faibussowitsch     // OK, returning PetscErrorCode
4453ba16761SJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
4463ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
4473ba16761SJacob Faibussowitsch   }
448667f096bSBarry Smith .ve
4493ba16761SJacob Faibussowitsch 
4503ba16761SJacob Faibussowitsch   Level: beginner
4513ba16761SJacob Faibussowitsch 
452667f096bSBarry Smith   Notes:
453667f096bSBarry Smith   Has identical usage to `PetscCall()`, except that it returns `void` on error instead of a
454667f096bSBarry Smith   `PetscErrorCode`. See `PetscCall()` for more detailed discussion.
455667f096bSBarry Smith 
456667f096bSBarry Smith   Note that users should prefer `PetscCallAbort()` to this routine. While this routine does
457667f096bSBarry Smith   "handle" errors by returning from the enclosing function, it effectively gobbles the
458667f096bSBarry Smith   error. Since the enclosing function itself returns `void`, its callers have no way of knowing
459667f096bSBarry Smith   that the routine returned early due to an error. `PetscCallAbort()` at least ensures that the
460667f096bSBarry Smith   program crashes gracefully.
461667f096bSBarry Smith 
4623ba16761SJacob Faibussowitsch .seealso: `PetscCall()`, `PetscErrorCode`
4633ba16761SJacob Faibussowitsch M*/
4649566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
4659566063dSJacob Faibussowitsch void PetscCall(PetscErrorCode);
466792fecdfSBarry Smith void PetscCallBack(const char *, PetscErrorCode);
4679566063dSJacob Faibussowitsch void PetscCallVoid(PetscErrorCode);
4689566063dSJacob Faibussowitsch #else
4699371c9d4SSatish Balay   #define PetscCall(...) \
4709371c9d4SSatish Balay     do { \
4713ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_q_; \
47237154d10SBarry Smith       PetscStackUpdateLine; \
4733ba16761SJacob Faibussowitsch       ierr_petsc_call_q_ = __VA_ARGS__; \
4743ba16761SJacob 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, " "); \
4759566063dSJacob Faibussowitsch     } while (0)
4769371c9d4SSatish Balay   #define PetscCallBack(function, ...) \
4779371c9d4SSatish Balay     do { \
4783ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_q_; \
479e33ced7fSLisandro Dalcin       PetscStackUpdateLine; \
480e33ced7fSLisandro Dalcin       PetscStackPushExternal(function); \
4813ba16761SJacob Faibussowitsch       ierr_petsc_call_q_ = __VA_ARGS__; \
482ef1023bdSBarry Smith       PetscStackPop; \
4833ba16761SJacob 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, " "); \
484ef1023bdSBarry Smith     } while (0)
4859371c9d4SSatish Balay   #define PetscCallVoid(...) \
4869371c9d4SSatish Balay     do { \
4873ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_void_; \
488e33ced7fSLisandro Dalcin       PetscStackUpdateLine; \
4893ba16761SJacob Faibussowitsch       ierr_petsc_call_void_ = __VA_ARGS__; \
4903ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_void_ != PETSC_SUCCESS)) { \
4913ba16761SJacob Faibussowitsch         ierr_petsc_call_void_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_void_, PETSC_ERROR_REPEAT, " "); \
4923ba16761SJacob Faibussowitsch         (void)ierr_petsc_call_void_; \
4939371c9d4SSatish Balay         return; \
4949371c9d4SSatish Balay       } \
4959566063dSJacob Faibussowitsch     } while (0)
4969566063dSJacob Faibussowitsch #endif
4979566063dSJacob Faibussowitsch 
4989566063dSJacob Faibussowitsch /*MC
4999566063dSJacob Faibussowitsch   CHKERRQ - Checks error code returned from PETSc function
50030de9b25SBarry Smith 
50130de9b25SBarry Smith   Synopsis:
502aaa7dc30SBarry Smith   #include <petscsys.h>
5039566063dSJacob Faibussowitsch   void CHKERRQ(PetscErrorCode ierr)
5049566063dSJacob Faibussowitsch 
5059566063dSJacob Faibussowitsch   Not Collective
5069566063dSJacob Faibussowitsch 
5072fe279fdSBarry Smith   Input Parameter:
5089566063dSJacob Faibussowitsch . ierr - nonzero error code
5099566063dSJacob Faibussowitsch 
5109566063dSJacob Faibussowitsch   Level: deprecated
5119566063dSJacob Faibussowitsch 
512667f096bSBarry Smith   Note:
513667f096bSBarry Smith   Deprecated in favor of `PetscCall()`. This routine behaves identically to it.
514667f096bSBarry Smith 
515db781477SPatrick Sanan .seealso: `PetscCall()`
5169566063dSJacob Faibussowitsch M*/
5179566063dSJacob Faibussowitsch #define CHKERRQ(...) PetscCall(__VA_ARGS__)
5189566063dSJacob Faibussowitsch #define CHKERRV(...) PetscCallVoid(__VA_ARGS__)
5199566063dSJacob Faibussowitsch 
520db9cea48SBarry Smith PETSC_EXTERN void PetscMPIErrorString(PetscMPIInt, char *);
521db9cea48SBarry Smith 
5229566063dSJacob Faibussowitsch /*MC
5239566063dSJacob Faibussowitsch   PetscCallMPI - Checks error code returned from MPI calls, if non-zero it calls the error
5249566063dSJacob Faibussowitsch   handler and then returns
5259566063dSJacob Faibussowitsch 
5269566063dSJacob Faibussowitsch   Synopsis:
5279566063dSJacob Faibussowitsch   #include <petscerror.h>
52849c86fc7SBarry Smith   void PetscCallMPI(MPI_Function(args))
52930de9b25SBarry Smith 
530eca87e8dSBarry Smith   Not Collective
53130de9b25SBarry Smith 
5322fe279fdSBarry Smith   Input Parameter:
53349c86fc7SBarry Smith . MPI_Function - an MPI function that returns an MPI error code
53430de9b25SBarry Smith 
535667f096bSBarry Smith   Level: beginner
536667f096bSBarry Smith 
5379566063dSJacob Faibussowitsch   Notes:
53887497f52SBarry Smith   Always returns the error code `PETSC_ERR_MPI`; the MPI error code and string are embedded in
5399566063dSJacob Faibussowitsch   the string error message. Do not use this to call any other routines (for example PETSc
5403ba16761SJacob Faibussowitsch   routines), it should only be used for direct MPI calls. The user may configure PETSc with the
5413ba16761SJacob Faibussowitsch   `--with-strict-petscerrorcode` option to check this at compile-time, otherwise they must
5429566063dSJacob Faibussowitsch   check this themselves.
5439566063dSJacob Faibussowitsch 
544*aaa8cc7dSPierre Jolivet   This routine can only be used in functions returning `PetscErrorCode` themselves. If the
5453ba16761SJacob Faibussowitsch   calling function returns a different type, use `PetscCallMPIAbort()` instead.
5463ba16761SJacob Faibussowitsch 
5479566063dSJacob Faibussowitsch   Example Usage:
5489566063dSJacob Faibussowitsch .vb
5499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(...)); // OK, calling MPI function
5509566063dSJacob Faibussowitsch 
5519566063dSJacob Faibussowitsch   PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead!
5529566063dSJacob Faibussowitsch .ve
5539566063dSJacob Faibussowitsch 
55449c86fc7SBarry Smith   Fortran Notes:
55587497f52SBarry Smith     The Fortran function from which this is used must declare a variable `PetscErrorCode` ierr and ierr must be
55649c86fc7SBarry Smith     the final argument to the MPI function being called.
55749c86fc7SBarry Smith 
55849c86fc7SBarry Smith     In the main program and in Fortran subroutines that do not have ierr as the final return parameter one
55987497f52SBarry Smith     should use `PetscCallMPIA()`
56049c86fc7SBarry Smith 
56149c86fc7SBarry Smith   Fortran Usage:
56249c86fc7SBarry Smith .vb
56349c86fc7SBarry Smith   PetscErrorCode ierr or integer ierr
56449c86fc7SBarry Smith   ...
56549c86fc7SBarry Smith   PetscCallMPI(MPI_Comm_size(...,ierr))
56649c86fc7SBarry Smith   PetscCallMPIA(MPI_Comm_size(...,ierr)) ! Will abort after calling error handler
56749c86fc7SBarry Smith 
56849c86fc7SBarry Smith   PetscCallMPI(MPI_Comm_size(...,eflag)) ! ERROR, final argument must be ierr
56949c86fc7SBarry Smith .ve
57049c86fc7SBarry Smith 
571db781477SPatrick Sanan .seealso: `SETERRMPI()`, `PetscCall()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`,
5723ba16761SJacob Faibussowitsch           `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
5733ba16761SJacob Faibussowitsch           `PetscError()`, `CHKMEMQ`
5743ba16761SJacob Faibussowitsch M*/
5753ba16761SJacob Faibussowitsch 
5763ba16761SJacob Faibussowitsch /*MC
5773ba16761SJacob Faibussowitsch   PetscCallMPIAbort - Like `PetscCallMPI()` but calls `MPI_Abort()` on error
5783ba16761SJacob Faibussowitsch 
5793ba16761SJacob Faibussowitsch   Synopsis:
5803ba16761SJacob Faibussowitsch   #include <petscerror.h>
5813ba16761SJacob Faibussowitsch   void PetscCallMPIAbort(MPI_Comm comm, MPI_Function(args))
5823ba16761SJacob Faibussowitsch 
5833ba16761SJacob Faibussowitsch   Not Collective
5843ba16761SJacob Faibussowitsch 
5853ba16761SJacob Faibussowitsch   Input Parameters:
5863ba16761SJacob Faibussowitsch + comm         - the MPI communicator to abort on
5873ba16761SJacob Faibussowitsch - MPI_Function - an MPI function that returns an MPI error code
5883ba16761SJacob Faibussowitsch 
589667f096bSBarry Smith   Level: beginner
590667f096bSBarry Smith 
5913ba16761SJacob Faibussowitsch   Notes:
5923ba16761SJacob Faibussowitsch   Usage is identical to `PetscCallMPI()`. See `PetscCallMPI()` for detailed discussion.
5933ba16761SJacob Faibussowitsch 
5943ba16761SJacob Faibussowitsch   This routine may be used in functions returning `void` or other non-`PetscErrorCode` types.
5953ba16761SJacob Faibussowitsch 
596989c49a2SBarry Smith   Fortran Note:
597989c49a2SBarry Smith   In Fortran this is called `PetscCallMPIA()` and is intended to be used in the main program while `PetscCallMPI()` is
598989c49a2SBarry Smith   used in Fortran subroutines.
599989c49a2SBarry Smith 
600989c49a2SBarry Smith   Developer Note:
601989c49a2SBarry Smith   This should have the same name in Fortran.
602989c49a2SBarry Smith 
6033ba16761SJacob Faibussowitsch .seealso: `PetscCallMPI()`, `PetscCallAbort()`, `SETERRABORT()`
60430de9b25SBarry Smith M*/
6053fcd9f07SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
60663a3b9bcSJacob Faibussowitsch void PetscCallMPI(PetscMPIInt);
6073ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm, PetscMPIInt);
608064a246eSJacob Faibussowitsch #else
6093ba16761SJacob Faibussowitsch   #define PetscCallMPI_Private(__PETSC_STACK_POP_FUNC__, __SETERR_FUNC__, __COMM__, ...) \
6109371c9d4SSatish Balay     do { \
6113ba16761SJacob Faibussowitsch       PetscMPIInt ierr_petsc_call_mpi_; \
612ef1023bdSBarry Smith       PetscStackUpdateLine; \
613792fecdfSBarry Smith       PetscStackPushExternal("MPI function"); \
614d71ae5a4SJacob Faibussowitsch       { \
6153ba16761SJacob Faibussowitsch         ierr_petsc_call_mpi_ = __VA_ARGS__; \
616d71ae5a4SJacob Faibussowitsch       } \
6173ba16761SJacob Faibussowitsch       __PETSC_STACK_POP_FUNC__; \
6183ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_mpi_ != MPI_SUCCESS)) { \
6193ba16761SJacob Faibussowitsch         char petsc_mpi_7_errorstring[2 * MPI_MAX_ERROR_STRING]; \
6203ba16761SJacob Faibussowitsch         PetscMPIErrorString(ierr_petsc_call_mpi_, (char *)petsc_mpi_7_errorstring); \
6213ba16761SJacob Faibussowitsch         __SETERR_FUNC__(__COMM__, PETSC_ERR_MPI, "MPI error %d %s", (int)ierr_petsc_call_mpi_, petsc_mpi_7_errorstring); \
6223fcd9f07SJacob Faibussowitsch       } \
6233fcd9f07SJacob Faibussowitsch     } while (0)
6243ba16761SJacob Faibussowitsch 
6253ba16761SJacob Faibussowitsch   #define PetscCallMPI(...)            PetscCallMPI_Private(PetscStackPop, SETERRQ, PETSC_COMM_SELF, __VA_ARGS__)
6263ba16761SJacob Faibussowitsch   #define PetscCallMPIAbort(comm, ...) PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRABORT, comm, __VA_ARGS__)
6279566063dSJacob Faibussowitsch #endif
628064a246eSJacob Faibussowitsch 
6297037b0edSPatrick Sanan /*MC
6309566063dSJacob Faibussowitsch   CHKERRMPI - Checks error code returned from MPI calls, if non-zero it calls the error
6319566063dSJacob Faibussowitsch   handler and then returns
6329566063dSJacob Faibussowitsch 
6339566063dSJacob Faibussowitsch   Synopsis:
6349566063dSJacob Faibussowitsch   #include <petscerror.h>
6359566063dSJacob Faibussowitsch   void CHKERRMPI(PetscErrorCode ierr)
6369566063dSJacob Faibussowitsch 
6379566063dSJacob Faibussowitsch   Not Collective
6389566063dSJacob Faibussowitsch 
6399566063dSJacob Faibussowitsch   Input Parameter:
6409566063dSJacob Faibussowitsch . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
6419566063dSJacob Faibussowitsch 
6429566063dSJacob Faibussowitsch   Level: deprecated
6439566063dSJacob Faibussowitsch 
644667f096bSBarry Smith   Note:
645667f096bSBarry Smith   Deprecated in favor of `PetscCallMPI()`. This routine behaves identically to it.
646667f096bSBarry Smith 
647db781477SPatrick Sanan .seealso: `PetscCallMPI()`
6489566063dSJacob Faibussowitsch M*/
6499566063dSJacob Faibussowitsch #define CHKERRMPI(...) PetscCallMPI(__VA_ARGS__)
6509566063dSJacob Faibussowitsch 
6519566063dSJacob Faibussowitsch /*MC
652989c49a2SBarry Smith   PetscCallAbort - Checks error code returned from PETSc function, if non-zero it aborts immediately by calling `MPI_Abort()`
6539566063dSJacob Faibussowitsch 
6549566063dSJacob Faibussowitsch   Synopsis:
6559566063dSJacob Faibussowitsch   #include <petscerror.h>
6569566063dSJacob Faibussowitsch   void PetscCallAbort(MPI_Comm comm, PetscErrorCode ierr)
6579566063dSJacob Faibussowitsch 
658c3339decSBarry Smith   Collective
6599566063dSJacob Faibussowitsch 
6609566063dSJacob Faibussowitsch   Input Parameters:
6619566063dSJacob Faibussowitsch + comm - the MPI communicator on which to abort
6629566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
6639566063dSJacob Faibussowitsch 
664667f096bSBarry Smith   Level: intermediate
665667f096bSBarry Smith 
6669566063dSJacob Faibussowitsch   Notes:
66787497f52SBarry Smith   This macro has identical type and usage semantics to `PetscCall()` with the important caveat
6689566063dSJacob Faibussowitsch   that this macro does not return. Instead, if ierr is nonzero it calls the PETSc error handler
66987497f52SBarry Smith   and then immediately calls `MPI_Abort()`. It can therefore be used anywhere.
6709566063dSJacob Faibussowitsch 
67187497f52SBarry Smith   As per `MPI_Abort()` semantics the communicator passed must be valid, although there is currently
67287497f52SBarry Smith   no attempt made at handling any potential errors from `MPI_Abort()`. Note that while
67387497f52SBarry Smith   `MPI_Abort()` is required to terminate only those processes which reside on comm, it is often
67487497f52SBarry Smith   the case that `MPI_Abort()` terminates *all* processes.
6759566063dSJacob Faibussowitsch 
6769566063dSJacob Faibussowitsch   Example Usage:
6779566063dSJacob Faibussowitsch .vb
6789566063dSJacob Faibussowitsch   PetscErrorCode boom(void) { return PETSC_ERR_MEM; }
6799566063dSJacob Faibussowitsch 
6809566063dSJacob Faibussowitsch   void foo(void)
6819566063dSJacob Faibussowitsch   {
6829566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type
6839566063dSJacob Faibussowitsch   }
6849566063dSJacob Faibussowitsch 
6859566063dSJacob Faibussowitsch   double bar(void)
6869566063dSJacob Faibussowitsch   {
6879566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type
6889566063dSJacob Faibussowitsch   }
6899566063dSJacob Faibussowitsch 
6909566063dSJacob Faibussowitsch   PetscCallAbort(MPI_COMM_NULL,boom()); // ERROR, communicator should be valid
6919566063dSJacob Faibussowitsch 
6929566063dSJacob Faibussowitsch   struct baz
6939566063dSJacob Faibussowitsch   {
6949566063dSJacob Faibussowitsch     baz()
6959566063dSJacob Faibussowitsch     {
6969566063dSJacob Faibussowitsch       PetscCallAbort(PETSC_COMM_SELF,boom()); // OK
6979566063dSJacob Faibussowitsch     }
6989566063dSJacob Faibussowitsch 
6999566063dSJacob Faibussowitsch     ~baz()
7009566063dSJacob Faibussowitsch     {
7019566063dSJacob Faibussowitsch       PetscCallAbort(PETSC_COMM_SELF,boom()); // OK (in fact the only way to handle PETSc errors)
7029566063dSJacob Faibussowitsch     }
7039566063dSJacob Faibussowitsch   };
7049566063dSJacob Faibussowitsch .ve
7059566063dSJacob Faibussowitsch 
706989c49a2SBarry Smith   Fortran Note:
707989c49a2SBarry Smith   Use `PetscCallA()`.
708989c49a2SBarry Smith 
709989c49a2SBarry Smith   Developer Note:
710989c49a2SBarry Smith   This should have the same name in Fortran as in C.
711989c49a2SBarry Smith 
712db781477SPatrick Sanan .seealso: `SETERRABORT()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`,
7135687f061SJacob Faibussowitsch           `SETERRQ()`, `CHKMEMQ`, `PetscCallMPI()`, `PetscCallCXXAbort()`
7149566063dSJacob Faibussowitsch M*/
7159566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
7169566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm, PetscErrorCode);
7179566063dSJacob Faibussowitsch void PetscCallContinue(PetscErrorCode);
7189566063dSJacob Faibussowitsch #else
7199371c9d4SSatish Balay   #define PetscCallAbort(comm, ...) \
7209371c9d4SSatish Balay     do { \
7217da6b3ccSLisandro Dalcin       PetscErrorCode ierr_petsc_call_abort_; \
7223ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
7237da6b3ccSLisandro Dalcin       ierr_petsc_call_abort_ = __VA_ARGS__; \
7243ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_abort_ != PETSC_SUCCESS)) { \
7253ba16761SJacob Faibussowitsch         ierr_petsc_call_abort_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_abort_, PETSC_ERROR_REPEAT, " "); \
7263ba16761SJacob Faibussowitsch         (void)MPI_Abort(comm, (PetscMPIInt)ierr_petsc_call_abort_); \
7279566063dSJacob Faibussowitsch       } \
7289566063dSJacob Faibussowitsch     } while (0)
7299371c9d4SSatish Balay   #define PetscCallContinue(...) \
7309371c9d4SSatish Balay     do { \
7317da6b3ccSLisandro Dalcin       PetscErrorCode ierr_petsc_call_continue_; \
7323ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
7337da6b3ccSLisandro Dalcin       ierr_petsc_call_continue_ = __VA_ARGS__; \
7343ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_continue_ != PETSC_SUCCESS)) { \
7353ba16761SJacob Faibussowitsch         ierr_petsc_call_continue_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_continue_, PETSC_ERROR_REPEAT, " "); \
7363ba16761SJacob Faibussowitsch         (void)ierr_petsc_call_continue_; \
7373ba16761SJacob Faibussowitsch       } \
7389566063dSJacob Faibussowitsch     } while (0)
7399566063dSJacob Faibussowitsch #endif
7409566063dSJacob Faibussowitsch 
7419566063dSJacob Faibussowitsch /*MC
7429566063dSJacob Faibussowitsch   CHKERRABORT - Checks error code returned from PETSc function. If non-zero it aborts immediately.
7439566063dSJacob Faibussowitsch 
7449566063dSJacob Faibussowitsch   Synopsis:
7459566063dSJacob Faibussowitsch   #include <petscerror.h>
7469566063dSJacob Faibussowitsch   void CHKERRABORT(MPI_Comm comm, PetscErrorCode ierr)
7479566063dSJacob Faibussowitsch 
7489566063dSJacob Faibussowitsch   Not Collective
7499566063dSJacob Faibussowitsch 
7509566063dSJacob Faibussowitsch   Input Parameters:
7519566063dSJacob Faibussowitsch + comm - the MPI communicator
7529566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
7539566063dSJacob Faibussowitsch 
7549566063dSJacob Faibussowitsch   Level: deprecated
7559566063dSJacob Faibussowitsch 
756667f096bSBarry Smith   Note:
757667f096bSBarry Smith   Deprecated in favor of `PetscCallAbort()`. This routine behaves identically to it.
758667f096bSBarry Smith 
759db781477SPatrick Sanan .seealso: `PetscCallAbort()`
7609566063dSJacob Faibussowitsch M*/
7619566063dSJacob Faibussowitsch #define CHKERRABORT(comm, ...) PetscCallAbort(comm, __VA_ARGS__)
7629566063dSJacob Faibussowitsch #define CHKERRCONTINUE(...)    PetscCallContinue(__VA_ARGS__)
7639566063dSJacob Faibussowitsch 
7649566063dSJacob Faibussowitsch /*MC
76587497f52SBarry Smith    CHKERRA - Fortran-only replacement for use of `CHKERRQ()` in the main program, which aborts immediately
766f388eb8bSPatrick Sanan 
767f388eb8bSPatrick Sanan    Synopsis:
768f388eb8bSPatrick Sanan    #include <petscsys.h>
769f388eb8bSPatrick Sanan    PetscErrorCode CHKERRA(PetscErrorCode ierr)
770f388eb8bSPatrick Sanan 
771f388eb8bSPatrick Sanan    Not Collective
772f388eb8bSPatrick Sanan 
7732fe279fdSBarry Smith    Input Parameter:
774f388eb8bSPatrick Sanan .  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
775f388eb8bSPatrick Sanan 
77687497f52SBarry Smith   Level: deprecated
777f388eb8bSPatrick Sanan 
77887497f52SBarry Smith    Note:
77987497f52SBarry Smith    This macro is rarely needed, normal usage is `PetscCallA()` in the main Fortran program.
780f388eb8bSPatrick Sanan 
781989c49a2SBarry Smith    Developer Note:
782989c49a2SBarry Smith    Why isn't this named `CHKERRABORT()` in Fortran?
783989c49a2SBarry Smith 
78487497f52SBarry Smith .seealso: `PetscCall()`, `PetscCallA()`, `PetscCallAbort()`, `CHKERRQ()`, `SETERRA()`, `SETERRQ()`, `SETERRABORT()`
785f388eb8bSPatrick Sanan M*/
786f388eb8bSPatrick Sanan 
7871e4f893aSSatish Balay PETSC_EXTERN PetscBool petscwaitonerrorflg;
7881e4f893aSSatish Balay PETSC_EXTERN PetscBool petscindebugger;
7897c66cc67SJunchao Zhang 
7907c66cc67SJunchao Zhang /*MC
791667f096bSBarry Smith    PETSCABORT - Call `MPI_Abort()` with an informative error code
7927c66cc67SJunchao Zhang 
793989c49a2SBarry Smith    No Fortran Support
794989c49a2SBarry Smith 
7957c66cc67SJunchao Zhang    Synopsis:
7967c66cc67SJunchao Zhang    #include <petscsys.h>
7977c66cc67SJunchao Zhang    PETSCABORT(MPI_Comm comm, PetscErrorCode ierr)
7987c66cc67SJunchao Zhang 
7997c66cc67SJunchao Zhang    Collective
8007c66cc67SJunchao Zhang 
8017c66cc67SJunchao Zhang    Input Parameters:
8027c66cc67SJunchao Zhang +  comm - A communicator, so that the error can be collective
8037c66cc67SJunchao Zhang -  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
8047c66cc67SJunchao Zhang 
805bf4d2887SBarry Smith    Level: advanced
8067c66cc67SJunchao Zhang 
807bf4d2887SBarry Smith    Notes:
808989c49a2SBarry Smith    If the option `-start_in_debugger` was used then this calls `abort()` to stop the program in the debugger.
809bf4d2887SBarry Smith 
810989c49a2SBarry Smith    if `PetscCIEnabledPortableErrorOutput` is set, which means the code is running in the PETSc test harness (make test),
811989c49a2SBarry Smith    and `comm` is `MPI_COMM_WORLD` it strives to exit cleanly without calling `MPI_Abort()` and instead calling `MPI_Finalize()`.
812660278c0SBarry Smith 
813989c49a2SBarry Smith    This is currently only used when an error propagates up to the C `main()` program and is detected by a `PetscCall()`, `PetscCallMPI()`,
814989c49a2SBarry Smith    or is set in `main()` with `SETERRQ()`. Abort calls such as `SETERRABORT()`,
815989c49a2SBarry Smith    `PetscCheckAbort()`, `PetscCallMPIAbort()`, and `PetscCallAbort()` always call `MPI_Abort()` and do not have any special
816989c49a2SBarry Smith    handling for the test harness.
817989c49a2SBarry Smith 
818989c49a2SBarry Smith    Developer Note:
819989c49a2SBarry Smith    Should the other abort calls also pass through this call instead of calling `MPI_Abort()` directly?
820989c49a2SBarry Smith 
821989c49a2SBarry Smith .seealso: `PetscError()`, `PetscCall()`, `SETERRABORT()`, `PetscCheckAbort()`, `PetscCallMPIAbort()`, `PetscCall()`, `PetscCallMPI()`,
822989c49a2SBarry Smith           `PetscCallAbort()`, `MPI_Abort()`
8237c66cc67SJunchao Zhang M*/
82429f88feaSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
82529f88feaSJacob Faibussowitsch void PETSCABORT(MPI_Comm, PetscErrorCode);
82629f88feaSJacob Faibussowitsch #else
8279371c9d4SSatish Balay   #define PETSCABORT(comm, ...) \
8289371c9d4SSatish Balay     do { \
8293ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_abort_; \
8303ba16761SJacob Faibussowitsch       if (petscwaitonerrorflg) { ierr_petsc_abort_ = PetscSleep(1000); } \
8313ba16761SJacob Faibussowitsch       if (petscindebugger) { \
8323ba16761SJacob Faibussowitsch         abort(); \
8333ba16761SJacob Faibussowitsch       } else { \
8347da6b3ccSLisandro Dalcin         PetscMPIInt size_; \
8353ba16761SJacob Faibussowitsch         ierr_petsc_abort_ = __VA_ARGS__; \
8367da6b3ccSLisandro Dalcin         MPI_Comm_size(comm, &size_); \
8377da6b3ccSLisandro Dalcin         if (PetscCIEnabledPortableErrorOutput && size_ == PetscGlobalSize && ierr_petsc_abort_ != PETSC_ERR_SIG) { \
8389371c9d4SSatish Balay           MPI_Finalize(); \
8399371c9d4SSatish Balay           exit(0); \
840660278c0SBarry Smith         } else if (PetscCIEnabledPortableErrorOutput && PetscGlobalSize == 1) { \
841660278c0SBarry Smith           exit(0); \
842660278c0SBarry Smith         } else { \
843660278c0SBarry Smith           MPI_Abort(comm, (PetscMPIInt)ierr_petsc_abort_); \
844660278c0SBarry Smith         } \
8453fcd9f07SJacob Faibussowitsch       } \
8467c66cc67SJunchao Zhang     } while (0)
84729f88feaSJacob Faibussowitsch #endif
848986eef2eSBarry Smith 
8499566063dSJacob Faibussowitsch #ifdef PETSC_CLANGUAGE_CXX
850986eef2eSBarry Smith   /*MC
8519566063dSJacob Faibussowitsch   PetscCallThrow - Checks error code, if non-zero it calls the C++ error handler which throws
8529566063dSJacob Faibussowitsch   an exception
853986eef2eSBarry Smith 
854986eef2eSBarry Smith   Synopsis:
8559566063dSJacob Faibussowitsch   #include <petscerror.h>
8569566063dSJacob Faibussowitsch   void PetscCallThrow(PetscErrorCode ierr)
857986eef2eSBarry Smith 
858986eef2eSBarry Smith   Not Collective
859986eef2eSBarry Smith 
8609566063dSJacob Faibussowitsch   Input Parameter:
861986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
862986eef2eSBarry Smith 
863667f096bSBarry Smith   Level: beginner
864667f096bSBarry Smith 
865986eef2eSBarry Smith   Notes:
8669566063dSJacob Faibussowitsch   Requires PETSc to be configured with clanguage = c++. Throws a std::runtime_error() on error.
867986eef2eSBarry Smith 
86887497f52SBarry Smith   Once the error handler throws the exception you can use `PetscCallVoid()` which returns without
86987497f52SBarry Smith   an error code (bad idea since the error is ignored) or `PetscCallAbort()` to have `MPI_Abort()`
8709566063dSJacob Faibussowitsch   called immediately.
8719566063dSJacob Faibussowitsch 
872db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscCall()`, `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`,
873db781477SPatrick Sanan           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`
874986eef2eSBarry Smith M*/
8759371c9d4SSatish Balay   #define PetscCallThrow(...) \
8769371c9d4SSatish Balay     do { \
8773ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
8783ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_throw_ = __VA_ARGS__; \
8793ba16761SJacob 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); \
880ffc4695bSBarry Smith     } while (0)
88185614651SBarry Smith 
882cc26af49SMatthew Knepley   /*MC
883cc26af49SMatthew Knepley   CHKERRXX - Checks error code, if non-zero it calls the C++ error handler which throws an exception
884cc26af49SMatthew Knepley 
885cc26af49SMatthew Knepley   Synopsis:
8869566063dSJacob Faibussowitsch   #include <petscerror.h>
8873af045c5SBarry Smith   void CHKERRXX(PetscErrorCode ierr)
888cc26af49SMatthew Knepley 
889eca87e8dSBarry Smith   Not Collective
890cc26af49SMatthew Knepley 
8919566063dSJacob Faibussowitsch   Input Parameter:
8923af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
893cc26af49SMatthew Knepley 
8949566063dSJacob Faibussowitsch   Level: deprecated
895cc26af49SMatthew Knepley 
896667f096bSBarry Smith   Note:
897667f096bSBarry Smith   Deprecated in favor of `PetscCallThrow()`. This routine behaves identically to it.
898667f096bSBarry Smith 
899db781477SPatrick Sanan .seealso: `PetscCallThrow()`
900cc26af49SMatthew Knepley M*/
9019566063dSJacob Faibussowitsch   #define CHKERRXX(...) PetscCallThrow(__VA_ARGS__)
902fd705b32SMatthew Knepley #endif
903fd705b32SMatthew Knepley 
9043ba16761SJacob Faibussowitsch #define PetscCallCXX_Private(__SETERR_FUNC__, __COMM__, ...) \
9053ba16761SJacob Faibussowitsch   do { \
9063ba16761SJacob Faibussowitsch     PetscStackUpdateLine; \
9073ba16761SJacob Faibussowitsch     try { \
9083ba16761SJacob Faibussowitsch       __VA_ARGS__; \
9093ba16761SJacob Faibussowitsch     } catch (const std::exception &e) { \
9103ba16761SJacob Faibussowitsch       __SETERR_FUNC__(__COMM__, PETSC_ERR_LIB, "%s", e.what()); \
9113ba16761SJacob Faibussowitsch     } \
9123ba16761SJacob Faibussowitsch   } while (0)
9133ba16761SJacob Faibussowitsch 
9143f520e80SJunchao Zhang /*MC
9159566063dSJacob Faibussowitsch   PetscCallCXX - Checks C++ function calls and if they throw an exception, catch it and then
9169566063dSJacob Faibussowitsch   return a PETSc error code
9173f520e80SJunchao Zhang 
9183f520e80SJunchao Zhang   Synopsis:
9199566063dSJacob Faibussowitsch   #include <petscerror.h>
9205687f061SJacob Faibussowitsch   void PetscCallCXX(...) noexcept;
9213f520e80SJunchao Zhang 
9223f520e80SJunchao Zhang   Not Collective
9233f520e80SJunchao Zhang 
9249566063dSJacob Faibussowitsch   Input Parameter:
9255687f061SJacob Faibussowitsch . __VA_ARGS__ - An arbitrary expression
9265687f061SJacob Faibussowitsch 
9275687f061SJacob Faibussowitsch   Level: beginner
9283f520e80SJunchao Zhang 
9293f520e80SJunchao Zhang   Notes:
9305687f061SJacob Faibussowitsch   `PetscCallCXX(...)` is a macro replacement for
9319566063dSJacob Faibussowitsch .vb
9329566063dSJacob Faibussowitsch   try {
9335687f061SJacob Faibussowitsch     __VA_ARGS__;
9349566063dSJacob Faibussowitsch   } catch (const std::exception& e) {
9359566063dSJacob Faibussowitsch     return ConvertToPetscErrorCode(e);
9369566063dSJacob Faibussowitsch   }
9379566063dSJacob Faibussowitsch .ve
9389566063dSJacob Faibussowitsch   Due to the fact that it catches any (reasonable) exception, it is essentially noexcept.
9393f520e80SJunchao Zhang 
9405687f061SJacob Faibussowitsch   If you cannot return a `PetscErrorCode` use `PetscCallCXXAbort()` instead.
9415687f061SJacob Faibussowitsch 
9429566063dSJacob Faibussowitsch   Example Usage:
9439566063dSJacob Faibussowitsch .vb
9449566063dSJacob Faibussowitsch   void foo(void) { throw std::runtime_error("error"); }
9453f520e80SJunchao Zhang 
9469566063dSJacob Faibussowitsch   void bar()
9479566063dSJacob Faibussowitsch   {
948e8952933SJacob Faibussowitsch     PetscCallCXX(foo()); // ERROR bar() does not return PetscErrorCode
9499566063dSJacob Faibussowitsch   }
9509566063dSJacob Faibussowitsch 
9519566063dSJacob Faibussowitsch   PetscErrorCode baz()
9529566063dSJacob Faibussowitsch   {
9539566063dSJacob Faibussowitsch     PetscCallCXX(foo()); // OK
9549566063dSJacob Faibussowitsch 
9559566063dSJacob Faibussowitsch     PetscCallCXX(
9569566063dSJacob Faibussowitsch       bar();
957d5b43468SJose E. Roman       foo(); // OK multiple statements allowed
9589566063dSJacob Faibussowitsch     );
9599566063dSJacob Faibussowitsch   }
9609566063dSJacob Faibussowitsch 
9619566063dSJacob Faibussowitsch   struct bop
9629566063dSJacob Faibussowitsch   {
9639566063dSJacob Faibussowitsch     bop()
9649566063dSJacob Faibussowitsch     {
965e8952933SJacob Faibussowitsch       PetscCallCXX(foo()); // ERROR returns PetscErrorCode, cannot be used in constructors
9669566063dSJacob Faibussowitsch     }
9679566063dSJacob Faibussowitsch   };
9689566063dSJacob Faibussowitsch 
969e8952933SJacob Faibussowitsch   // ERROR contains do-while, cannot be used as function-try block
9709566063dSJacob Faibussowitsch   PetscErrorCode qux() PetscCallCXX(
9719566063dSJacob Faibussowitsch     bar();
9729566063dSJacob Faibussowitsch     baz();
9739566063dSJacob Faibussowitsch     foo();
9749566063dSJacob Faibussowitsch     return 0;
9759566063dSJacob Faibussowitsch   )
9769566063dSJacob Faibussowitsch .ve
9779566063dSJacob Faibussowitsch 
9785687f061SJacob Faibussowitsch .seealso: `PetscCallCXXAbort()`, `PetscCallThrow()`, `SETERRQ()`, `PetscCall()`,
9795687f061SJacob Faibussowitsch           `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
9805687f061SJacob Faibussowitsch           `PetscError()`, `CHKMEMQ`
9813f520e80SJunchao Zhang M*/
9823ba16761SJacob Faibussowitsch #define PetscCallCXX(...) PetscCallCXX_Private(SETERRQ, PETSC_COMM_SELF, __VA_ARGS__)
9835687f061SJacob Faibussowitsch 
9845687f061SJacob Faibussowitsch /*MC
9855687f061SJacob Faibussowitsch   PetscCallCXXAbort - Like `PetscCallCXX()` but calls `MPI_Abort()` instead of returning an
9865687f061SJacob Faibussowitsch   error-code
9875687f061SJacob Faibussowitsch 
9885687f061SJacob Faibussowitsch   Synopsis:
9895687f061SJacob Faibussowitsch   #include <petscerror.h>
9905687f061SJacob Faibussowitsch   void PetscCallCXXAbort(MPI_Comm comm, ...) noexcept;
9915687f061SJacob Faibussowitsch 
99220f4b53cSBarry Smith   Collective; No Fortran Support
9935687f061SJacob Faibussowitsch 
9945687f061SJacob Faibussowitsch   Input Parameters:
9955687f061SJacob Faibussowitsch + comm        - The MPI communicator to abort on
9965687f061SJacob Faibussowitsch - __VA_ARGS__ - An arbitrary expression
9975687f061SJacob Faibussowitsch 
9985687f061SJacob Faibussowitsch   Level: beginner
9995687f061SJacob Faibussowitsch 
10005687f061SJacob Faibussowitsch   Notes:
10015687f061SJacob Faibussowitsch   This macro may be used to check C++ expressions for exceptions in cases where you cannot
10025687f061SJacob Faibussowitsch   return an error code. This includes constructors, destructors, copy/move assignment functions
10035687f061SJacob Faibussowitsch   or constructors among others.
10045687f061SJacob Faibussowitsch 
10055687f061SJacob Faibussowitsch   If an exception is caught, the macro calls `SETERRABORT()` on `comm`. The exception must
10065687f061SJacob Faibussowitsch   derive from `std::exception` in order to be caught.
10075687f061SJacob Faibussowitsch 
10085687f061SJacob Faibussowitsch   If the routine _can_ return an error-code it is highly advised to use `PetscCallCXX()`
10095687f061SJacob Faibussowitsch   instead.
10105687f061SJacob Faibussowitsch 
10115687f061SJacob Faibussowitsch   See `PetscCallCXX()` for additional discussion.
10125687f061SJacob Faibussowitsch 
10135687f061SJacob Faibussowitsch   Example Usage:
10145687f061SJacob Faibussowitsch .vb
10155687f061SJacob Faibussowitsch   class Foo
10165687f061SJacob Faibussowitsch   {
10175687f061SJacob Faibussowitsch     std::vector<int> data_;
10185687f061SJacob Faibussowitsch 
10195687f061SJacob Faibussowitsch   public:
10205687f061SJacob Faibussowitsch     // normally std::vector::reserve() may raise an exception, but since we handle it with
10215687f061SJacob Faibussowitsch     // PetscCallCXXAbort() we may mark this routine as noexcept!
10225687f061SJacob Faibussowitsch     Foo() noexcept
10235687f061SJacob Faibussowitsch     {
10245687f061SJacob Faibussowitsch       PetscCallCXXAbort(PETSC_COMM_SELF, data_.reserve(10));
10255687f061SJacob Faibussowitsch     }
10265687f061SJacob Faibussowitsch   };
10275687f061SJacob Faibussowitsch 
10285687f061SJacob Faibussowitsch   std::vector<int> bar()
10295687f061SJacob Faibussowitsch   {
10305687f061SJacob Faibussowitsch     std::vector<int> v;
10315687f061SJacob Faibussowitsch 
10325687f061SJacob Faibussowitsch     PetscFunctionBegin;
10335687f061SJacob Faibussowitsch     // OK!
10345687f061SJacob Faibussowitsch     PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1));
10355687f061SJacob Faibussowitsch     PetscFunctionReturn(v);
10365687f061SJacob Faibussowitsch   }
10375687f061SJacob Faibussowitsch 
10385687f061SJacob Faibussowitsch   PetscErrorCode baz()
10395687f061SJacob Faibussowitsch   {
10405687f061SJacob Faibussowitsch     std::vector<int> v;
10415687f061SJacob Faibussowitsch 
10425687f061SJacob Faibussowitsch     PetscFunctionBegin;
10435687f061SJacob Faibussowitsch     // WRONG! baz() returns a PetscErrorCode, prefer PetscCallCXX() instead
10445687f061SJacob Faibussowitsch     PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1));
10453ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
10465687f061SJacob Faibussowitsch   }
10475687f061SJacob Faibussowitsch .ve
10485687f061SJacob Faibussowitsch 
10495687f061SJacob Faibussowitsch .seealso: `PetscCallCXX()`, `SETERRABORT()`, `PetscCallAbort()`
10505687f061SJacob Faibussowitsch M*/
10513ba16761SJacob Faibussowitsch #define PetscCallCXXAbort(comm, ...) PetscCallCXX_Private(SETERRABORT, comm, __VA_ARGS__)
10523f520e80SJunchao Zhang 
105330de9b25SBarry Smith /*MC
10549566063dSJacob Faibussowitsch   CHKERRCXX - Checks C++ function calls and if they throw an exception, catch it and then
10559566063dSJacob Faibussowitsch   return a PETSc error code
10569566063dSJacob Faibussowitsch 
10579566063dSJacob Faibussowitsch   Synopsis:
10589566063dSJacob Faibussowitsch   #include <petscerror.h>
10599566063dSJacob Faibussowitsch   void CHKERRCXX(func) noexcept;
10609566063dSJacob Faibussowitsch 
10619566063dSJacob Faibussowitsch   Not Collective
10629566063dSJacob Faibussowitsch 
10639566063dSJacob Faibussowitsch   Input Parameter:
10649566063dSJacob Faibussowitsch . func - C++ function calls
10659566063dSJacob Faibussowitsch 
10669566063dSJacob Faibussowitsch   Level: deprecated
10679566063dSJacob Faibussowitsch 
1068667f096bSBarry Smith   Note:
1069667f096bSBarry Smith   Deprecated in favor of `PetscCallCXX()`. This routine behaves identically to it.
1070667f096bSBarry Smith 
1071db781477SPatrick Sanan .seealso: `PetscCallCXX()`
10729566063dSJacob Faibussowitsch M*/
10739566063dSJacob Faibussowitsch #define CHKERRCXX(...) PetscCallCXX(__VA_ARGS__)
10749566063dSJacob Faibussowitsch 
10759566063dSJacob Faibussowitsch /*MC
107630de9b25SBarry Smith    CHKMEMQ - Checks the memory for corruption, calls error handler if any is detected
107730de9b25SBarry Smith 
107830de9b25SBarry Smith    Synopsis:
1079aaa7dc30SBarry Smith    #include <petscsys.h>
108091d3bdf4SKris Buschelman    CHKMEMQ;
108130de9b25SBarry Smith 
1082eca87e8dSBarry Smith    Not Collective
1083eca87e8dSBarry Smith 
108430de9b25SBarry Smith   Level: beginner
108530de9b25SBarry Smith 
108630de9b25SBarry Smith    Notes:
1087a17b96a8SKyle Gerard Felker     We highly recommend using Valgrind https://petsc.org/release/faq/#valgrind or for NVIDIA CUDA systems
10885ed36255SBarry Smith     https://docs.nvidia.com/cuda/cuda-memcheck/index.html for finding memory problems. The ``CHKMEMQ`` macro is useful on systems that
10895ed36255SBarry Smith     do not have valgrind, but is not as good as valgrind or cuda-memcheck.
10901957e957SBarry Smith 
1091667f096bSBarry Smith     Must run with the option `-malloc_debug` (`-malloc_test` in debug mode; or if `PetscMallocSetDebug()` called) to enable this option
109230de9b25SBarry Smith 
109330de9b25SBarry Smith     Once the error handler is called the calling function is then returned from with the given error code.
109430de9b25SBarry Smith 
109530de9b25SBarry Smith     By defaults prints location where memory that is corrupted was allocated.
109630de9b25SBarry Smith 
109787497f52SBarry Smith     Use `CHKMEMA` for functions that return void
1098f621e05eSBarry Smith 
1099db781477SPatrick Sanan .seealso: `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `SETERRQ()`, `PetscMallocValidate()`
110030de9b25SBarry Smith M*/
11016d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
1102064a246eSJacob Faibussowitsch   #define CHKMEMQ
1103064a246eSJacob Faibussowitsch   #define CHKMEMA
11046d210af2SJacob Faibussowitsch #else
11059371c9d4SSatish Balay   #define CHKMEMQ \
11069371c9d4SSatish Balay     do { \
11073ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_memq_ = PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__); \
11083ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_memq_ != PETSC_SUCCESS)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_memq_, PETSC_ERROR_REPEAT, " "); \
110986d09637SLisandro Dalcin     } while (0)
11106d210af2SJacob Faibussowitsch   #define CHKMEMA PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__)
1111064a246eSJacob Faibussowitsch #endif
11129566063dSJacob Faibussowitsch 
1113668f157eSBarry Smith /*E
1114668f157eSBarry Smith   PetscErrorType - passed to the PETSc error handling routines indicating if this is the first or a later call to the error handlers
1115668f157eSBarry Smith 
1116668f157eSBarry Smith   Level: advanced
1117668f157eSBarry Smith 
1118667f096bSBarry Smith   Note:
111987497f52SBarry Smith   `PETSC_ERROR_IN_CXX` indicates the error was detected in C++ and an exception should be generated
1120d736bfebSBarry Smith 
1121667f096bSBarry Smith   Developer Note:
1122667f096bSBarry Smith     This is currently used to decide when to print the detailed information about the run in `PetscTraceBackErrorHandler()`
1123668f157eSBarry Smith 
112487497f52SBarry Smith .seealso: `PetscError()`, `SETERRQ()`
1125668f157eSBarry Smith E*/
11269371c9d4SSatish Balay typedef enum {
11279371c9d4SSatish Balay   PETSC_ERROR_INITIAL = 0,
11289371c9d4SSatish Balay   PETSC_ERROR_REPEAT  = 1,
11299371c9d4SSatish Balay   PETSC_ERROR_IN_CXX  = 2
11309371c9d4SSatish Balay } PetscErrorType;
11314b209cf6SBarry Smith 
1132eb9e708aSLisandro Dalcin #if defined(__clang_analyzer__)
1133eb9e708aSLisandro Dalcin __attribute__((analyzer_noreturn))
1134eb9e708aSLisandro Dalcin #endif
11359371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode
11369371c9d4SSatish Balay PetscError(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, ...) PETSC_ATTRIBUTE_COLD PETSC_ATTRIBUTE_FORMAT(7, 8);
1137eb9e708aSLisandro Dalcin 
1138014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscErrorPrintfInitialize(void);
11393ba16761SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscErrorMessage(PetscErrorCode, const char *[], char **);
1140d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscTraceBackErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1141d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscIgnoreErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1142d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscEmacsClientErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1143d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscMPIAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1144d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1145d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAttachDebuggerErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1146d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscReturnErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1147efca3c55SSatish Balay PETSC_EXTERN PetscErrorCode PetscPushErrorHandler(PetscErrorCode (*handler)(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *), void *);
1148014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopErrorHandler(void);
11498d359177SBarry Smith PETSC_EXTERN PetscErrorCode PetscSignalHandlerDefault(int, void *);
1150014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPushSignalHandler(PetscErrorCode (*)(int, void *), void *);
1151014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopSignalHandler(void);
115228559dc8SJed Brown PETSC_EXTERN PetscErrorCode PetscCheckPointerSetIntensity(PetscInt);
1153c2a741eeSJunchao Zhang PETSC_EXTERN void           PetscSignalSegvCheckPointerOrMpi(void);
1154d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use PetscSignalSegvCheckPointerOrMpi() (since version 3.13)") static inline void PetscSignalSegvCheckPointer(void)
1155d71ae5a4SJacob Faibussowitsch {
11569371c9d4SSatish Balay   PetscSignalSegvCheckPointerOrMpi();
11579371c9d4SSatish Balay }
1158329f5518SBarry Smith 
1159639ff905SBarry Smith /*MC
1160639ff905SBarry Smith     PetscErrorPrintf - Prints error messages.
1161639ff905SBarry Smith 
1162cf53795eSBarry Smith     Not Collective; No Fortran Support
1163cf53795eSBarry Smith 
1164639ff905SBarry Smith    Synopsis:
1165aaa7dc30SBarry Smith     #include <petscsys.h>
1166639ff905SBarry Smith      PetscErrorCode (*PetscErrorPrintf)(const char format[],...);
1167639ff905SBarry Smith 
1168f899ff85SJose E. Roman     Input Parameter:
1169cd05f99aSJacob Faibussowitsch .   format - the usual `printf()` format string
1170639ff905SBarry Smith 
1171639ff905SBarry Smith    Options Database Keys:
11721957e957SBarry Smith +    -error_output_stdout - cause error messages to be printed to stdout instead of the (default) stderr
1173e1bc860dSBarry Smith -    -error_output_none - to turn off all printing of error messages (does not change the way the error is handled.)
1174639ff905SBarry Smith 
1175cf53795eSBarry Smith    Level: developer
1176cf53795eSBarry Smith 
117795452b02SPatrick Sanan    Notes:
117895452b02SPatrick Sanan     Use
1179667f096bSBarry Smith .vb
1180667f096bSBarry Smith      PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the
1181667f096bSBarry Smith                         error is handled.) and
1182667f096bSBarry Smith      PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on or you can use your own function
1183667f096bSBarry Smith .ve
1184639ff905SBarry Smith      Use
1185667f096bSBarry Smith .vb
118687497f52SBarry Smith      `PETSC_STDERR` = FILE* obtained from a file open etc. to have stderr printed to the file.
118787497f52SBarry Smith      `PETSC_STDOUT` = FILE* obtained from a file open etc. to have stdout printed to the file.
1188667f096bSBarry Smith .ve
1189639ff905SBarry Smith 
1190639ff905SBarry Smith        Use
119187497f52SBarry Smith       `PetscPushErrorHandler()` to provide your own error handler that determines what kind of messages to print
1192639ff905SBarry Smith 
1193db781477SPatrick Sanan .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscHelpPrintf()`, `PetscPrintf()`, `PetscPushErrorHandler()`, `PetscVFPrintf()`, `PetscHelpPrintf()`
1194639ff905SBarry Smith M*/
11953ca90d2dSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode (*PetscErrorPrintf)(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1196639ff905SBarry Smith 
1197cf0818bdSBarry Smith /*E
1198cf0818bdSBarry Smith      PetscFPTrap - types of floating point exceptions that may be trapped
1199cf0818bdSBarry Smith 
120087497f52SBarry Smith      Currently only `PETSC_FP_TRAP_OFF` and `PETSC_FP_TRAP_ON` are handled. All others are treated as `PETSC_FP_TRAP_ON`.
1201cf0818bdSBarry Smith 
1202cf0818bdSBarry Smith      Level: intermediate
1203cf0818bdSBarry Smith 
1204cf0818bdSBarry Smith .seealso: `PetscSetFPTrap()`, `PetscPushFPTrap()`
1205cf0818bdSBarry Smith  E*/
12069371c9d4SSatish Balay typedef enum {
12079371c9d4SSatish Balay   PETSC_FP_TRAP_OFF      = 0,
12089371c9d4SSatish Balay   PETSC_FP_TRAP_INDIV    = 1,
12099371c9d4SSatish Balay   PETSC_FP_TRAP_FLTOPERR = 2,
12109371c9d4SSatish Balay   PETSC_FP_TRAP_FLTOVF   = 4,
12119371c9d4SSatish Balay   PETSC_FP_TRAP_FLTUND   = 8,
12129371c9d4SSatish Balay   PETSC_FP_TRAP_FLTDIV   = 16,
12139371c9d4SSatish Balay   PETSC_FP_TRAP_FLTINEX  = 32
12149371c9d4SSatish Balay } PetscFPTrap;
1215bd2b07b1SBarry 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)
1216014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSetFPTrap(PetscFPTrap);
1217014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPush(PetscFPTrap);
1218014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPop(void);
1219aba4c478SBarry Smith PETSC_EXTERN PetscErrorCode PetscDetermineInitialFPTrap(void);
122054a8ef01SBarry Smith 
12213a40ed3dSBarry Smith /*
12223a40ed3dSBarry Smith       Allows the code to build a stack frame as it runs
12233a40ed3dSBarry Smith */
12243a40ed3dSBarry Smith 
122599cd645aSJed Brown #define PETSCSTACKSIZE 64
12263a40ed3dSBarry Smith typedef struct {
12270e33f6ddSBarry Smith   const char *function[PETSCSTACKSIZE];
12280e33f6ddSBarry Smith   const char *file[PETSCSTACKSIZE];
1229184914b5SBarry Smith   int         line[PETSCSTACKSIZE];
1230362febeeSStefano Zampini   int         petscroutine[PETSCSTACKSIZE]; /* 0 external called from petsc, 1 petsc functions, 2 petsc user functions */
1231184914b5SBarry Smith   int         currentsize;
1232a2f94806SJed Brown   int         hotdepth;
12334be741a6SBarry Smith   PetscBool   check; /* option to check for correct Push/Pop semantics, true for default petscstack but not other stacks */
12343a40ed3dSBarry Smith } PetscStack;
1235dfb7d7afSStefano Zampini #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
123627104ee2SJacob Faibussowitsch PETSC_EXTERN PetscStack petscstack;
123727104ee2SJacob Faibussowitsch #endif
12383a40ed3dSBarry Smith 
12395d12eec7SSatish Balay #if defined(PETSC_SERIALIZE_FUNCTIONS)
12405d12eec7SSatish Balay   #include <petsc/private/petscfptimpl.h>
12415d12eec7SSatish Balay   /*
12425d12eec7SSatish Balay    Registers the current function into the global function pointer to function name table
12435d12eec7SSatish Balay 
12445d12eec7SSatish Balay    Have to fix this to handle errors but cannot return error since used in PETSC_VIEWER_DRAW_() etc
12455d12eec7SSatish Balay */
12469371c9d4SSatish Balay   #define PetscRegister__FUNCT__() \
12479371c9d4SSatish Balay     do { \
12485d12eec7SSatish Balay       static PetscBool __chked = PETSC_FALSE; \
12495d12eec7SSatish Balay       if (!__chked) { \
12509371c9d4SSatish Balay         void *ptr; \
12513ba16761SJacob Faibussowitsch         PetscCallAbort(PETSC_COMM_SELF, PetscDLSym(NULL, PETSC_FUNCTION_NAME, &ptr)); \
12525d12eec7SSatish Balay         __chked = PETSC_TRUE; \
12539371c9d4SSatish Balay       } \
12549371c9d4SSatish Balay     } while (0)
12555d12eec7SSatish Balay #else
12565d12eec7SSatish Balay   #define PetscRegister__FUNCT__()
12575d12eec7SSatish Balay #endif
12585d12eec7SSatish Balay 
1259eae3dc7dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) || defined(__clang_analyzer__)
12606d210af2SJacob Faibussowitsch   #define PetscStackPushNoCheck(funct, petsc_routine, hot)
126137154d10SBarry Smith   #define PetscStackUpdateLine
1262792fecdfSBarry Smith   #define PetscStackPushExternal(funct)
12636d210af2SJacob Faibussowitsch   #define PetscStackPopNoCheck
12646d210af2SJacob Faibussowitsch   #define PetscStackClearTop
12656d210af2SJacob Faibussowitsch   #define PetscFunctionBegin
12666d210af2SJacob Faibussowitsch   #define PetscFunctionBeginUser
12676d210af2SJacob Faibussowitsch   #define PetscFunctionBeginHot
12680a57284eSJacob Faibussowitsch   #define PetscFunctionReturn(...)  return __VA_ARGS__
12696d210af2SJacob Faibussowitsch   #define PetscFunctionReturnVoid() return
12706d210af2SJacob Faibussowitsch   #define PetscStackPop
12716d210af2SJacob Faibussowitsch   #define PetscStackPush(f)
1272dfb7d7afSStefano Zampini #elif defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1273660278c0SBarry Smith 
12749371c9d4SSatish Balay   #define PetscStackPush_Private(stack__, file__, func__, line__, petsc_routine__, hot__) \
12759371c9d4SSatish Balay     do { \
12765a96b57dSJacob Faibussowitsch       if (stack__.currentsize < PETSCSTACKSIZE) { \
12775a96b57dSJacob Faibussowitsch         stack__.function[stack__.currentsize] = func__; \
1278ef1023bdSBarry Smith         if (petsc_routine__) { \
1279ef1023bdSBarry Smith           stack__.file[stack__.currentsize] = file__; \
12805a96b57dSJacob Faibussowitsch           stack__.line[stack__.currentsize] = line__; \
1281ef1023bdSBarry Smith         } else { \
1282648bc8c4SBarry Smith           stack__.file[stack__.currentsize] = PETSC_NULLPTR; \
1283ef1023bdSBarry Smith           stack__.line[stack__.currentsize] = 0; \
1284ef1023bdSBarry Smith         } \
12855a96b57dSJacob Faibussowitsch         stack__.petscroutine[stack__.currentsize] = petsc_routine__; \
12865a96b57dSJacob Faibussowitsch       } \
12875a96b57dSJacob Faibussowitsch       ++stack__.currentsize; \
12885a96b57dSJacob Faibussowitsch       stack__.hotdepth += (hot__ || stack__.hotdepth); \
12895a96b57dSJacob Faibussowitsch     } while (0)
12905a96b57dSJacob Faibussowitsch 
12914be741a6SBarry Smith   /* uses PetscCheckAbort() because may be used in a function that does not return an error code */
12929371c9d4SSatish Balay   #define PetscStackPop_Private(stack__, func__) \
12939371c9d4SSatish Balay     do { \
12944be741a6SBarry 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__); \
12955a96b57dSJacob Faibussowitsch       if (--stack__.currentsize < PETSCSTACKSIZE) { \
12969371c9d4SSatish 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", \
12979371c9d4SSatish Balay                         stack__.function[stack__.currentsize], stack__.file[stack__.currentsize], stack__.line[stack__.currentsize], func__, __FILE__, __LINE__); \
12985a96b57dSJacob Faibussowitsch         stack__.function[stack__.currentsize]     = PETSC_NULLPTR; \
12995a96b57dSJacob Faibussowitsch         stack__.file[stack__.currentsize]         = PETSC_NULLPTR; \
13005a96b57dSJacob Faibussowitsch         stack__.line[stack__.currentsize]         = 0; \
13015a96b57dSJacob Faibussowitsch         stack__.petscroutine[stack__.currentsize] = 0; \
13025a96b57dSJacob Faibussowitsch       } \
13035a96b57dSJacob Faibussowitsch       stack__.hotdepth = PetscMax(stack__.hotdepth - 1, 0); \
13045a96b57dSJacob Faibussowitsch     } while (0)
13055a96b57dSJacob Faibussowitsch 
1306586f9135SBarry Smith   /*MC
1307586f9135SBarry Smith    PetscStackPushNoCheck - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is
1308586f9135SBarry Smith    currently in the source code.
1309586f9135SBarry Smith 
1310586f9135SBarry Smith    Not Collective
1311586f9135SBarry Smith 
1312586f9135SBarry Smith    Synopsis:
1313586f9135SBarry Smith    #include <petscsys.h>
1314586f9135SBarry Smith    void PetscStackPushNoCheck(char *funct,int petsc_routine,PetscBool hot);
1315586f9135SBarry Smith 
1316586f9135SBarry Smith    Input Parameters:
1317586f9135SBarry Smith +  funct - the function name
1318586f9135SBarry Smith .  petsc_routine - 2 user function, 1 PETSc function, 0 some other function
1319586f9135SBarry Smith -  hot - indicates that the function may be called often so expensive error checking should be turned off inside the function
1320586f9135SBarry Smith 
1321586f9135SBarry Smith    Level: developer
1322586f9135SBarry Smith 
1323586f9135SBarry Smith    Notes:
1324586f9135SBarry 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
132587497f52SBarry 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
1326586f9135SBarry Smith    help debug the problem.
1327586f9135SBarry Smith 
1328ef1023bdSBarry Smith    This version does not check the memory corruption (an expensive operation), use `PetscStackPush()` to check the memory.
1329ef1023bdSBarry Smith 
1330792fecdfSBarry Smith    Use `PetscStackPushExternal()` for a function call that is about to be made to a non-PETSc or user function (such as BLAS etc).
1331ef1023bdSBarry Smith 
1332586f9135SBarry Smith    The default stack is a global variable called `petscstack`.
1333586f9135SBarry Smith 
1334586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1335ef1023bdSBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPush()`, `PetscStackPop`,
1336792fecdfSBarry Smith           `PetscStackPushExternal()`
1337586f9135SBarry Smith M*/
13389371c9d4SSatish Balay   #define PetscStackPushNoCheck(funct, petsc_routine, hot) \
13399371c9d4SSatish Balay     do { \
1340e04113cfSBarry Smith       PetscStackSAWsTakeAccess(); \
13415a96b57dSJacob Faibussowitsch       PetscStackPush_Private(petscstack, __FILE__, funct, __LINE__, petsc_routine, hot); \
1342e04113cfSBarry Smith       PetscStackSAWsGrantAccess(); \
1343441dd030SJed Brown     } while (0)
1344441dd030SJed Brown 
1345586f9135SBarry Smith   /*MC
134687497f52SBarry Smith    PetscStackUpdateLine - in a function that has a `PetscFunctionBegin` or `PetscFunctionBeginUser` updates the stack line number to the
134737154d10SBarry Smith    current line number.
134837154d10SBarry Smith 
134937154d10SBarry Smith    Not Collective
135037154d10SBarry Smith 
135137154d10SBarry Smith    Synopsis:
135237154d10SBarry Smith    #include <petscsys.h>
135337154d10SBarry Smith    void PetscStackUpdateLine
135437154d10SBarry Smith 
135537154d10SBarry Smith    Level: developer
135637154d10SBarry Smith 
135737154d10SBarry Smith    Notes:
135887497f52SBarry Smith    Using `PetscCall()` and friends automatically handles this process
135987497f52SBarry Smith 
136037154d10SBarry 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
136137154d10SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
136237154d10SBarry Smith    help debug the problem.
136337154d10SBarry Smith 
136437154d10SBarry Smith    The default stack is a global variable called petscstack.
136537154d10SBarry Smith 
136637154d10SBarry Smith    This is used by `PetscCall()` and is otherwise not like to be needed
136737154d10SBarry Smith 
136837154d10SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`, `PetscCall()`
136937154d10SBarry Smith M*/
137037154d10SBarry Smith   #define PetscStackUpdateLine \
13713ba16761SJacob Faibussowitsch     do { \
13723ba16761SJacob Faibussowitsch       if (petscstack.currentsize > 0 && petscstack.function[petscstack.currentsize - 1] == PETSC_FUNCTION_NAME) { petscstack.line[petscstack.currentsize - 1] = __LINE__; } \
13733ba16761SJacob Faibussowitsch     } while (0)
137437154d10SBarry Smith 
137537154d10SBarry Smith   /*MC
1376792fecdfSBarry Smith    PetscStackPushExternal - Pushes a new function name onto the PETSc default stack that tracks where the running program is
1377ef1023bdSBarry Smith    currently in the source code. Does not include the filename or line number since this is called by the calling routine
1378ef1023bdSBarry Smith    for non-PETSc or user functions.
1379ef1023bdSBarry Smith 
1380ef1023bdSBarry Smith    Not Collective
1381ef1023bdSBarry Smith 
1382ef1023bdSBarry Smith    Synopsis:
1383ef1023bdSBarry Smith    #include <petscsys.h>
1384792fecdfSBarry Smith    void PetscStackPushExternal(char *funct);
1385ef1023bdSBarry Smith 
13862fe279fdSBarry Smith    Input Parameter:
1387ef1023bdSBarry Smith .  funct - the function name
1388ef1023bdSBarry Smith 
1389ef1023bdSBarry Smith    Level: developer
1390ef1023bdSBarry Smith 
1391ef1023bdSBarry Smith    Notes:
139287497f52SBarry Smith    Using `PetscCallExternal()` and friends automatically handles this process
139387497f52SBarry Smith 
1394ef1023bdSBarry 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
1395ef1023bdSBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1396ef1023bdSBarry Smith    help debug the problem.
1397ef1023bdSBarry Smith 
1398ef1023bdSBarry Smith    The default stack is a global variable called `petscstack`.
1399ef1023bdSBarry Smith 
1400ef1023bdSBarry Smith    This is to be used when calling an external package function such as a BLAS function.
1401ef1023bdSBarry Smith 
1402ef1023bdSBarry Smith    This also updates the stack line number for the current stack function.
1403ef1023bdSBarry Smith 
1404ef1023bdSBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1405ef1023bdSBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1406ef1023bdSBarry Smith M*/
14079371c9d4SSatish Balay   #define PetscStackPushExternal(funct) \
14089371c9d4SSatish Balay     do { \
14099371c9d4SSatish Balay       PetscStackUpdateLine; \
14109371c9d4SSatish Balay       PetscStackPushNoCheck(funct, 0, PETSC_TRUE); \
14119371c9d4SSatish Balay     } while (0);
1412ef1023bdSBarry Smith 
1413ef1023bdSBarry Smith   /*MC
1414586f9135SBarry Smith    PetscStackPopNoCheck - Pops a function name from the PETSc default stack that tracks where the running program is
1415586f9135SBarry Smith    currently in the source code.
1416586f9135SBarry Smith 
1417586f9135SBarry Smith    Not Collective
1418586f9135SBarry Smith 
1419586f9135SBarry Smith    Synopsis:
1420586f9135SBarry Smith    #include <petscsys.h>
1421586f9135SBarry Smith    void PetscStackPopNoCheck(char *funct);
1422586f9135SBarry Smith 
1423586f9135SBarry Smith    Input Parameter:
1424586f9135SBarry Smith .   funct - the function name
1425586f9135SBarry Smith 
1426586f9135SBarry Smith    Level: developer
1427586f9135SBarry Smith 
1428586f9135SBarry Smith    Notes:
142987497f52SBarry Smith    Using `PetscCall()`, `PetscCallExternal()`, `PetscCallBack()` and friends negates the need to call this
143087497f52SBarry Smith 
1431586f9135SBarry 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
1432586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1433586f9135SBarry Smith    help debug the problem.
1434586f9135SBarry Smith 
1435586f9135SBarry Smith    The default stack is a global variable called petscstack.
1436586f9135SBarry Smith 
1437586f9135SBarry Smith    Developer Note:
1438586f9135SBarry Smith    `PetscStackPopNoCheck()` takes a function argument while  `PetscStackPop` does not, this difference is likely just historical.
1439586f9135SBarry Smith 
1440586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1441586f9135SBarry Smith M*/
14429371c9d4SSatish Balay   #define PetscStackPopNoCheck(funct) \
14439371c9d4SSatish Balay     do { \
1444362febeeSStefano Zampini       PetscStackSAWsTakeAccess(); \
14455a96b57dSJacob Faibussowitsch       PetscStackPop_Private(petscstack, funct); \
1446362febeeSStefano Zampini       PetscStackSAWsGrantAccess(); \
1447362febeeSStefano Zampini     } while (0)
1448362febeeSStefano Zampini 
14499371c9d4SSatish Balay   #define PetscStackClearTop \
14509371c9d4SSatish Balay     do { \
1451e04113cfSBarry Smith       PetscStackSAWsTakeAccess(); \
14529371c9d4SSatish Balay       if (petscstack.currentsize > 0 && --petscstack.currentsize < PETSCSTACKSIZE) { \
145327104ee2SJacob Faibussowitsch         petscstack.function[petscstack.currentsize]     = PETSC_NULLPTR; \
145427104ee2SJacob Faibussowitsch         petscstack.file[petscstack.currentsize]         = PETSC_NULLPTR; \
145527104ee2SJacob Faibussowitsch         petscstack.line[petscstack.currentsize]         = 0; \
145627104ee2SJacob Faibussowitsch         petscstack.petscroutine[petscstack.currentsize] = 0; \
1457441dd030SJed Brown       } \
145827104ee2SJacob Faibussowitsch       petscstack.hotdepth = PetscMax(petscstack.hotdepth - 1, 0); \
1459e04113cfSBarry Smith       PetscStackSAWsGrantAccess(); \
1460441dd030SJed Brown     } while (0)
1461441dd030SJed Brown 
146230de9b25SBarry Smith   /*MC
14631957e957SBarry Smith    PetscFunctionBegin - First executable line of each PETSc function,  used for error handling. Final
146487497f52SBarry Smith       line of PETSc functions should be `PetscFunctionReturn`(0);
146530de9b25SBarry Smith 
146630de9b25SBarry Smith    Synopsis:
1467aaa7dc30SBarry Smith    #include <petscsys.h>
146830de9b25SBarry Smith    void PetscFunctionBegin;
146930de9b25SBarry Smith 
147020f4b53cSBarry Smith    Not Collective; No Fortran Support
1471eca87e8dSBarry Smith 
147230de9b25SBarry Smith    Usage:
147330de9b25SBarry Smith .vb
147430de9b25SBarry Smith      int something;
147530de9b25SBarry Smith 
147630de9b25SBarry Smith      PetscFunctionBegin;
147730de9b25SBarry Smith .ve
147830de9b25SBarry Smith 
147930de9b25SBarry Smith    Level: developer
148030de9b25SBarry Smith 
148120f4b53cSBarry Smith    Note:
148220f4b53cSBarry Smith      Use `PetscFunctionBeginUser` for application codes.
148320f4b53cSBarry Smith 
1484586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`
148530de9b25SBarry Smith 
148630de9b25SBarry Smith M*/
14879371c9d4SSatish Balay   #define PetscFunctionBegin \
14889371c9d4SSatish Balay     do { \
1489362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_FALSE); \
1490a2f94806SJed Brown       PetscRegister__FUNCT__(); \
1491a2f94806SJed Brown     } while (0)
1492a2f94806SJed Brown 
1493a2f94806SJed Brown   /*MC
149487497f52SBarry Smith    PetscFunctionBeginHot - Substitute for `PetscFunctionBegin` to be used in functions that are called in
1495a2f94806SJed Brown    performance-critical circumstances.  Use of this function allows for lighter profiling by default.
1496a2f94806SJed Brown 
1497a2f94806SJed Brown    Synopsis:
1498aaa7dc30SBarry Smith    #include <petscsys.h>
1499a2f94806SJed Brown    void PetscFunctionBeginHot;
1500a2f94806SJed Brown 
150120f4b53cSBarry Smith    Not Collective; No Fortran Support
1502a2f94806SJed Brown 
1503a2f94806SJed Brown    Usage:
1504a2f94806SJed Brown .vb
1505a2f94806SJed Brown      int something;
1506a2f94806SJed Brown 
1507a2f94806SJed Brown      PetscFunctionBeginHot;
1508a2f94806SJed Brown .ve
1509a2f94806SJed Brown 
1510a2f94806SJed Brown    Level: developer
1511a2f94806SJed Brown 
1512586f9135SBarry Smith .seealso: `PetscFunctionBegin`, `PetscFunctionReturn()`, `PetscStackPushNoCheck()`
1513a2f94806SJed Brown 
1514a2f94806SJed Brown M*/
15159371c9d4SSatish Balay   #define PetscFunctionBeginHot \
15169371c9d4SSatish Balay     do { \
1517362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_TRUE); \
15182d53ad75SBarry Smith       PetscRegister__FUNCT__(); \
151953c77d0aSJed Brown     } while (0)
152053c77d0aSJed Brown 
1521a8d2bbe5SBarry Smith   /*MC
1522530d85c1SBarry Smith    PetscFunctionBeginUser - First executable line of user provided routines
1523a8d2bbe5SBarry Smith 
1524a8d2bbe5SBarry Smith    Synopsis:
1525aaa7dc30SBarry Smith    #include <petscsys.h>
1526a8d2bbe5SBarry Smith    void PetscFunctionBeginUser;
1527a8d2bbe5SBarry Smith 
152820f4b53cSBarry Smith    Not Collective; No Fortran Support
1529a8d2bbe5SBarry Smith 
1530a8d2bbe5SBarry Smith    Usage:
1531a8d2bbe5SBarry Smith .vb
1532a8d2bbe5SBarry Smith      int something;
1533a8d2bbe5SBarry Smith 
1534ac285190SBarry Smith      PetscFunctionBeginUser;
1535a8d2bbe5SBarry Smith .ve
1536a8d2bbe5SBarry Smith 
153720f4b53cSBarry Smith    Level: intermediate
153820f4b53cSBarry Smith 
1539a8d2bbe5SBarry Smith    Notes:
1540530d85c1SBarry Smith       Functions that incorporate this must call `PetscFunctionReturn()` instead of return except for main().
1541530d85c1SBarry Smith 
1542530d85c1SBarry Smith       May be used before `PetscInitialize()`
15431957e957SBarry Smith 
1544530d85c1SBarry Smith       This is identical to `PetscFunctionBegin` except it labels the routine as a user
1545ac285190SBarry Smith       routine instead of as a PETSc library routine.
1546ac285190SBarry Smith 
1547586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, `PetscFunctionBeginHot`, `PetscStackPushNoCheck()`
1548a8d2bbe5SBarry Smith 
1549a8d2bbe5SBarry Smith M*/
15509371c9d4SSatish Balay   #define PetscFunctionBeginUser \
15519371c9d4SSatish Balay     do { \
1552362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 2, PETSC_FALSE); \
1553a8d2bbe5SBarry Smith       PetscRegister__FUNCT__(); \
1554a8d2bbe5SBarry Smith     } while (0)
1555a8d2bbe5SBarry Smith 
1556586f9135SBarry Smith   /*MC
1557586f9135SBarry Smith    PetscStackPush - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is
1558586f9135SBarry Smith    currently in the source code and verifies the memory is not corrupted.
1559586f9135SBarry Smith 
1560586f9135SBarry Smith    Not Collective
1561586f9135SBarry Smith 
1562586f9135SBarry Smith    Synopsis:
1563586f9135SBarry Smith    #include <petscsys.h>
1564586f9135SBarry Smith    void PetscStackPush(char *funct)
1565586f9135SBarry Smith 
1566586f9135SBarry Smith    Input Parameter:
1567586f9135SBarry Smith .  funct - the function name
1568586f9135SBarry Smith 
1569586f9135SBarry Smith    Level: developer
1570586f9135SBarry Smith 
1571586f9135SBarry Smith    Notes:
1572586f9135SBarry 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
1573586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1574586f9135SBarry Smith    help debug the problem.
1575586f9135SBarry Smith 
1576586f9135SBarry Smith    The default stack is a global variable called petscstack.
1577586f9135SBarry Smith 
1578586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1579586f9135SBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1580586f9135SBarry Smith M*/
15819371c9d4SSatish Balay   #define PetscStackPush(n) \
15829371c9d4SSatish Balay     do { \
1583362febeeSStefano Zampini       PetscStackPushNoCheck(n, 0, PETSC_FALSE); \
158415681b3cSBarry Smith       CHKMEMQ; \
158515681b3cSBarry Smith     } while (0)
15863a40ed3dSBarry Smith 
1587586f9135SBarry Smith   /*MC
1588586f9135SBarry Smith    PetscStackPop - Pops a function name from the PETSc default stack that tracks where the running program is
1589586f9135SBarry Smith    currently in the source code and verifies the memory is not corrupted.
1590586f9135SBarry Smith 
1591586f9135SBarry Smith    Not Collective
1592586f9135SBarry Smith 
1593586f9135SBarry Smith    Synopsis:
1594586f9135SBarry Smith    #include <petscsys.h>
1595586f9135SBarry Smith    void PetscStackPop
1596586f9135SBarry Smith 
1597586f9135SBarry Smith    Level: developer
1598586f9135SBarry Smith 
1599586f9135SBarry Smith    Notes:
1600586f9135SBarry 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
1601586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1602586f9135SBarry Smith    help debug the problem.
1603586f9135SBarry Smith 
1604586f9135SBarry Smith    The default stack is a global variable called petscstack.
1605586f9135SBarry Smith 
1606586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPopNoCheck()`, `PetscStackPush()`
1607586f9135SBarry Smith M*/
16089371c9d4SSatish Balay   #define PetscStackPop \
16099371c9d4SSatish Balay     do { \
1610441dd030SJed Brown       CHKMEMQ; \
1611362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
161215681b3cSBarry Smith     } while (0)
1613d64ed03dSBarry Smith 
161430de9b25SBarry Smith   /*MC
16150a57284eSJacob Faibussowitsch    PetscFunctionReturn - Last executable line of each PETSc function used for error
16160a57284eSJacob Faibussowitsch    handling. Replaces `return()`.
161730de9b25SBarry Smith 
161830de9b25SBarry Smith    Synopsis:
16190a57284eSJacob Faibussowitsch    #include <petscerror.h>
16200a57284eSJacob Faibussowitsch    void PetscFunctionReturn(...)
162130de9b25SBarry Smith 
1622667f096bSBarry Smith    Not Collective; No Fortran Support
1623eca87e8dSBarry Smith 
16245687f061SJacob Faibussowitsch    Level: beginner
16255687f061SJacob Faibussowitsch 
16260a57284eSJacob Faibussowitsch    Notes:
16270a57284eSJacob Faibussowitsch    This routine is a macro, so while it does not "return" anything itself, it does return from
16280a57284eSJacob Faibussowitsch    the function in the literal sense.
16290a57284eSJacob Faibussowitsch 
16300a57284eSJacob Faibussowitsch    Usually the return value is the integer literal `0` (for example in any function returning
16310a57284eSJacob Faibussowitsch    `PetscErrorCode`), however it is possible to return any arbitrary type. The arguments of
16320a57284eSJacob Faibussowitsch    this macro are placed before the `return` statement as-is.
16330a57284eSJacob Faibussowitsch 
16340a57284eSJacob Faibussowitsch    Any routine which returns via `PetscFunctionReturn()` must begin with a corresponding
16350a57284eSJacob Faibussowitsch    `PetscFunctionBegin`.
16360a57284eSJacob Faibussowitsch 
16370a57284eSJacob Faibussowitsch    For routines which return `void` use `PetscFunctionReturnVoid()` instead.
16380a57284eSJacob Faibussowitsch 
16390a57284eSJacob Faibussowitsch    Example Usage:
164030de9b25SBarry Smith .vb
16410a57284eSJacob Faibussowitsch    PetscErrorCode foo(int *x)
16420a57284eSJacob Faibussowitsch    {
16430a57284eSJacob Faibussowitsch      PetscFunctionBegin; // don't forget the begin!
16440a57284eSJacob Faibussowitsch      *x = 10;
16453ba16761SJacob Faibussowitsch      PetscFunctionReturn(PETSC_SUCCESS);
164630de9b25SBarry Smith    }
164730de9b25SBarry Smith .ve
164830de9b25SBarry Smith 
16490a57284eSJacob Faibussowitsch    May return any arbitrary type\:
16500a57284eSJacob Faibussowitsch .vb
16510a57284eSJacob Faibussowitsch   struct Foo
16520a57284eSJacob Faibussowitsch   {
16530a57284eSJacob Faibussowitsch     int x;
16540a57284eSJacob Faibussowitsch   };
16550a57284eSJacob Faibussowitsch 
16560a57284eSJacob Faibussowitsch   struct Foo make_foo(int value)
16570a57284eSJacob Faibussowitsch   {
16580a57284eSJacob Faibussowitsch     struct Foo f;
16590a57284eSJacob Faibussowitsch 
16600a57284eSJacob Faibussowitsch     PetscFunctionBegin;
16610a57284eSJacob Faibussowitsch     f.x = value;
16620a57284eSJacob Faibussowitsch     PetscFunctionReturn(f);
16630a57284eSJacob Faibussowitsch   }
16640a57284eSJacob Faibussowitsch .ve
16650a57284eSJacob Faibussowitsch 
16660a57284eSJacob Faibussowitsch .seealso: `PetscFunctionBegin`, `PetscFunctionBeginUser`, `PetscFunctionReturnVoid()`,
16670a57284eSJacob Faibussowitsch           `PetscStackPopNoCheck()`
166830de9b25SBarry Smith M*/
16695687f061SJacob Faibussowitsch   #define PetscFunctionReturn(...) \
16709371c9d4SSatish Balay     do { \
1671362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
16725687f061SJacob Faibussowitsch       return __VA_ARGS__; \
167327104ee2SJacob Faibussowitsch     } while (0)
1674d64ed03dSBarry Smith 
16750a57284eSJacob Faibussowitsch   /*MC
16760a57284eSJacob Faibussowitsch   PetscFunctionReturnVoid - Like `PetscFunctionReturn()` but returns `void`
16770a57284eSJacob Faibussowitsch 
16780a57284eSJacob Faibussowitsch   Synopsis:
16790a57284eSJacob Faibussowitsch   #include <petscerror.h>
16800a57284eSJacob Faibussowitsch   void PetscFunctionReturnVoid()
16810a57284eSJacob Faibussowitsch 
16820a57284eSJacob Faibussowitsch   Not Collective
16830a57284eSJacob Faibussowitsch 
16845687f061SJacob Faibussowitsch   Level: beginner
16855687f061SJacob Faibussowitsch 
16865687f061SJacob Faibussowitsch   Note:
16870a57284eSJacob Faibussowitsch   Behaves identically to `PetscFunctionReturn()` except that it returns `void`. That is, this
16885687f061SJacob Faibussowitsch   macro culminates with `return`.
16890a57284eSJacob Faibussowitsch 
16900a57284eSJacob Faibussowitsch   Example Usage:
16910a57284eSJacob Faibussowitsch .vb
16920a57284eSJacob Faibussowitsch   void foo()
16930a57284eSJacob Faibussowitsch   {
16940a57284eSJacob Faibussowitsch     PetscFunctionBegin; // must start with PetscFunctionBegin!
16950a57284eSJacob Faibussowitsch     bar();
16960a57284eSJacob Faibussowitsch     baz();
16970a57284eSJacob Faibussowitsch     PetscFunctionReturnVoid();
16980a57284eSJacob Faibussowitsch   }
16990a57284eSJacob Faibussowitsch .ve
17000a57284eSJacob Faibussowitsch 
17010a57284eSJacob Faibussowitsch .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, PetscFunctionBeginUser`
17020a57284eSJacob Faibussowitsch M*/
17039371c9d4SSatish Balay   #define PetscFunctionReturnVoid() \
17049371c9d4SSatish Balay     do { \
1705362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
170627104ee2SJacob Faibussowitsch       return; \
170727104ee2SJacob Faibussowitsch     } while (0)
170827104ee2SJacob Faibussowitsch #else /* PETSC_USE_DEBUG */
170927104ee2SJacob Faibussowitsch   #define PetscStackPushNoCheck(funct, petsc_routine, hot)
171037154d10SBarry Smith   #define PetscStackUpdateLine
1711792fecdfSBarry Smith   #define PetscStackPushExternal(funct)
17123ba16761SJacob Faibussowitsch   #define PetscStackPopNoCheck(...)
171327104ee2SJacob Faibussowitsch   #define PetscStackClearTop
17143a40ed3dSBarry Smith   #define PetscFunctionBegin
17150bdf7c52SPeter Brune   #define PetscFunctionBeginUser
1716a2f94806SJed Brown   #define PetscFunctionBeginHot
17170a57284eSJacob Faibussowitsch   #define PetscFunctionReturn(...)  return __VA_ARGS__
17185665465eSBarry Smith   #define PetscFunctionReturnVoid() return
1719812af9f3SBarry Smith   #define PetscStackPop             CHKMEMQ
1720812af9f3SBarry Smith   #define PetscStackPush(f)         CHKMEMQ
172127104ee2SJacob Faibussowitsch #endif /* PETSC_USE_DEBUG */
17223a40ed3dSBarry Smith 
17236d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
17243ba16761SJacob Faibussowitsch   #define PetscStackCallExternalVoid(...)
17253ba16761SJacob Faibussowitsch template <typename F, typename... Args>
17263ba16761SJacob Faibussowitsch void PetscCallExternal(F, Args...);
17276d210af2SJacob Faibussowitsch #else
1728586f9135SBarry Smith   /*MC
1729e77caa6dSBarry Smith     PetscStackCallExternalVoid - Calls an external library routine or user function after pushing the name of the routine on the stack.
1730eb6b5d47SBarry Smith 
1731eb6b5d47SBarry Smith    Input Parameters:
1732eb6b5d47SBarry Smith +   name - string that gives the name of the function being called
1733586f9135SBarry Smith -   routine - actual call to the routine, for example, functionname(a,b)
1734fd3f9acdSBarry Smith 
1735586f9135SBarry Smith    Level: developer
1736eb6b5d47SBarry Smith 
1737586f9135SBarry Smith    Note:
1738792fecdfSBarry Smith    Often one should use `PetscCallExternal()` instead. This routine is intended for external library routines that DO NOT return error codes
1739eb6b5d47SBarry Smith 
1740586f9135SBarry Smith    In debug mode this also checks the memory for corruption at the end of the function call.
1741586f9135SBarry Smith 
1742586f9135SBarry Smith    Certain external packages, such as BLAS/LAPACK may have their own macros for managing the call, error checking, etc.
1743586f9135SBarry Smith 
1744586f9135SBarry Smith    Developer Note:
1745586f9135SBarry 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.
1746586f9135SBarry Smith 
1747792fecdfSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscCallExternal()`, `PetscCallBLAS()`
1748586f9135SBarry Smith @*/
17493ba16761SJacob Faibussowitsch   #define PetscStackCallExternalVoid(name, ...) \
17509371c9d4SSatish Balay     do { \
17519371c9d4SSatish Balay       PetscStackPush(name); \
17523ba16761SJacob Faibussowitsch       __VA_ARGS__; \
17539371c9d4SSatish Balay       PetscStackPop; \
17549371c9d4SSatish Balay     } while (0)
1755eb6b5d47SBarry Smith 
1756586f9135SBarry Smith   /*MC
1757792fecdfSBarry Smith     PetscCallExternal - Calls an external library routine that returns an error code after pushing the name of the routine on the stack.
1758fd3f9acdSBarry Smith 
1759fd3f9acdSBarry Smith    Input Parameters:
1760fd3f9acdSBarry Smith +   func-  name of the routine
1761586f9135SBarry Smith -   args - arguments to the routine
1762586f9135SBarry Smith 
1763586f9135SBarry Smith    Level: developer
1764fd3f9acdSBarry Smith 
176595452b02SPatrick Sanan    Notes:
1766e77caa6dSBarry Smith    This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not.
1767dbf62e16SBarry Smith 
1768586f9135SBarry Smith    In debug mode this also checks the memory for corruption at the end of the function call.
1769fd3f9acdSBarry Smith 
177087497f52SBarry Smith    Assumes the error return code of the function is an integer and that a value of 0 indicates success
177187497f52SBarry Smith 
1772586f9135SBarry Smith    Developer Note:
1773d5b43468SJose E. Roman    This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc.
1774586f9135SBarry Smith 
1775e77caa6dSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()`
1776586f9135SBarry Smith M*/
17779371c9d4SSatish Balay   #define PetscCallExternal(func, ...) \
17789371c9d4SSatish Balay     do { \
1779a74df02fSJacob Faibussowitsch       PetscStackPush(PetscStringize(func)); \
17803ba16761SJacob Faibussowitsch       int ierr_petsc_call_external_ = func(__VA_ARGS__); \
17811d4906efSStefano Zampini       PetscStackPop; \
17823ba16761SJacob 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_); \
1783fd3f9acdSBarry Smith     } while (0)
17846d210af2SJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */
178506d1fe2cSBarry Smith 
178606d1fe2cSBarry Smith #endif
1787