xref: /petsc/include/petscerror.h (revision a1fd7ae3ce1a484b2c7fddb18ad0962aba6eb1a3)
154a8ef01SBarry Smith /*
2f621e05eSBarry Smith     Contains all error handling interfaces for PETSc.
354a8ef01SBarry Smith */
4a4963045SJacob Faibussowitsch #pragma once
56c7e564aSBarry Smith 
6e1bf4ed2SJacob Faibussowitsch #include <petscmacros.h>
7e1bf4ed2SJacob Faibussowitsch #include <petscsystypes.h>
8e1bf4ed2SJacob Faibussowitsch 
9af75f258SJacob Faibussowitsch #if defined(__cplusplus)
10af75f258SJacob Faibussowitsch   #include <exception> // std::exception
11af75f258SJacob Faibussowitsch #endif
12af75f258SJacob Faibussowitsch 
13ac09b921SBarry Smith /* SUBMANSEC = Sys */
14ac09b921SBarry Smith 
15edd03b47SJacob Faibussowitsch #define SETERRQ1(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
16edd03b47SJacob Faibussowitsch #define SETERRQ2(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
17edd03b47SJacob Faibussowitsch #define SETERRQ3(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
18edd03b47SJacob Faibussowitsch #define SETERRQ4(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
19edd03b47SJacob Faibussowitsch #define SETERRQ5(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
20edd03b47SJacob Faibussowitsch #define SETERRQ6(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
21edd03b47SJacob Faibussowitsch #define SETERRQ7(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
22edd03b47SJacob Faibussowitsch #define SETERRQ8(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
23edd03b47SJacob Faibussowitsch #define SETERRQ9(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
2498921bdaSJacob Faibussowitsch 
2530de9b25SBarry Smith /*MC
261957e957SBarry Smith    SETERRQ - Macro to be called when an error has been detected,
2730de9b25SBarry Smith 
2830de9b25SBarry Smith    Synopsis:
29aaa7dc30SBarry Smith    #include <petscsys.h>
3098921bdaSJacob Faibussowitsch    PetscErrorCode SETERRQ(MPI_Comm comm, PetscErrorCode ierr, char *message, ...)
3130de9b25SBarry Smith 
32d083f849SBarry Smith    Collective
3330de9b25SBarry Smith 
3430de9b25SBarry Smith    Input Parameters:
3595bd0b28SBarry Smith +  comm    - An MPI communicator, use `PETSC_COMM_SELF` unless you know all ranks of another communicator will detect the error
363af045c5SBarry Smith .  ierr    - nonzero error code, see the list of standard error codes in include/petscerror.h
3730de9b25SBarry Smith -  message - error message
3830de9b25SBarry Smith 
3930de9b25SBarry Smith   Level: beginner
4030de9b25SBarry Smith 
4130de9b25SBarry Smith    Notes:
4287497f52SBarry Smith    This is rarely needed, one should use `PetscCheck()` and `PetscCall()` and friends to automatically handle error conditions.
4330de9b25SBarry Smith    Once the error handler is called the calling function is then returned from with the given error code.
4430de9b25SBarry Smith 
4587497f52SBarry Smith    Experienced users can set the error handler with `PetscPushErrorHandler()`.
4630de9b25SBarry Smith 
4716a05f60SBarry Smith    Fortran Note:
48989c49a2SBarry Smith    `SETERRQ()` may be called from Fortran subroutines but `SETERRA()` must be called from the
493b1008b8SBarry Smith    Fortran main program.
503b1008b8SBarry Smith 
51db781477SPatrick Sanan .seealso: `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
52af27ebaaSBarry Smith           `PetscError()`, `PetscCall()`, `CHKMEMQ`, `CHKERRA()`, `PetscCallMPI()`, `PetscErrorCode`
5330de9b25SBarry Smith M*/
54e809461dSJacob Faibussowitsch #define SETERRQ(comm, ierr, ...) \
55e809461dSJacob Faibussowitsch   do { \
56e809461dSJacob Faibussowitsch     PetscErrorCode ierr_seterrq_petsc_ = PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \
57e809461dSJacob Faibussowitsch     return ierr_seterrq_petsc_ ? ierr_seterrq_petsc_ : PETSC_ERR_RETURN; \
58e809461dSJacob Faibussowitsch   } while (0)
59986eef2eSBarry Smith 
60648c30bcSBarry Smith #define SETERRQNULL(comm, ierr, ...) \
61648c30bcSBarry Smith   do { \
62648c30bcSBarry Smith     (void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \
63648c30bcSBarry Smith     return NULL; \
64648c30bcSBarry Smith   } while (0)
65648c30bcSBarry Smith 
66ffc4695bSBarry Smith /*
67ffc4695bSBarry Smith     Returned from PETSc functions that are called from MPI, such as related to attributes
68ffc4695bSBarry Smith       Do not confuse PETSC_MPI_ERROR_CODE and PETSC_ERR_MPI, the first is registered with MPI and returned to MPI as
69888a1fb3SBarry Smith       an error code, the latter is a regular PETSc error code passed within PETSc code indicating an error was detected in an MPI call.
70ffc4695bSBarry Smith */
71ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CLASS;
72ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CODE;
73ffc4695bSBarry Smith 
74986eef2eSBarry Smith /*MC
75986eef2eSBarry Smith    SETERRMPI - Macro to be called when an error has been detected within an MPI callback function
76986eef2eSBarry Smith 
77989c49a2SBarry Smith    No Fortran Support
78989c49a2SBarry Smith 
79986eef2eSBarry Smith    Synopsis:
80986eef2eSBarry Smith    #include <petscsys.h>
8198921bdaSJacob Faibussowitsch    PetscErrorCode SETERRMPI(MPI_Comm comm, PetscErrorCode ierr, char *message, ...)
82986eef2eSBarry Smith 
83d083f849SBarry Smith    Collective
84986eef2eSBarry Smith 
85986eef2eSBarry Smith    Input Parameters:
8695bd0b28SBarry Smith +  comm    - An MPI communicator, use `PETSC_COMM_SELF` unless you know all ranks of another communicator will detect the error
87986eef2eSBarry Smith .  ierr    - nonzero error code, see the list of standard error codes in include/petscerror.h
88986eef2eSBarry Smith -  message - error message
89986eef2eSBarry Smith 
90986eef2eSBarry Smith   Level: developer
91986eef2eSBarry Smith 
9295bd0b28SBarry Smith    Note:
9387497f52SBarry Smith    This macro is FOR USE IN MPI CALLBACK FUNCTIONS ONLY, such as those passed to `MPI_Comm_create_keyval()`. It always returns the error code `PETSC_MPI_ERROR_CODE`
9487497f52SBarry Smith   which is registered with `MPI_Add_error_code()` when PETSc is initialized.
95986eef2eSBarry Smith 
96af27ebaaSBarry Smith .seealso: `SETERRQ()`, `PetscCall()`, `PetscCallMPI()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `PetscErrorCode`
97986eef2eSBarry Smith M*/
983ba16761SJacob Faibussowitsch #define SETERRMPI(comm, ierr, ...) return ((void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__), PETSC_MPI_ERROR_CODE)
99ee8199e6SMatthew G. Knepley 
100ee8199e6SMatthew G. Knepley /*MC
101f388eb8bSPatrick Sanan    SETERRA - Fortran-only macro that can be called when an error has been detected from the main program
102f388eb8bSPatrick Sanan 
103f388eb8bSPatrick Sanan    Synopsis:
104f388eb8bSPatrick Sanan    #include <petscsys.h>
105f388eb8bSPatrick Sanan    PetscErrorCode SETERRA(MPI_Comm comm, PetscErrorCode ierr, char *message)
106f388eb8bSPatrick Sanan 
107f388eb8bSPatrick Sanan    Collective
108f388eb8bSPatrick Sanan 
109f388eb8bSPatrick Sanan    Input Parameters:
11095bd0b28SBarry Smith +  comm    - An MPI communicator, so that the error can be collective
111f388eb8bSPatrick Sanan .  ierr    - nonzero error code, see the list of standard error codes in include/petscerror.h
112f388eb8bSPatrick Sanan -  message - error message in the printf format
113f388eb8bSPatrick Sanan 
114f388eb8bSPatrick Sanan    Level: beginner
115f388eb8bSPatrick Sanan 
116f388eb8bSPatrick Sanan    Notes:
11787497f52SBarry Smith    This should only be used with Fortran. With C/C++, use `SETERRQ()`.
118f388eb8bSPatrick Sanan 
11987497f52SBarry Smith    `SETERRQ()` may be called from Fortran subroutines but `SETERRA()` must be called from the
120f388eb8bSPatrick Sanan     Fortran main program.
121f388eb8bSPatrick Sanan 
122af27ebaaSBarry Smith .seealso: `SETERRQ()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()`, `PetscErrorCode`
123f388eb8bSPatrick Sanan M*/
124f388eb8bSPatrick Sanan 
125f388eb8bSPatrick Sanan /*MC
126fa190f98SMatthew G. Knepley    SETERRABORT - Macro that can be called when an error has been detected,
127fa190f98SMatthew G. Knepley 
128fa190f98SMatthew G. Knepley    Synopsis:
129fa190f98SMatthew G. Knepley    #include <petscsys.h>
13098921bdaSJacob Faibussowitsch    PetscErrorCode SETERRABORT(MPI_Comm comm, PetscErrorCode ierr, char *message, ...)
131fa190f98SMatthew G. Knepley 
132d083f849SBarry Smith    Collective
133fa190f98SMatthew G. Knepley 
134fa190f98SMatthew G. Knepley    Input Parameters:
13595bd0b28SBarry Smith +  comm    - An MPI communicator, so that the error can be collective
1363af045c5SBarry Smith .  ierr    - nonzero error code, see the list of standard error codes in include/petscerror.h
137fa190f98SMatthew G. Knepley -  message - error message in the printf format
138fa190f98SMatthew G. Knepley 
139fa190f98SMatthew G. Knepley    Level: beginner
140fa190f98SMatthew G. Knepley 
141fa190f98SMatthew G. Knepley    Notes:
14287497f52SBarry Smith    This function just calls `MPI_Abort()`.
14387497f52SBarry Smith 
14487497f52SBarry Smith    This should only be called in routines that cannot return an error code, such as in C++ constructors.
145fa190f98SMatthew G. Knepley 
146989c49a2SBarry Smith    Fortran Note:
147989c49a2SBarry Smith    Use `SETERRA()` in Fortran main program and `SETERRQ()` in Fortran subroutines
148989c49a2SBarry Smith 
149989c49a2SBarry Smith    Developer Note:
150989c49a2SBarry Smith    In Fortran `SETERRA()` could be called `SETERRABORT()` since they serve the same purpose
151989c49a2SBarry Smith 
152af27ebaaSBarry Smith .seealso: `SETERRQ()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `PetscCall()`, `CHKMEMQ`, `PetscErrorCode`
153fa190f98SMatthew G. Knepley M*/
1549371c9d4SSatish Balay #define SETERRABORT(comm, ierr, ...) \
1559371c9d4SSatish Balay   do { \
1563ba16761SJacob Faibussowitsch     (void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \
15798921bdaSJacob Faibussowitsch     MPI_Abort(comm, ierr); \
15898921bdaSJacob Faibussowitsch   } while (0)
1599a00fa46SSatish Balay 
16030de9b25SBarry Smith /*MC
161139c8099SSatish Balay   PetscCheck - Checks that a particular condition is true; if not true, then returns the provided error code
1622c71b3e2SJacob Faibussowitsch 
1632c71b3e2SJacob Faibussowitsch   Synopsis:
1642c71b3e2SJacob Faibussowitsch   #include <petscerror.h>
1652c71b3e2SJacob Faibussowitsch   void PetscCheck(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
1662c71b3e2SJacob Faibussowitsch 
1677cdbe19fSJose E. Roman   Collective; No Fortran Support
1682c71b3e2SJacob Faibussowitsch 
1692c71b3e2SJacob Faibussowitsch   Input Parameters:
1702c71b3e2SJacob Faibussowitsch + cond    - The boolean condition
1712c71b3e2SJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
1722c71b3e2SJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
1732c71b3e2SJacob Faibussowitsch - message - Error message in printf format
1742c71b3e2SJacob Faibussowitsch 
175667f096bSBarry Smith   Level: beginner
176667f096bSBarry Smith 
1772c71b3e2SJacob Faibussowitsch   Notes:
1782c71b3e2SJacob Faibussowitsch   Enabled in both optimized and debug builds.
1792c71b3e2SJacob Faibussowitsch 
180a2454668SBarry Smith   As a general rule, `PetscCheck()` is used to check "usage error" (for example, passing an incorrect value as a function argument),
181a2454668SBarry Smith   `PetscAssert()` is used to "check for bugs in PETSc" (for example, is a value in a PETSc data structure nonsensical).
182a2454668SBarry Smith   However, for functions that are called in a "hot spot", for example, thousands of times in a loop, `PetscAssert()` should be used instead
183a2454668SBarry Smith   of `PetscCheck()` since the former is compiled out in PETSc's optimization code.
184a2454668SBarry Smith 
18587497f52SBarry Smith   Calls `SETERRQ()` if the assertion fails, so can only be called from functions returning a
18687497f52SBarry Smith   `PetscErrorCode` (or equivalent type after conversion).
1872c71b3e2SJacob Faibussowitsch 
188af27ebaaSBarry Smith .seealso: `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheckAbort()`, `PetscErrorCode`
1899566063dSJacob Faibussowitsch M*/
1909371c9d4SSatish Balay #define PetscCheck(cond, comm, ierr, ...) \
1919371c9d4SSatish Balay   do { \
1929371c9d4SSatish Balay     if (PetscUnlikely(!(cond))) SETERRQ(comm, ierr, __VA_ARGS__); \
1939371c9d4SSatish Balay   } while (0)
1942c71b3e2SJacob Faibussowitsch 
1952c71b3e2SJacob Faibussowitsch /*MC
1964be741a6SBarry Smith   PetscCheckAbort - Check that a particular condition is true, otherwise prints error and aborts
1974be741a6SBarry Smith 
1984be741a6SBarry Smith   Synopsis:
1994be741a6SBarry Smith   #include <petscerror.h>
2004be741a6SBarry Smith   void PetscCheckAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
2014be741a6SBarry Smith 
2027cdbe19fSJose E. Roman   Collective; No Fortran Support
2034be741a6SBarry Smith 
2044be741a6SBarry Smith   Input Parameters:
2054be741a6SBarry Smith + cond    - The boolean condition
2064be741a6SBarry Smith . comm    - The communicator on which the check can be collective on
2074be741a6SBarry Smith . ierr    - A nonzero error code, see include/petscerror.h for the complete list
2084be741a6SBarry Smith - message - Error message in printf format
2094be741a6SBarry Smith 
210667f096bSBarry Smith   Level: developer
211667f096bSBarry Smith 
2124be741a6SBarry Smith   Notes:
2134be741a6SBarry Smith   Enabled in both optimized and debug builds.
2144be741a6SBarry Smith 
21587497f52SBarry Smith   Calls `SETERRABORT()` if the assertion fails, can be called from a function that does not return an
21687497f52SBarry Smith   error code, such as a C++ constructor. usually `PetscCheck()` should be used.
2174be741a6SBarry Smith 
218af27ebaaSBarry Smith .seealso: `PetscAssertAbort()`, `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheck()`, `SETERRABORT()`, `PetscErrorCode`
2194be741a6SBarry Smith M*/
2209371c9d4SSatish Balay #define PetscCheckAbort(cond, comm, ierr, ...) \
2210e6b6b59SJacob Faibussowitsch   do { \
2220e6b6b59SJacob Faibussowitsch     if (PetscUnlikely(!(cond))) SETERRABORT(comm, ierr, __VA_ARGS__); \
223c7481402SJacob Faibussowitsch   } while (0)
2244be741a6SBarry Smith 
2254be741a6SBarry Smith /*MC
2269ace16cdSJacob Faibussowitsch   PetscAssert - Assert that a particular condition is true
2279ace16cdSJacob Faibussowitsch 
2289ace16cdSJacob Faibussowitsch   Synopsis:
2299ace16cdSJacob Faibussowitsch   #include <petscerror.h>
2309ace16cdSJacob Faibussowitsch   void PetscAssert(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
2319ace16cdSJacob Faibussowitsch 
2327cdbe19fSJose E. Roman   Collective; No Fortran Support
2339ace16cdSJacob Faibussowitsch 
2349ace16cdSJacob Faibussowitsch   Input Parameters:
2359ace16cdSJacob Faibussowitsch + cond    - The boolean condition
2369ace16cdSJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
2379ace16cdSJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
238af27ebaaSBarry Smith - message - Error message in `printf()` format
2399ace16cdSJacob Faibussowitsch 
240667f096bSBarry Smith   Level: beginner
241667f096bSBarry Smith 
2429ace16cdSJacob Faibussowitsch   Notes:
243c7481402SJacob Faibussowitsch   Equivalent to `PetscCheck()` if debugging is enabled, and `PetscAssume(cond)` otherwise.
2442c71b3e2SJacob Faibussowitsch 
24587497f52SBarry Smith   See `PetscCheck()` for usage and behaviour.
24687497f52SBarry Smith 
24787497f52SBarry Smith   This is needed instead of simply using `assert()` because this correctly handles the collective nature of errors under MPI
2489ace16cdSJacob Faibussowitsch 
249af27ebaaSBarry Smith .seealso: `PetscCheck()`, `SETERRQ()`, `PetscError()`, `PetscAssertAbort()`, `PetscErrorCode`
2509566063dSJacob Faibussowitsch M*/
251c7481402SJacob Faibussowitsch #if PetscDefined(USE_DEBUG)
252c7481402SJacob Faibussowitsch   #define PetscAssert(cond, comm, ierr, ...) PetscCheck(cond, comm, ierr, __VA_ARGS__)
253c7481402SJacob Faibussowitsch #else
254c7481402SJacob Faibussowitsch   #define PetscAssert(cond, ...) PetscAssume(cond)
255c7481402SJacob Faibussowitsch #endif
2569ace16cdSJacob Faibussowitsch 
2579ace16cdSJacob Faibussowitsch /*MC
2580e6b6b59SJacob Faibussowitsch   PetscAssertAbort - Assert that a particular condition is true, otherwise prints error and aborts
2590e6b6b59SJacob Faibussowitsch 
2600e6b6b59SJacob Faibussowitsch   Synopsis:
2610e6b6b59SJacob Faibussowitsch   #include <petscerror.h>
2620e6b6b59SJacob Faibussowitsch   void PetscAssertAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
2630e6b6b59SJacob Faibussowitsch 
2647cdbe19fSJose E. Roman   Collective; No Fortran Support
2650e6b6b59SJacob Faibussowitsch 
2660e6b6b59SJacob Faibussowitsch   Input Parameters:
2670e6b6b59SJacob Faibussowitsch + cond    - The boolean condition
2680e6b6b59SJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
2690e6b6b59SJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
2700e6b6b59SJacob Faibussowitsch - message - Error message in printf format
2710e6b6b59SJacob Faibussowitsch 
272667f096bSBarry Smith   Level: beginner
273667f096bSBarry Smith 
27495bd0b28SBarry Smith   Note:
2750e6b6b59SJacob Faibussowitsch   Enabled only in debug builds. See `PetscCheckAbort()` for usage.
2760e6b6b59SJacob Faibussowitsch 
2770e6b6b59SJacob Faibussowitsch .seealso: `PetscCheckAbort()`, `PetscAssert()`, `PetscCheck()`, `SETERRABORT()`, `PetscError()`
2780e6b6b59SJacob Faibussowitsch M*/
2793ba16761SJacob Faibussowitsch #if PetscDefined(USE_DEBUG)
2803ba16761SJacob Faibussowitsch   #define PetscAssertAbort(cond, comm, ierr, ...) PetscCheckAbort(cond, comm, ierr, __VA_ARGS__)
2813ba16761SJacob Faibussowitsch #else
2823ba16761SJacob Faibussowitsch   #define PetscAssertAbort(cond, comm, ierr, ...) PetscAssume(cond)
2833ba16761SJacob Faibussowitsch #endif
2840e6b6b59SJacob Faibussowitsch 
2850e6b6b59SJacob Faibussowitsch /*MC
2863ba16761SJacob Faibussowitsch   PetscCall - Calls a PETSc function and then checks the resulting error code, if it is
2873ba16761SJacob Faibussowitsch   non-zero it calls the error handler and returns from the current function with the error
2883ba16761SJacob Faibussowitsch   code.
2899566063dSJacob Faibussowitsch 
2909566063dSJacob Faibussowitsch   Synopsis:
2919566063dSJacob Faibussowitsch   #include <petscerror.h>
29249c86fc7SBarry Smith   void PetscCall(PetscFunction(args))
2939566063dSJacob Faibussowitsch 
2949566063dSJacob Faibussowitsch   Not Collective
2959566063dSJacob Faibussowitsch 
2969566063dSJacob Faibussowitsch   Input Parameter:
29749c86fc7SBarry Smith . PetscFunction - any PETSc function that returns an error code
2989566063dSJacob Faibussowitsch 
299667f096bSBarry Smith   Level: beginner
300667f096bSBarry Smith 
3019566063dSJacob Faibussowitsch   Notes:
3029566063dSJacob Faibussowitsch   Once the error handler is called the calling function is then returned from with the given
30387497f52SBarry Smith   error code. Experienced users can set the error handler with `PetscPushErrorHandler()`.
3049566063dSJacob Faibussowitsch 
30587497f52SBarry Smith   `PetscCall()` cannot be used in functions returning a datatype not convertible to
3068af9f190SBarry Smith   `PetscErrorCode`. For example, `PetscCall()` may not be used in functions returning `void`, use
3073ba16761SJacob Faibussowitsch   `PetscCallAbort()` or `PetscCallVoid()` in this case.
3089566063dSJacob Faibussowitsch 
3099566063dSJacob Faibussowitsch   Example Usage:
3109566063dSJacob Faibussowitsch .vb
3119566063dSJacob Faibussowitsch   PetscCall(PetscInitiailize(...)); // OK to call even when PETSc is not yet initialized!
3129566063dSJacob Faibussowitsch 
3139566063dSJacob Faibussowitsch   struct my_struct
3149566063dSJacob Faibussowitsch   {
3159566063dSJacob Faibussowitsch     void *data;
3169566063dSJacob Faibussowitsch   } my_complex_type;
3179566063dSJacob Faibussowitsch 
3189566063dSJacob Faibussowitsch   struct my_struct bar(void)
3199566063dSJacob Faibussowitsch   {
3206aad120cSJose E. Roman     PetscCall(foo(15)); // ERROR PetscErrorCode not convertible to struct my_struct!
3219566063dSJacob Faibussowitsch   }
3229566063dSJacob Faibussowitsch 
3236aad120cSJose E. Roman   PetscCall(bar()) // ERROR input not convertible to PetscErrorCode
3249566063dSJacob Faibussowitsch .ve
3259566063dSJacob Faibussowitsch 
32687497f52SBarry Smith   It is also possible to call this directly on a `PetscErrorCode` variable
32749c86fc7SBarry Smith .vb
32849c86fc7SBarry Smith   PetscCall(ierr);  // check if ierr is nonzero
32949c86fc7SBarry Smith .ve
33049c86fc7SBarry Smith 
331792fecdfSBarry Smith   Should not be used to call callback functions provided by users, `PetscCallBack()` should be used in that situation.
332ef1023bdSBarry Smith 
3336a8be23eSBarry Smith   `PetscUseTypeMethod()` or `PetscTryTypeMethod()` should be used when calling functions pointers contained in a PETSc object's `ops` array
3346a8be23eSBarry Smith 
33549c86fc7SBarry Smith   Fortran Notes:
33698d50953SPierre Jolivet     The Fortran function in which this is used must declare a `PetscErrorCode` variable necessarily named `ierr`, and `ierr` must be
33787497f52SBarry Smith     the final argument to the PETSc function being called.
33849c86fc7SBarry Smith 
33998d50953SPierre Jolivet     In the main program and in Fortran subroutines that do not have `ierr` as the final return parameter, one
34087497f52SBarry Smith     should use `PetscCallA()`
34149c86fc7SBarry Smith 
34249c86fc7SBarry Smith   Example Fortran Usage:
34349c86fc7SBarry Smith .vb
34449c86fc7SBarry Smith   PetscErrorCode ierr
34549c86fc7SBarry Smith   Vec v
34649c86fc7SBarry Smith 
34749c86fc7SBarry Smith   ...
34849c86fc7SBarry Smith   PetscCall(VecShift(v, 1.0, ierr))
34949c86fc7SBarry Smith   PetscCallA(VecShift(v, 1.0, ierr))
35049c86fc7SBarry Smith .ve
35149c86fc7SBarry Smith 
352667f096bSBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`,
353667f096bSBarry Smith           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`,
354648c30bcSBarry Smith           `CHKERRMPI()`, `PetscCallBack()`, `PetscCallAbort()`, `PetscCallVoid()`, `PetscCallNull()`
355648c30bcSBarry Smith M*/
356648c30bcSBarry Smith 
357648c30bcSBarry Smith /*MC
358648c30bcSBarry Smith   PetscCallNull - Calls a PETSc function and then checks the resulting error code, if it is
359648c30bcSBarry Smith   non-zero it calls the error handler and returns a `NULL`
360648c30bcSBarry Smith 
361648c30bcSBarry Smith   Synopsis:
362648c30bcSBarry Smith   #include <petscerror.h>
363648c30bcSBarry Smith   void PetscCallNull(PetscFunction(args))
364648c30bcSBarry Smith 
365648c30bcSBarry Smith   Not Collective; No Fortran Support
366648c30bcSBarry Smith 
367648c30bcSBarry Smith   Input Parameter:
368648c30bcSBarry Smith . PetscFunction - any PETSc function that returns something that can be returned as a `NULL`
369648c30bcSBarry Smith 
370648c30bcSBarry Smith   Level: developer
371648c30bcSBarry Smith 
372648c30bcSBarry Smith .seealso: `PetscCall()`, `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`,
373648c30bcSBarry Smith           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`,
374648c30bcSBarry Smith           `CHKERRMPI()`, `PetscCallBack()`, `PetscCallAbort()`, `PetscCallVoid()`, `PetscCall()`
3759566063dSJacob Faibussowitsch M*/
376ef1023bdSBarry Smith 
377ef1023bdSBarry Smith /*MC
37898d50953SPierre Jolivet    PetscCallA - Fortran-only macro that should be used in the main program and subroutines that do not have `ierr` as the final return parameter, to call PETSc functions instead of using
37998d50953SPierre Jolivet    `PetscCall()` which should be used in other Fortran subroutines
380989c49a2SBarry Smith 
381989c49a2SBarry Smith    Synopsis:
382989c49a2SBarry Smith    #include <petscsys.h>
383989c49a2SBarry Smith    PetscErrorCode PetscCallA(PetscFunction(arguments, ierr))
384989c49a2SBarry Smith 
385989c49a2SBarry Smith    Collective
386989c49a2SBarry Smith 
387989c49a2SBarry Smith    Input Parameter:
388989c49a2SBarry Smith .  PetscFunction(arguments,ierr) - the call to the function
389989c49a2SBarry Smith 
390989c49a2SBarry Smith   Level: beginner
391989c49a2SBarry Smith 
392989c49a2SBarry Smith    Notes:
393989c49a2SBarry Smith    This should only be used with Fortran. With C/C++, use `PetscCall()` always.
394989c49a2SBarry Smith 
39598d50953SPierre Jolivet    The Fortran function in which this is used must declare a `PetscErrorCode` variable necessarily named `ierr`
396989c49a2SBarry Smith    Use `SETERRA()` to set an error in a Fortran main program and `SETERRQ()` in Fortran subroutines
397989c49a2SBarry Smith 
398989c49a2SBarry Smith .seealso: `SETERRQ()`, `SETERRA()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()`
399989c49a2SBarry Smith M*/
400989c49a2SBarry Smith 
401989c49a2SBarry Smith /*MC
402792fecdfSBarry 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
403ef1023bdSBarry Smith   handler and returns from the current function with the error code.
404ef1023bdSBarry Smith 
405ef1023bdSBarry Smith   Synopsis:
406ef1023bdSBarry Smith   #include <petscerror.h>
407792fecdfSBarry Smith   void PetscCallBack(const char *functionname, PetscFunction(args))
408ef1023bdSBarry Smith 
4097cdbe19fSJose E. Roman   Not Collective; No Fortran Support
410ef1023bdSBarry Smith 
411ef1023bdSBarry Smith   Input Parameters:
412ef1023bdSBarry Smith + functionname - the name of the function being called, this can be a string with spaces that describes the meaning of the callback
41387497f52SBarry Smith - PetscFunction - user provided callback function that returns an error code
414ef1023bdSBarry Smith 
415ef1023bdSBarry Smith   Example Usage:
416ef1023bdSBarry Smith .vb
417792fecdfSBarry Smith   PetscCallBack("XXX callback to do something", a->callback(...));
418ef1023bdSBarry Smith .ve
419ef1023bdSBarry Smith 
420ef1023bdSBarry Smith   Level: developer
421ef1023bdSBarry Smith 
422667f096bSBarry Smith   Notes:
4237e6f8dd6SBarry Smith   `PetscUseTypeMethod()` and ` PetscTryTypeMethod()` are the preferred API for this functionality. But when the callback functions are associated with a
4247e6f8dd6SBarry Smith   `DMSNES` or `DMTS` this API must be used.
4257e6f8dd6SBarry Smith 
426667f096bSBarry Smith   Once the error handler is called the calling function is then returned from with the given
427667f096bSBarry Smith   error code. Experienced users can set the error handler with `PetscPushErrorHandler()`.
428667f096bSBarry Smith 
429667f096bSBarry Smith   `PetscCallBack()` should only be called in PETSc when a call is being made to a user provided call-back routine.
430667f096bSBarry Smith 
4317e6f8dd6SBarry Smith   Developer Note:
4327e6f8dd6SBarry Smith   It would be good to provide a new API for when the callbacks are associated with `DMSNES` or `DMTS` so this routine could be used less
4337e6f8dd6SBarry Smith 
43487497f52SBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscCall()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`
4357e6f8dd6SBarry Smith           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, `CHKERRMPI()`, `PetscCall()`,  `PetscUseTypeMethod()`, `PetscTryTypeMethod()`
436ef1023bdSBarry Smith M*/
437ef1023bdSBarry Smith 
4383ba16761SJacob Faibussowitsch /*MC
4398af9f190SBarry Smith   PetscCallVoid - Like `PetscCall()` but for use in functions that return `void`
4403ba16761SJacob Faibussowitsch 
4413ba16761SJacob Faibussowitsch   Synopsis:
4423ba16761SJacob Faibussowitsch   #include <petscerror.h>
4438af9f190SBarry Smith   void PetscCallVoid(PetscFunction(args))
4443ba16761SJacob Faibussowitsch 
4457cdbe19fSJose E. Roman   Not Collective; No Fortran Support
4463ba16761SJacob Faibussowitsch 
4473ba16761SJacob Faibussowitsch   Input Parameter:
4483ba16761SJacob Faibussowitsch . PetscFunction - any PETSc function that returns an error code
4493ba16761SJacob Faibussowitsch 
4503ba16761SJacob Faibussowitsch   Example Usage:
4513ba16761SJacob Faibussowitsch .vb
4523ba16761SJacob Faibussowitsch   void foo()
4533ba16761SJacob Faibussowitsch   {
4543ba16761SJacob Faibussowitsch     KSP ksp;
4553ba16761SJacob Faibussowitsch 
4563ba16761SJacob Faibussowitsch     PetscFunctionBeginUser;
4573ba16761SJacob Faibussowitsch     // OK, properly handles PETSc error codes
4583ba16761SJacob Faibussowitsch     PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp));
4598af9f190SBarry Smith     PetscFunctionReturnVoid();
4603ba16761SJacob Faibussowitsch   }
4613ba16761SJacob Faibussowitsch 
4623ba16761SJacob Faibussowitsch   PetscErrorCode bar()
4633ba16761SJacob Faibussowitsch   {
4643ba16761SJacob Faibussowitsch     KSP ksp;
4653ba16761SJacob Faibussowitsch 
4663ba16761SJacob Faibussowitsch     PetscFunctionBeginUser;
4673ba16761SJacob Faibussowitsch     // ERROR, Non-void function 'bar' should return a value
4683ba16761SJacob Faibussowitsch     PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp));
4693ba16761SJacob Faibussowitsch     // OK, returning PetscErrorCode
4703ba16761SJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
4713ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
4723ba16761SJacob Faibussowitsch   }
473667f096bSBarry Smith .ve
4743ba16761SJacob Faibussowitsch 
4753ba16761SJacob Faibussowitsch   Level: beginner
4763ba16761SJacob Faibussowitsch 
477667f096bSBarry Smith   Notes:
478667f096bSBarry Smith   Has identical usage to `PetscCall()`, except that it returns `void` on error instead of a
479667f096bSBarry Smith   `PetscErrorCode`. See `PetscCall()` for more detailed discussion.
480667f096bSBarry Smith 
481667f096bSBarry Smith   Note that users should prefer `PetscCallAbort()` to this routine. While this routine does
482667f096bSBarry Smith   "handle" errors by returning from the enclosing function, it effectively gobbles the
483667f096bSBarry Smith   error. Since the enclosing function itself returns `void`, its callers have no way of knowing
484667f096bSBarry Smith   that the routine returned early due to an error. `PetscCallAbort()` at least ensures that the
485667f096bSBarry Smith   program crashes gracefully.
486667f096bSBarry Smith 
487648c30bcSBarry Smith .seealso: `PetscCall()`, `PetscErrorCode`, `PetscCallAbort()`, `PetscCallNull()`
4883ba16761SJacob Faibussowitsch M*/
4899566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
4909566063dSJacob Faibussowitsch void PetscCall(PetscErrorCode);
491792fecdfSBarry Smith void PetscCallBack(const char *, PetscErrorCode);
4929566063dSJacob Faibussowitsch void PetscCallVoid(PetscErrorCode);
493648c30bcSBarry Smith void PetscCallNull(PetscErrorCode);
4949566063dSJacob Faibussowitsch #else
4959371c9d4SSatish Balay   #define PetscCall(...) \
4969371c9d4SSatish Balay     do { \
4973ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_q_; \
49837154d10SBarry Smith       PetscStackUpdateLine; \
4993ba16761SJacob Faibussowitsch       ierr_petsc_call_q_ = __VA_ARGS__; \
5003ba16761SJacob 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, " "); \
5019566063dSJacob Faibussowitsch     } while (0)
502648c30bcSBarry Smith   #define PetscCallNull(...) \
503648c30bcSBarry Smith     do { \
504648c30bcSBarry Smith       PetscErrorCode ierr_petsc_call_q_; \
505648c30bcSBarry Smith       PetscStackUpdateLine; \
506648c30bcSBarry Smith       ierr_petsc_call_q_ = __VA_ARGS__; \
507648c30bcSBarry Smith       if (PetscUnlikely(ierr_petsc_call_q_ != PETSC_SUCCESS)) { \
508648c30bcSBarry Smith         (void)PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); \
509648c30bcSBarry Smith         PetscFunctionReturn(NULL); \
510648c30bcSBarry Smith       } \
511648c30bcSBarry Smith     } while (0)
5129371c9d4SSatish Balay   #define PetscCallBack(function, ...) \
5139371c9d4SSatish Balay     do { \
5143ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_q_; \
515e33ced7fSLisandro Dalcin       PetscStackUpdateLine; \
516e33ced7fSLisandro Dalcin       PetscStackPushExternal(function); \
5173ba16761SJacob Faibussowitsch       ierr_petsc_call_q_ = __VA_ARGS__; \
518ef1023bdSBarry Smith       PetscStackPop; \
5193ba16761SJacob 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, " "); \
520ef1023bdSBarry Smith     } while (0)
5219371c9d4SSatish Balay   #define PetscCallVoid(...) \
5229371c9d4SSatish Balay     do { \
5233ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_void_; \
524e33ced7fSLisandro Dalcin       PetscStackUpdateLine; \
5253ba16761SJacob Faibussowitsch       ierr_petsc_call_void_ = __VA_ARGS__; \
5263ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_void_ != PETSC_SUCCESS)) { \
5273ba16761SJacob Faibussowitsch         ierr_petsc_call_void_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_void_, PETSC_ERROR_REPEAT, " "); \
5283ba16761SJacob Faibussowitsch         (void)ierr_petsc_call_void_; \
5299371c9d4SSatish Balay         return; \
5309371c9d4SSatish Balay       } \
5319566063dSJacob Faibussowitsch     } while (0)
5329566063dSJacob Faibussowitsch #endif
5339566063dSJacob Faibussowitsch 
5349566063dSJacob Faibussowitsch /*MC
5359566063dSJacob Faibussowitsch   CHKERRQ - Checks error code returned from PETSc function
53630de9b25SBarry Smith 
53730de9b25SBarry Smith   Synopsis:
538aaa7dc30SBarry Smith   #include <petscsys.h>
5399566063dSJacob Faibussowitsch   void CHKERRQ(PetscErrorCode ierr)
5409566063dSJacob Faibussowitsch 
5419566063dSJacob Faibussowitsch   Not Collective
5429566063dSJacob Faibussowitsch 
5432fe279fdSBarry Smith   Input Parameter:
5449566063dSJacob Faibussowitsch . ierr - nonzero error code
5459566063dSJacob Faibussowitsch 
5469566063dSJacob Faibussowitsch   Level: deprecated
5479566063dSJacob Faibussowitsch 
548667f096bSBarry Smith   Note:
549667f096bSBarry Smith   Deprecated in favor of `PetscCall()`. This routine behaves identically to it.
550667f096bSBarry Smith 
551db781477SPatrick Sanan .seealso: `PetscCall()`
5529566063dSJacob Faibussowitsch M*/
5539566063dSJacob Faibussowitsch #define CHKERRQ(...) PetscCall(__VA_ARGS__)
5549566063dSJacob Faibussowitsch #define CHKERRV(...) PetscCallVoid(__VA_ARGS__)
5559566063dSJacob Faibussowitsch 
556*a1fd7ae3SBarry Smith PETSC_EXTERN void PetscMPIErrorString(PetscMPIInt, size_t, char *);
557db9cea48SBarry Smith 
5589566063dSJacob Faibussowitsch /*MC
5599566063dSJacob Faibussowitsch   PetscCallMPI - Checks error code returned from MPI calls, if non-zero it calls the error
560648c30bcSBarry Smith   handler and then returns a `PetscErrorCode`
5619566063dSJacob Faibussowitsch 
5629566063dSJacob Faibussowitsch   Synopsis:
5639566063dSJacob Faibussowitsch   #include <petscerror.h>
56449c86fc7SBarry Smith   void PetscCallMPI(MPI_Function(args))
56530de9b25SBarry Smith 
566eca87e8dSBarry Smith   Not Collective
56730de9b25SBarry Smith 
5682fe279fdSBarry Smith   Input Parameter:
56949c86fc7SBarry Smith . MPI_Function - an MPI function that returns an MPI error code
57030de9b25SBarry Smith 
571667f096bSBarry Smith   Level: beginner
572667f096bSBarry Smith 
5739566063dSJacob Faibussowitsch   Notes:
57487497f52SBarry Smith   Always returns the error code `PETSC_ERR_MPI`; the MPI error code and string are embedded in
5759566063dSJacob Faibussowitsch   the string error message. Do not use this to call any other routines (for example PETSc
5763ba16761SJacob Faibussowitsch   routines), it should only be used for direct MPI calls. The user may configure PETSc with the
5773ba16761SJacob Faibussowitsch   `--with-strict-petscerrorcode` option to check this at compile-time, otherwise they must
5789566063dSJacob Faibussowitsch   check this themselves.
5799566063dSJacob Faibussowitsch 
580aaa8cc7dSPierre Jolivet   This routine can only be used in functions returning `PetscErrorCode` themselves. If the
5813ba16761SJacob Faibussowitsch   calling function returns a different type, use `PetscCallMPIAbort()` instead.
5823ba16761SJacob Faibussowitsch 
5839566063dSJacob Faibussowitsch   Example Usage:
5849566063dSJacob Faibussowitsch .vb
5859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(...)); // OK, calling MPI function
5869566063dSJacob Faibussowitsch 
5879566063dSJacob Faibussowitsch   PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead!
5889566063dSJacob Faibussowitsch .ve
5899566063dSJacob Faibussowitsch 
59049c86fc7SBarry Smith   Fortran Notes:
59187497f52SBarry Smith     The Fortran function from which this is used must declare a variable `PetscErrorCode` ierr and ierr must be
59249c86fc7SBarry Smith     the final argument to the MPI function being called.
59349c86fc7SBarry Smith 
59449c86fc7SBarry Smith     In the main program and in Fortran subroutines that do not have ierr as the final return parameter one
59587497f52SBarry Smith     should use `PetscCallMPIA()`
59649c86fc7SBarry Smith 
59749c86fc7SBarry Smith   Fortran Usage:
59849c86fc7SBarry Smith .vb
59949c86fc7SBarry Smith   PetscErrorCode ierr or integer ierr
60049c86fc7SBarry Smith   ...
60149c86fc7SBarry Smith   PetscCallMPI(MPI_Comm_size(...,ierr))
60249c86fc7SBarry Smith   PetscCallMPIA(MPI_Comm_size(...,ierr)) ! Will abort after calling error handler
60349c86fc7SBarry Smith 
60449c86fc7SBarry Smith   PetscCallMPI(MPI_Comm_size(...,eflag)) ! ERROR, final argument must be ierr
60549c86fc7SBarry Smith .ve
60649c86fc7SBarry Smith 
607db781477SPatrick Sanan .seealso: `SETERRMPI()`, `PetscCall()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`,
6083ba16761SJacob Faibussowitsch           `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
609648c30bcSBarry Smith           `PetscError()`, `CHKMEMQ`, `PetscCallMPINull()`
610648c30bcSBarry Smith M*/
611648c30bcSBarry Smith 
612648c30bcSBarry Smith /*MC
613648c30bcSBarry Smith   PetscCallMPINull - Checks error code returned from MPI calls, if non-zero it calls the error
614648c30bcSBarry Smith   handler and then returns a `NULL`
615648c30bcSBarry Smith 
616648c30bcSBarry Smith   Synopsis:
617648c30bcSBarry Smith   #include <petscerror.h>
618648c30bcSBarry Smith   void PetscCallMPINull(MPI_Function(args))
619648c30bcSBarry Smith 
620648c30bcSBarry Smith   Not Collective; No Fortran Support
621648c30bcSBarry Smith 
622648c30bcSBarry Smith   Input Parameter:
623648c30bcSBarry Smith . MPI_Function - an MPI function that returns an MPI error code
624648c30bcSBarry Smith 
625648c30bcSBarry Smith   Level: beginner
626648c30bcSBarry Smith 
627648c30bcSBarry Smith   Notes:
628648c30bcSBarry Smith   Always passes the error code `PETSC_ERR_MPI` to the error handler `PetscError()`; the MPI error code and string are embedded in
629648c30bcSBarry Smith   the string error message. Do not use this to call any other routines (for example PETSc
630648c30bcSBarry Smith   routines), it should only be used for direct MPI calls.
631648c30bcSBarry Smith 
632648c30bcSBarry Smith   This routine can only be used in functions returning anything that can be returned as a `NULL` themselves. If the
633648c30bcSBarry Smith   calling function returns a different type, use `PetscCallMPIAbort()` instead.
634648c30bcSBarry Smith 
635648c30bcSBarry Smith   Example Usage:
636648c30bcSBarry Smith .vb
637648c30bcSBarry Smith   PetscCallMPINull(MPI_Comm_size(...)); // OK, calling MPI function
638648c30bcSBarry Smith 
639648c30bcSBarry Smith   PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead!
640648c30bcSBarry Smith .ve
641648c30bcSBarry Smith 
642648c30bcSBarry Smith .seealso: `SETERRMPI()`, `PetscCall()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`,
643648c30bcSBarry Smith           `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
644648c30bcSBarry Smith           `PetscError()`, `CHKMEMQ`, `PetscCallMPI()`
6453ba16761SJacob Faibussowitsch M*/
6463ba16761SJacob Faibussowitsch 
6473ba16761SJacob Faibussowitsch /*MC
6483ba16761SJacob Faibussowitsch   PetscCallMPIAbort - Like `PetscCallMPI()` but calls `MPI_Abort()` on error
6493ba16761SJacob Faibussowitsch 
6503ba16761SJacob Faibussowitsch   Synopsis:
6513ba16761SJacob Faibussowitsch   #include <petscerror.h>
6523ba16761SJacob Faibussowitsch   void PetscCallMPIAbort(MPI_Comm comm, MPI_Function(args))
6533ba16761SJacob Faibussowitsch 
6543ba16761SJacob Faibussowitsch   Not Collective
6553ba16761SJacob Faibussowitsch 
6563ba16761SJacob Faibussowitsch   Input Parameters:
6573ba16761SJacob Faibussowitsch + comm         - the MPI communicator to abort on
6583ba16761SJacob Faibussowitsch - MPI_Function - an MPI function that returns an MPI error code
6593ba16761SJacob Faibussowitsch 
660667f096bSBarry Smith   Level: beginner
661667f096bSBarry Smith 
6623ba16761SJacob Faibussowitsch   Notes:
6633ba16761SJacob Faibussowitsch   Usage is identical to `PetscCallMPI()`. See `PetscCallMPI()` for detailed discussion.
6643ba16761SJacob Faibussowitsch 
6653ba16761SJacob Faibussowitsch   This routine may be used in functions returning `void` or other non-`PetscErrorCode` types.
6663ba16761SJacob Faibussowitsch 
667989c49a2SBarry Smith   Fortran Note:
668989c49a2SBarry Smith   In Fortran this is called `PetscCallMPIA()` and is intended to be used in the main program while `PetscCallMPI()` is
669989c49a2SBarry Smith   used in Fortran subroutines.
670989c49a2SBarry Smith 
671989c49a2SBarry Smith   Developer Note:
672989c49a2SBarry Smith   This should have the same name in Fortran.
673989c49a2SBarry Smith 
6743ba16761SJacob Faibussowitsch .seealso: `PetscCallMPI()`, `PetscCallAbort()`, `SETERRABORT()`
67530de9b25SBarry Smith M*/
6763fcd9f07SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
67763a3b9bcSJacob Faibussowitsch void PetscCallMPI(PetscMPIInt);
6783ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm, PetscMPIInt);
679648c30bcSBarry Smith void PetscCallMPINull(PetscMPIInt);
680064a246eSJacob Faibussowitsch #else
6813ba16761SJacob Faibussowitsch   #define PetscCallMPI_Private(__PETSC_STACK_POP_FUNC__, __SETERR_FUNC__, __COMM__, ...) \
6829371c9d4SSatish Balay     do { \
6833ba16761SJacob Faibussowitsch       PetscMPIInt ierr_petsc_call_mpi_; \
684ef1023bdSBarry Smith       PetscStackUpdateLine; \
685792fecdfSBarry Smith       PetscStackPushExternal("MPI function"); \
686d71ae5a4SJacob Faibussowitsch       { \
6873ba16761SJacob Faibussowitsch         ierr_petsc_call_mpi_ = __VA_ARGS__; \
688d71ae5a4SJacob Faibussowitsch       } \
6893ba16761SJacob Faibussowitsch       __PETSC_STACK_POP_FUNC__; \
6903ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_mpi_ != MPI_SUCCESS)) { \
6913ba16761SJacob Faibussowitsch         char petsc_mpi_7_errorstring[2 * MPI_MAX_ERROR_STRING]; \
692*a1fd7ae3SBarry Smith         PetscMPIErrorString(ierr_petsc_call_mpi_, 2 * MPI_MAX_ERROR_STRING, (char *)petsc_mpi_7_errorstring); \
693*a1fd7ae3SBarry Smith         __SETERR_FUNC__(__COMM__, PETSC_ERR_MPI, "MPI error %d %s", ierr_petsc_call_mpi_, petsc_mpi_7_errorstring); \
6943fcd9f07SJacob Faibussowitsch       } \
6953fcd9f07SJacob Faibussowitsch     } while (0)
6963ba16761SJacob Faibussowitsch 
6973ba16761SJacob Faibussowitsch   #define PetscCallMPI(...)            PetscCallMPI_Private(PetscStackPop, SETERRQ, PETSC_COMM_SELF, __VA_ARGS__)
6983ba16761SJacob Faibussowitsch   #define PetscCallMPIAbort(comm, ...) PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRABORT, comm, __VA_ARGS__)
699648c30bcSBarry Smith   #define PetscCallMPINull(...)        PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRQNULL, PETSC_COMM_SELF, __VA_ARGS__)
7009566063dSJacob Faibussowitsch #endif
701064a246eSJacob Faibussowitsch 
7027037b0edSPatrick Sanan /*MC
7039566063dSJacob Faibussowitsch   CHKERRMPI - Checks error code returned from MPI calls, if non-zero it calls the error
7049566063dSJacob Faibussowitsch   handler and then returns
7059566063dSJacob Faibussowitsch 
7069566063dSJacob Faibussowitsch   Synopsis:
7079566063dSJacob Faibussowitsch   #include <petscerror.h>
7089566063dSJacob Faibussowitsch   void CHKERRMPI(PetscErrorCode ierr)
7099566063dSJacob Faibussowitsch 
7109566063dSJacob Faibussowitsch   Not Collective
7119566063dSJacob Faibussowitsch 
7129566063dSJacob Faibussowitsch   Input Parameter:
7139566063dSJacob Faibussowitsch . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
7149566063dSJacob Faibussowitsch 
7159566063dSJacob Faibussowitsch   Level: deprecated
7169566063dSJacob Faibussowitsch 
717667f096bSBarry Smith   Note:
718667f096bSBarry Smith   Deprecated in favor of `PetscCallMPI()`. This routine behaves identically to it.
719667f096bSBarry Smith 
720db781477SPatrick Sanan .seealso: `PetscCallMPI()`
7219566063dSJacob Faibussowitsch M*/
7229566063dSJacob Faibussowitsch #define CHKERRMPI(...) PetscCallMPI(__VA_ARGS__)
7239566063dSJacob Faibussowitsch 
7249566063dSJacob Faibussowitsch /*MC
725989c49a2SBarry Smith   PetscCallAbort - Checks error code returned from PETSc function, if non-zero it aborts immediately by calling `MPI_Abort()`
7269566063dSJacob Faibussowitsch 
7279566063dSJacob Faibussowitsch   Synopsis:
7289566063dSJacob Faibussowitsch   #include <petscerror.h>
7299566063dSJacob Faibussowitsch   void PetscCallAbort(MPI_Comm comm, PetscErrorCode ierr)
7309566063dSJacob Faibussowitsch 
731c3339decSBarry Smith   Collective
7329566063dSJacob Faibussowitsch 
7339566063dSJacob Faibussowitsch   Input Parameters:
7349566063dSJacob Faibussowitsch + comm - the MPI communicator on which to abort
7359566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
7369566063dSJacob Faibussowitsch 
737667f096bSBarry Smith   Level: intermediate
738667f096bSBarry Smith 
7399566063dSJacob Faibussowitsch   Notes:
74087497f52SBarry Smith   This macro has identical type and usage semantics to `PetscCall()` with the important caveat
7419566063dSJacob Faibussowitsch   that this macro does not return. Instead, if ierr is nonzero it calls the PETSc error handler
74287497f52SBarry Smith   and then immediately calls `MPI_Abort()`. It can therefore be used anywhere.
7439566063dSJacob Faibussowitsch 
74487497f52SBarry Smith   As per `MPI_Abort()` semantics the communicator passed must be valid, although there is currently
74587497f52SBarry Smith   no attempt made at handling any potential errors from `MPI_Abort()`. Note that while
74687497f52SBarry Smith   `MPI_Abort()` is required to terminate only those processes which reside on comm, it is often
74787497f52SBarry Smith   the case that `MPI_Abort()` terminates *all* processes.
7489566063dSJacob Faibussowitsch 
7499566063dSJacob Faibussowitsch   Example Usage:
7509566063dSJacob Faibussowitsch .vb
7519566063dSJacob Faibussowitsch   PetscErrorCode boom(void) { return PETSC_ERR_MEM; }
7529566063dSJacob Faibussowitsch 
7539566063dSJacob Faibussowitsch   void foo(void)
7549566063dSJacob Faibussowitsch   {
7559566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type
7569566063dSJacob Faibussowitsch   }
7579566063dSJacob Faibussowitsch 
7589566063dSJacob Faibussowitsch   double bar(void)
7599566063dSJacob Faibussowitsch   {
7609566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type
7619566063dSJacob Faibussowitsch   }
7629566063dSJacob Faibussowitsch 
7639566063dSJacob Faibussowitsch   PetscCallAbort(MPI_COMM_NULL,boom()); // ERROR, communicator should be valid
7649566063dSJacob Faibussowitsch 
7659566063dSJacob Faibussowitsch   struct baz
7669566063dSJacob Faibussowitsch   {
7679566063dSJacob Faibussowitsch     baz()
7689566063dSJacob Faibussowitsch     {
7699566063dSJacob Faibussowitsch       PetscCallAbort(PETSC_COMM_SELF,boom()); // OK
7709566063dSJacob Faibussowitsch     }
7719566063dSJacob Faibussowitsch 
7729566063dSJacob Faibussowitsch     ~baz()
7739566063dSJacob Faibussowitsch     {
7749566063dSJacob Faibussowitsch       PetscCallAbort(PETSC_COMM_SELF,boom()); // OK (in fact the only way to handle PETSc errors)
7759566063dSJacob Faibussowitsch     }
7769566063dSJacob Faibussowitsch   };
7779566063dSJacob Faibussowitsch .ve
7789566063dSJacob Faibussowitsch 
779989c49a2SBarry Smith   Fortran Note:
780989c49a2SBarry Smith   Use `PetscCallA()`.
781989c49a2SBarry Smith 
782989c49a2SBarry Smith   Developer Note:
783989c49a2SBarry Smith   This should have the same name in Fortran as in C.
784989c49a2SBarry Smith 
785db781477SPatrick Sanan .seealso: `SETERRABORT()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`,
7865687f061SJacob Faibussowitsch           `SETERRQ()`, `CHKMEMQ`, `PetscCallMPI()`, `PetscCallCXXAbort()`
7879566063dSJacob Faibussowitsch M*/
7889566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
7899566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm, PetscErrorCode);
7909566063dSJacob Faibussowitsch void PetscCallContinue(PetscErrorCode);
7919566063dSJacob Faibussowitsch #else
7929371c9d4SSatish Balay   #define PetscCallAbort(comm, ...) \
7939371c9d4SSatish Balay     do { \
7947da6b3ccSLisandro Dalcin       PetscErrorCode ierr_petsc_call_abort_; \
7953ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
7967da6b3ccSLisandro Dalcin       ierr_petsc_call_abort_ = __VA_ARGS__; \
7973ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_abort_ != PETSC_SUCCESS)) { \
7983ba16761SJacob Faibussowitsch         ierr_petsc_call_abort_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_abort_, PETSC_ERROR_REPEAT, " "); \
7993ba16761SJacob Faibussowitsch         (void)MPI_Abort(comm, (PetscMPIInt)ierr_petsc_call_abort_); \
8009566063dSJacob Faibussowitsch       } \
8019566063dSJacob Faibussowitsch     } while (0)
8029371c9d4SSatish Balay   #define PetscCallContinue(...) \
8039371c9d4SSatish Balay     do { \
8047da6b3ccSLisandro Dalcin       PetscErrorCode ierr_petsc_call_continue_; \
8053ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
8067da6b3ccSLisandro Dalcin       ierr_petsc_call_continue_ = __VA_ARGS__; \
8073ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_continue_ != PETSC_SUCCESS)) { \
8083ba16761SJacob Faibussowitsch         ierr_petsc_call_continue_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_continue_, PETSC_ERROR_REPEAT, " "); \
8093ba16761SJacob Faibussowitsch         (void)ierr_petsc_call_continue_; \
8103ba16761SJacob Faibussowitsch       } \
8119566063dSJacob Faibussowitsch     } while (0)
8129566063dSJacob Faibussowitsch #endif
8139566063dSJacob Faibussowitsch 
8149566063dSJacob Faibussowitsch /*MC
8159566063dSJacob Faibussowitsch   CHKERRABORT - Checks error code returned from PETSc function. If non-zero it aborts immediately.
8169566063dSJacob Faibussowitsch 
8179566063dSJacob Faibussowitsch   Synopsis:
8189566063dSJacob Faibussowitsch   #include <petscerror.h>
8199566063dSJacob Faibussowitsch   void CHKERRABORT(MPI_Comm comm, PetscErrorCode ierr)
8209566063dSJacob Faibussowitsch 
8219566063dSJacob Faibussowitsch   Not Collective
8229566063dSJacob Faibussowitsch 
8239566063dSJacob Faibussowitsch   Input Parameters:
8249566063dSJacob Faibussowitsch + comm - the MPI communicator
8259566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
8269566063dSJacob Faibussowitsch 
8279566063dSJacob Faibussowitsch   Level: deprecated
8289566063dSJacob Faibussowitsch 
829667f096bSBarry Smith   Note:
830667f096bSBarry Smith   Deprecated in favor of `PetscCallAbort()`. This routine behaves identically to it.
831667f096bSBarry Smith 
832af27ebaaSBarry Smith .seealso: `PetscCallAbort()`, `PetscErrorCode`
8339566063dSJacob Faibussowitsch M*/
8349566063dSJacob Faibussowitsch #define CHKERRABORT(comm, ...) PetscCallAbort(comm, __VA_ARGS__)
8359566063dSJacob Faibussowitsch #define CHKERRCONTINUE(...)    PetscCallContinue(__VA_ARGS__)
8369566063dSJacob Faibussowitsch 
8379566063dSJacob Faibussowitsch /*MC
83887497f52SBarry Smith    CHKERRA - Fortran-only replacement for use of `CHKERRQ()` in the main program, which aborts immediately
839f388eb8bSPatrick Sanan 
840f388eb8bSPatrick Sanan    Synopsis:
841f388eb8bSPatrick Sanan    #include <petscsys.h>
842f388eb8bSPatrick Sanan    PetscErrorCode CHKERRA(PetscErrorCode ierr)
843f388eb8bSPatrick Sanan 
844f388eb8bSPatrick Sanan    Not Collective
845f388eb8bSPatrick Sanan 
8462fe279fdSBarry Smith    Input Parameter:
847f388eb8bSPatrick Sanan .  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
848f388eb8bSPatrick Sanan 
84987497f52SBarry Smith    Level: deprecated
850f388eb8bSPatrick Sanan 
85187497f52SBarry Smith    Note:
85287497f52SBarry Smith    This macro is rarely needed, normal usage is `PetscCallA()` in the main Fortran program.
853f388eb8bSPatrick Sanan 
854989c49a2SBarry Smith    Developer Note:
855989c49a2SBarry Smith    Why isn't this named `CHKERRABORT()` in Fortran?
856989c49a2SBarry Smith 
85787497f52SBarry Smith .seealso: `PetscCall()`, `PetscCallA()`, `PetscCallAbort()`, `CHKERRQ()`, `SETERRA()`, `SETERRQ()`, `SETERRABORT()`
858f388eb8bSPatrick Sanan M*/
859f388eb8bSPatrick Sanan 
8601e4f893aSSatish Balay PETSC_EXTERN PetscBool petscwaitonerrorflg;
8611e4f893aSSatish Balay PETSC_EXTERN PetscBool petscindebugger;
86235f00c14SToby Isaac PETSC_EXTERN PetscBool petscabortmpifinalize;
8637c66cc67SJunchao Zhang 
8647c66cc67SJunchao Zhang /*MC
865667f096bSBarry Smith    PETSCABORT - Call `MPI_Abort()` with an informative error code
8667c66cc67SJunchao Zhang 
8677c66cc67SJunchao Zhang    Synopsis:
8687c66cc67SJunchao Zhang    #include <petscsys.h>
8697c66cc67SJunchao Zhang    PETSCABORT(MPI_Comm comm, PetscErrorCode ierr)
8707c66cc67SJunchao Zhang 
8717cdbe19fSJose E. Roman    Collective; No Fortran Support
8727c66cc67SJunchao Zhang 
8737c66cc67SJunchao Zhang    Input Parameters:
87495bd0b28SBarry Smith +  comm - An MPI communicator, so that the error can be collective
8757c66cc67SJunchao Zhang -  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
8767c66cc67SJunchao Zhang 
877bf4d2887SBarry Smith    Level: advanced
8787c66cc67SJunchao Zhang 
879bf4d2887SBarry Smith    Notes:
880989c49a2SBarry Smith    If the option `-start_in_debugger` was used then this calls `abort()` to stop the program in the debugger.
881bf4d2887SBarry Smith 
882989c49a2SBarry Smith    if `PetscCIEnabledPortableErrorOutput` is set, which means the code is running in the PETSc test harness (make test),
883989c49a2SBarry Smith    and `comm` is `MPI_COMM_WORLD` it strives to exit cleanly without calling `MPI_Abort()` and instead calling `MPI_Finalize()`.
884660278c0SBarry Smith 
885989c49a2SBarry Smith    This is currently only used when an error propagates up to the C `main()` program and is detected by a `PetscCall()`, `PetscCallMPI()`,
886989c49a2SBarry Smith    or is set in `main()` with `SETERRQ()`. Abort calls such as `SETERRABORT()`,
887989c49a2SBarry Smith    `PetscCheckAbort()`, `PetscCallMPIAbort()`, and `PetscCallAbort()` always call `MPI_Abort()` and do not have any special
888989c49a2SBarry Smith    handling for the test harness.
889989c49a2SBarry Smith 
890989c49a2SBarry Smith    Developer Note:
891989c49a2SBarry Smith    Should the other abort calls also pass through this call instead of calling `MPI_Abort()` directly?
892989c49a2SBarry Smith 
893989c49a2SBarry Smith .seealso: `PetscError()`, `PetscCall()`, `SETERRABORT()`, `PetscCheckAbort()`, `PetscCallMPIAbort()`, `PetscCall()`, `PetscCallMPI()`,
894af27ebaaSBarry Smith           `PetscCallAbort()`, `MPI_Abort()`, `PetscErrorCode`
8957c66cc67SJunchao Zhang M*/
89629f88feaSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
89729f88feaSJacob Faibussowitsch void PETSCABORT(MPI_Comm, PetscErrorCode);
89829f88feaSJacob Faibussowitsch #else
8999371c9d4SSatish Balay   #define PETSCABORT(comm, ...) \
9009371c9d4SSatish Balay     do { \
9013ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_abort_; \
9023ba16761SJacob Faibussowitsch       if (petscwaitonerrorflg) { ierr_petsc_abort_ = PetscSleep(1000); } \
9033ba16761SJacob Faibussowitsch       if (petscindebugger) { \
9043ba16761SJacob Faibussowitsch         abort(); \
9053ba16761SJacob Faibussowitsch       } else { \
9067da6b3ccSLisandro Dalcin         PetscMPIInt size_; \
9073ba16761SJacob Faibussowitsch         ierr_petsc_abort_ = __VA_ARGS__; \
9087da6b3ccSLisandro Dalcin         MPI_Comm_size(comm, &size_); \
90935f00c14SToby Isaac         if (PetscCIEnabledPortableErrorOutput && (size_ == PetscGlobalSize || petscabortmpifinalize) && ierr_petsc_abort_ != PETSC_ERR_SIG) { \
9109371c9d4SSatish Balay           MPI_Finalize(); \
9119371c9d4SSatish Balay           exit(0); \
912660278c0SBarry Smith         } else if (PetscCIEnabledPortableErrorOutput && PetscGlobalSize == 1) { \
913660278c0SBarry Smith           exit(0); \
914660278c0SBarry Smith         } else { \
915660278c0SBarry Smith           MPI_Abort(comm, (PetscMPIInt)ierr_petsc_abort_); \
916660278c0SBarry Smith         } \
9173fcd9f07SJacob Faibussowitsch       } \
9187c66cc67SJunchao Zhang     } while (0)
91929f88feaSJacob Faibussowitsch #endif
920986eef2eSBarry Smith 
9219566063dSJacob Faibussowitsch #ifdef PETSC_CLANGUAGE_CXX
922986eef2eSBarry Smith   /*MC
9239566063dSJacob Faibussowitsch   PetscCallThrow - Checks error code, if non-zero it calls the C++ error handler which throws
9249566063dSJacob Faibussowitsch   an exception
925986eef2eSBarry Smith 
926986eef2eSBarry Smith   Synopsis:
9279566063dSJacob Faibussowitsch   #include <petscerror.h>
9289566063dSJacob Faibussowitsch   void PetscCallThrow(PetscErrorCode ierr)
929986eef2eSBarry Smith 
930986eef2eSBarry Smith   Not Collective
931986eef2eSBarry Smith 
9329566063dSJacob Faibussowitsch   Input Parameter:
933986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
934986eef2eSBarry Smith 
935667f096bSBarry Smith   Level: beginner
936667f096bSBarry Smith 
937986eef2eSBarry Smith   Notes:
938af27ebaaSBarry Smith   Requires PETSc to be configured with clanguage of c++. Throws a std::runtime_error() on error.
939986eef2eSBarry Smith 
94087497f52SBarry Smith   Once the error handler throws the exception you can use `PetscCallVoid()` which returns without
94187497f52SBarry Smith   an error code (bad idea since the error is ignored) or `PetscCallAbort()` to have `MPI_Abort()`
9429566063dSJacob Faibussowitsch   called immediately.
9439566063dSJacob Faibussowitsch 
944db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscCall()`, `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`,
945db781477SPatrick Sanan           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`
946986eef2eSBarry Smith M*/
9479371c9d4SSatish Balay   #define PetscCallThrow(...) \
9489371c9d4SSatish Balay     do { \
9493ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
9503ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_throw_ = __VA_ARGS__; \
9513ba16761SJacob 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); \
952ffc4695bSBarry Smith     } while (0)
95385614651SBarry Smith 
954cc26af49SMatthew Knepley   /*MC
955cc26af49SMatthew Knepley   CHKERRXX - Checks error code, if non-zero it calls the C++ error handler which throws an exception
956cc26af49SMatthew Knepley 
957cc26af49SMatthew Knepley   Synopsis:
9589566063dSJacob Faibussowitsch   #include <petscerror.h>
9593af045c5SBarry Smith   void CHKERRXX(PetscErrorCode ierr)
960cc26af49SMatthew Knepley 
961eca87e8dSBarry Smith   Not Collective
962cc26af49SMatthew Knepley 
9639566063dSJacob Faibussowitsch   Input Parameter:
9643af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
965cc26af49SMatthew Knepley 
9669566063dSJacob Faibussowitsch   Level: deprecated
967cc26af49SMatthew Knepley 
968667f096bSBarry Smith   Note:
969667f096bSBarry Smith   Deprecated in favor of `PetscCallThrow()`. This routine behaves identically to it.
970667f096bSBarry Smith 
971db781477SPatrick Sanan .seealso: `PetscCallThrow()`
972cc26af49SMatthew Knepley M*/
9739566063dSJacob Faibussowitsch   #define CHKERRXX(...) PetscCallThrow(__VA_ARGS__)
974fd705b32SMatthew Knepley #endif
975fd705b32SMatthew Knepley 
9763ba16761SJacob Faibussowitsch #define PetscCallCXX_Private(__SETERR_FUNC__, __COMM__, ...) \
9773ba16761SJacob Faibussowitsch   do { \
9783ba16761SJacob Faibussowitsch     PetscStackUpdateLine; \
9793ba16761SJacob Faibussowitsch     try { \
9803ba16761SJacob Faibussowitsch       __VA_ARGS__; \
9813ba16761SJacob Faibussowitsch     } catch (const std::exception &e) { \
9823ba16761SJacob Faibussowitsch       __SETERR_FUNC__(__COMM__, PETSC_ERR_LIB, "%s", e.what()); \
9833ba16761SJacob Faibussowitsch     } \
9843ba16761SJacob Faibussowitsch   } while (0)
9853ba16761SJacob Faibussowitsch 
9863f520e80SJunchao Zhang /*MC
9879566063dSJacob Faibussowitsch   PetscCallCXX - Checks C++ function calls and if they throw an exception, catch it and then
9889566063dSJacob Faibussowitsch   return a PETSc error code
9893f520e80SJunchao Zhang 
9903f520e80SJunchao Zhang   Synopsis:
9919566063dSJacob Faibussowitsch   #include <petscerror.h>
9925687f061SJacob Faibussowitsch   void PetscCallCXX(...) noexcept;
9933f520e80SJunchao Zhang 
9943f520e80SJunchao Zhang   Not Collective
9953f520e80SJunchao Zhang 
9969566063dSJacob Faibussowitsch   Input Parameter:
9975687f061SJacob Faibussowitsch . __VA_ARGS__ - An arbitrary expression
9985687f061SJacob Faibussowitsch 
9995687f061SJacob Faibussowitsch   Level: beginner
10003f520e80SJunchao Zhang 
10013f520e80SJunchao Zhang   Notes:
10025687f061SJacob Faibussowitsch   `PetscCallCXX(...)` is a macro replacement for
10039566063dSJacob Faibussowitsch .vb
10049566063dSJacob Faibussowitsch   try {
10055687f061SJacob Faibussowitsch     __VA_ARGS__;
10069566063dSJacob Faibussowitsch   } catch (const std::exception& e) {
10079566063dSJacob Faibussowitsch     return ConvertToPetscErrorCode(e);
10089566063dSJacob Faibussowitsch   }
10099566063dSJacob Faibussowitsch .ve
10109566063dSJacob Faibussowitsch   Due to the fact that it catches any (reasonable) exception, it is essentially noexcept.
10113f520e80SJunchao Zhang 
10125687f061SJacob Faibussowitsch   If you cannot return a `PetscErrorCode` use `PetscCallCXXAbort()` instead.
10135687f061SJacob Faibussowitsch 
10149566063dSJacob Faibussowitsch   Example Usage:
10159566063dSJacob Faibussowitsch .vb
10169566063dSJacob Faibussowitsch   void foo(void) { throw std::runtime_error("error"); }
10173f520e80SJunchao Zhang 
10189566063dSJacob Faibussowitsch   void bar()
10199566063dSJacob Faibussowitsch   {
1020e8952933SJacob Faibussowitsch     PetscCallCXX(foo()); // ERROR bar() does not return PetscErrorCode
10219566063dSJacob Faibussowitsch   }
10229566063dSJacob Faibussowitsch 
10239566063dSJacob Faibussowitsch   PetscErrorCode baz()
10249566063dSJacob Faibussowitsch   {
10259566063dSJacob Faibussowitsch     PetscCallCXX(foo()); // OK
10269566063dSJacob Faibussowitsch 
10279566063dSJacob Faibussowitsch     PetscCallCXX(
10289566063dSJacob Faibussowitsch       bar();
1029d5b43468SJose E. Roman       foo(); // OK multiple statements allowed
10309566063dSJacob Faibussowitsch     );
10319566063dSJacob Faibussowitsch   }
10329566063dSJacob Faibussowitsch 
10339566063dSJacob Faibussowitsch   struct bop
10349566063dSJacob Faibussowitsch   {
10359566063dSJacob Faibussowitsch     bop()
10369566063dSJacob Faibussowitsch     {
1037e8952933SJacob Faibussowitsch       PetscCallCXX(foo()); // ERROR returns PetscErrorCode, cannot be used in constructors
10389566063dSJacob Faibussowitsch     }
10399566063dSJacob Faibussowitsch   };
10409566063dSJacob Faibussowitsch 
1041e8952933SJacob Faibussowitsch   // ERROR contains do-while, cannot be used as function-try block
10429566063dSJacob Faibussowitsch   PetscErrorCode qux() PetscCallCXX(
10439566063dSJacob Faibussowitsch     bar();
10449566063dSJacob Faibussowitsch     baz();
10459566063dSJacob Faibussowitsch     foo();
10469566063dSJacob Faibussowitsch     return 0;
10479566063dSJacob Faibussowitsch   )
10489566063dSJacob Faibussowitsch .ve
10499566063dSJacob Faibussowitsch 
10505687f061SJacob Faibussowitsch .seealso: `PetscCallCXXAbort()`, `PetscCallThrow()`, `SETERRQ()`, `PetscCall()`,
10515687f061SJacob Faibussowitsch           `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
10525687f061SJacob Faibussowitsch           `PetscError()`, `CHKMEMQ`
10533f520e80SJunchao Zhang M*/
10543ba16761SJacob Faibussowitsch #define PetscCallCXX(...) PetscCallCXX_Private(SETERRQ, PETSC_COMM_SELF, __VA_ARGS__)
10555687f061SJacob Faibussowitsch 
10565687f061SJacob Faibussowitsch /*MC
10575687f061SJacob Faibussowitsch   PetscCallCXXAbort - Like `PetscCallCXX()` but calls `MPI_Abort()` instead of returning an
10585687f061SJacob Faibussowitsch   error-code
10595687f061SJacob Faibussowitsch 
10605687f061SJacob Faibussowitsch   Synopsis:
10615687f061SJacob Faibussowitsch   #include <petscerror.h>
10625687f061SJacob Faibussowitsch   void PetscCallCXXAbort(MPI_Comm comm, ...) noexcept;
10635687f061SJacob Faibussowitsch 
106420f4b53cSBarry Smith   Collective; No Fortran Support
10655687f061SJacob Faibussowitsch 
10665687f061SJacob Faibussowitsch   Input Parameters:
10675687f061SJacob Faibussowitsch + comm        - The MPI communicator to abort on
10685687f061SJacob Faibussowitsch - __VA_ARGS__ - An arbitrary expression
10695687f061SJacob Faibussowitsch 
10705687f061SJacob Faibussowitsch   Level: beginner
10715687f061SJacob Faibussowitsch 
10725687f061SJacob Faibussowitsch   Notes:
10735687f061SJacob Faibussowitsch   This macro may be used to check C++ expressions for exceptions in cases where you cannot
10745687f061SJacob Faibussowitsch   return an error code. This includes constructors, destructors, copy/move assignment functions
10755687f061SJacob Faibussowitsch   or constructors among others.
10765687f061SJacob Faibussowitsch 
10775687f061SJacob Faibussowitsch   If an exception is caught, the macro calls `SETERRABORT()` on `comm`. The exception must
10785687f061SJacob Faibussowitsch   derive from `std::exception` in order to be caught.
10795687f061SJacob Faibussowitsch 
10805687f061SJacob Faibussowitsch   If the routine _can_ return an error-code it is highly advised to use `PetscCallCXX()`
10815687f061SJacob Faibussowitsch   instead.
10825687f061SJacob Faibussowitsch 
10835687f061SJacob Faibussowitsch   See `PetscCallCXX()` for additional discussion.
10845687f061SJacob Faibussowitsch 
10855687f061SJacob Faibussowitsch   Example Usage:
10865687f061SJacob Faibussowitsch .vb
10875687f061SJacob Faibussowitsch   class Foo
10885687f061SJacob Faibussowitsch   {
10895687f061SJacob Faibussowitsch     std::vector<int> data_;
10905687f061SJacob Faibussowitsch 
10915687f061SJacob Faibussowitsch   public:
10925687f061SJacob Faibussowitsch     // normally std::vector::reserve() may raise an exception, but since we handle it with
10935687f061SJacob Faibussowitsch     // PetscCallCXXAbort() we may mark this routine as noexcept!
10945687f061SJacob Faibussowitsch     Foo() noexcept
10955687f061SJacob Faibussowitsch     {
10965687f061SJacob Faibussowitsch       PetscCallCXXAbort(PETSC_COMM_SELF, data_.reserve(10));
10975687f061SJacob Faibussowitsch     }
10985687f061SJacob Faibussowitsch   };
10995687f061SJacob Faibussowitsch 
11005687f061SJacob Faibussowitsch   std::vector<int> bar()
11015687f061SJacob Faibussowitsch   {
11025687f061SJacob Faibussowitsch     std::vector<int> v;
11035687f061SJacob Faibussowitsch 
11045687f061SJacob Faibussowitsch     PetscFunctionBegin;
11055687f061SJacob Faibussowitsch     // OK!
11065687f061SJacob Faibussowitsch     PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1));
11075687f061SJacob Faibussowitsch     PetscFunctionReturn(v);
11085687f061SJacob Faibussowitsch   }
11095687f061SJacob Faibussowitsch 
11105687f061SJacob Faibussowitsch   PetscErrorCode baz()
11115687f061SJacob Faibussowitsch   {
11125687f061SJacob Faibussowitsch     std::vector<int> v;
11135687f061SJacob Faibussowitsch 
11145687f061SJacob Faibussowitsch     PetscFunctionBegin;
11155687f061SJacob Faibussowitsch     // WRONG! baz() returns a PetscErrorCode, prefer PetscCallCXX() instead
11165687f061SJacob Faibussowitsch     PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1));
11173ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
11185687f061SJacob Faibussowitsch   }
11195687f061SJacob Faibussowitsch .ve
11205687f061SJacob Faibussowitsch 
11215687f061SJacob Faibussowitsch .seealso: `PetscCallCXX()`, `SETERRABORT()`, `PetscCallAbort()`
11225687f061SJacob Faibussowitsch M*/
11233ba16761SJacob Faibussowitsch #define PetscCallCXXAbort(comm, ...) PetscCallCXX_Private(SETERRABORT, comm, __VA_ARGS__)
11243f520e80SJunchao Zhang 
112530de9b25SBarry Smith /*MC
11269566063dSJacob Faibussowitsch   CHKERRCXX - Checks C++ function calls and if they throw an exception, catch it and then
11279566063dSJacob Faibussowitsch   return a PETSc error code
11289566063dSJacob Faibussowitsch 
11299566063dSJacob Faibussowitsch   Synopsis:
11309566063dSJacob Faibussowitsch   #include <petscerror.h>
11319566063dSJacob Faibussowitsch   void CHKERRCXX(func) noexcept;
11329566063dSJacob Faibussowitsch 
11339566063dSJacob Faibussowitsch   Not Collective
11349566063dSJacob Faibussowitsch 
11359566063dSJacob Faibussowitsch   Input Parameter:
11369566063dSJacob Faibussowitsch . func - C++ function calls
11379566063dSJacob Faibussowitsch 
11389566063dSJacob Faibussowitsch   Level: deprecated
11399566063dSJacob Faibussowitsch 
1140667f096bSBarry Smith   Note:
1141667f096bSBarry Smith   Deprecated in favor of `PetscCallCXX()`. This routine behaves identically to it.
1142667f096bSBarry Smith 
1143db781477SPatrick Sanan .seealso: `PetscCallCXX()`
11449566063dSJacob Faibussowitsch M*/
11459566063dSJacob Faibussowitsch #define CHKERRCXX(...) PetscCallCXX(__VA_ARGS__)
11469566063dSJacob Faibussowitsch 
11479566063dSJacob Faibussowitsch /*MC
114830de9b25SBarry Smith    CHKMEMQ - Checks the memory for corruption, calls error handler if any is detected
114930de9b25SBarry Smith 
115030de9b25SBarry Smith    Synopsis:
1151aaa7dc30SBarry Smith    #include <petscsys.h>
115291d3bdf4SKris Buschelman    CHKMEMQ;
115330de9b25SBarry Smith 
1154eca87e8dSBarry Smith    Not Collective
1155eca87e8dSBarry Smith 
115630de9b25SBarry Smith    Level: beginner
115730de9b25SBarry Smith 
115830de9b25SBarry Smith    Notes:
1159af27ebaaSBarry Smith    We recommend using Valgrind <https://petsc.org/release/faq/#valgrind> or for NVIDIA CUDA systems
1160af27ebaaSBarry Smith    <https://docs.nvidia.com/cuda/cuda-memcheck/index.html> for finding memory problems. The ``CHKMEMQ`` macro is useful on systems that
11615ed36255SBarry Smith    do not have valgrind, but is not as good as valgrind or cuda-memcheck.
11621957e957SBarry Smith 
1163667f096bSBarry Smith    Must run with the option `-malloc_debug` (`-malloc_test` in debug mode; or if `PetscMallocSetDebug()` called) to enable this option
116430de9b25SBarry Smith 
116530de9b25SBarry Smith    Once the error handler is called the calling function is then returned from with the given error code.
116630de9b25SBarry Smith 
116730de9b25SBarry Smith    By defaults prints location where memory that is corrupted was allocated.
116830de9b25SBarry Smith 
11698af9f190SBarry Smith    Use `CHKMEMA` for functions that return `void`
1170f621e05eSBarry Smith 
1171db781477SPatrick Sanan .seealso: `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `SETERRQ()`, `PetscMallocValidate()`
117230de9b25SBarry Smith M*/
11736d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
1174064a246eSJacob Faibussowitsch   #define CHKMEMQ
1175064a246eSJacob Faibussowitsch   #define CHKMEMA
11766d210af2SJacob Faibussowitsch #else
11779371c9d4SSatish Balay   #define CHKMEMQ \
11789371c9d4SSatish Balay     do { \
11793ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_memq_ = PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__); \
11803ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_memq_ != PETSC_SUCCESS)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_memq_, PETSC_ERROR_REPEAT, " "); \
118186d09637SLisandro Dalcin     } while (0)
11826d210af2SJacob Faibussowitsch   #define CHKMEMA PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__)
1183064a246eSJacob Faibussowitsch #endif
11849566063dSJacob Faibussowitsch 
1185668f157eSBarry Smith /*E
1186668f157eSBarry Smith   PetscErrorType - passed to the PETSc error handling routines indicating if this is the first or a later call to the error handlers
1187668f157eSBarry Smith 
1188668f157eSBarry Smith   Level: advanced
1189668f157eSBarry Smith 
1190667f096bSBarry Smith   Note:
119187497f52SBarry Smith   `PETSC_ERROR_IN_CXX` indicates the error was detected in C++ and an exception should be generated
1192d736bfebSBarry Smith 
1193667f096bSBarry Smith   Developer Note:
1194667f096bSBarry Smith   This is currently used to decide when to print the detailed information about the run in `PetscTraceBackErrorHandler()`
1195668f157eSBarry Smith 
119687497f52SBarry Smith .seealso: `PetscError()`, `SETERRQ()`
1197668f157eSBarry Smith E*/
11989371c9d4SSatish Balay typedef enum {
11999371c9d4SSatish Balay   PETSC_ERROR_INITIAL = 0,
12009371c9d4SSatish Balay   PETSC_ERROR_REPEAT  = 1,
12019371c9d4SSatish Balay   PETSC_ERROR_IN_CXX  = 2
12029371c9d4SSatish Balay } PetscErrorType;
12034b209cf6SBarry Smith 
1204eb9e708aSLisandro Dalcin #if defined(__clang_analyzer__)
1205eb9e708aSLisandro Dalcin __attribute__((analyzer_noreturn))
1206eb9e708aSLisandro Dalcin #endif
12079371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode
12089371c9d4SSatish Balay PetscError(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, ...) PETSC_ATTRIBUTE_COLD PETSC_ATTRIBUTE_FORMAT(7, 8);
1209eb9e708aSLisandro Dalcin 
1210014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscErrorPrintfInitialize(void);
12113ba16761SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscErrorMessage(PetscErrorCode, const char *[], char **);
1212d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscTraceBackErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1213d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscIgnoreErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1214d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscEmacsClientErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1215d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscMPIAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1216d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1217d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAttachDebuggerErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1218d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscReturnErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1219efca3c55SSatish Balay PETSC_EXTERN PetscErrorCode PetscPushErrorHandler(PetscErrorCode (*handler)(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *), void *);
1220014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopErrorHandler(void);
12218d359177SBarry Smith PETSC_EXTERN PetscErrorCode PetscSignalHandlerDefault(int, void *);
1222014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPushSignalHandler(PetscErrorCode (*)(int, void *), void *);
1223014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopSignalHandler(void);
122428559dc8SJed Brown PETSC_EXTERN PetscErrorCode PetscCheckPointerSetIntensity(PetscInt);
1225c2a741eeSJunchao Zhang PETSC_EXTERN void           PetscSignalSegvCheckPointerOrMpi(void);
1226edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 13, 0, "PetscSignalSegvCheckPointerOrMpi()", ) static inline void PetscSignalSegvCheckPointer(void)
1227d71ae5a4SJacob Faibussowitsch {
12289371c9d4SSatish Balay   PetscSignalSegvCheckPointerOrMpi();
12299371c9d4SSatish Balay }
1230329f5518SBarry Smith 
1231639ff905SBarry Smith /*MC
1232639ff905SBarry Smith     PetscErrorPrintf - Prints error messages.
1233639ff905SBarry Smith 
1234639ff905SBarry Smith    Synopsis:
1235aaa7dc30SBarry Smith     #include <petscsys.h>
1236639ff905SBarry Smith      PetscErrorCode (*PetscErrorPrintf)(const char format[], ...);
1237639ff905SBarry Smith 
12387cdbe19fSJose E. Roman     Not Collective; No Fortran Support
12397cdbe19fSJose E. Roman 
1240f899ff85SJose E. Roman     Input Parameter:
1241cd05f99aSJacob Faibussowitsch .   format - the usual `printf()` format string
1242639ff905SBarry Smith 
1243639ff905SBarry Smith    Options Database Keys:
12441957e957SBarry Smith +  -error_output_stdout - cause error messages to be printed to stdout instead of the (default) stderr
1245e1bc860dSBarry Smith -  -error_output_none   - to turn off all printing of error messages (does not change the way the error is handled.)
1246639ff905SBarry Smith 
1247cf53795eSBarry Smith    Level: developer
1248cf53795eSBarry Smith 
124995452b02SPatrick Sanan    Notes:
125095452b02SPatrick Sanan    Use
1251667f096bSBarry Smith .vb
125295bd0b28SBarry Smith      PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the error is handled) and
1253667f096bSBarry Smith      PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on or you can use your own function
1254667f096bSBarry Smith .ve
1255639ff905SBarry Smith    Use
1256667f096bSBarry Smith .vb
125787497f52SBarry Smith      `PETSC_STDERR` = FILE* obtained from a file open etc. to have stderr printed to the file.
125887497f52SBarry Smith      `PETSC_STDOUT` = FILE* obtained from a file open etc. to have stdout printed to the file.
1259667f096bSBarry Smith .ve
1260639ff905SBarry Smith    Use
1261af27ebaaSBarry Smith .vb
126287497f52SBarry Smith       `PetscPushErrorHandler()` to provide your own error handler that determines what kind of messages to print
1263af27ebaaSBarry Smith .ve
1264639ff905SBarry Smith 
1265db781477SPatrick Sanan .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscHelpPrintf()`, `PetscPrintf()`, `PetscPushErrorHandler()`, `PetscVFPrintf()`, `PetscHelpPrintf()`
1266639ff905SBarry Smith M*/
12673ca90d2dSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode (*PetscErrorPrintf)(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1268639ff905SBarry Smith 
1269cf0818bdSBarry Smith /*E
1270cf0818bdSBarry Smith      PetscFPTrap - types of floating point exceptions that may be trapped
1271cf0818bdSBarry Smith 
127287497f52SBarry Smith      Currently only `PETSC_FP_TRAP_OFF` and `PETSC_FP_TRAP_ON` are handled. All others are treated as `PETSC_FP_TRAP_ON`.
1273cf0818bdSBarry Smith 
1274cf0818bdSBarry Smith      Level: intermediate
1275cf0818bdSBarry Smith 
12767de69702SBarry Smith .seealso: `PetscSetFPTrap()`, `PetscFPTrapPush()`
1277cf0818bdSBarry Smith  E*/
12789371c9d4SSatish Balay typedef enum {
12799371c9d4SSatish Balay   PETSC_FP_TRAP_OFF      = 0,
12809371c9d4SSatish Balay   PETSC_FP_TRAP_INDIV    = 1,
12819371c9d4SSatish Balay   PETSC_FP_TRAP_FLTOPERR = 2,
12829371c9d4SSatish Balay   PETSC_FP_TRAP_FLTOVF   = 4,
12839371c9d4SSatish Balay   PETSC_FP_TRAP_FLTUND   = 8,
12849371c9d4SSatish Balay   PETSC_FP_TRAP_FLTDIV   = 16,
12859371c9d4SSatish Balay   PETSC_FP_TRAP_FLTINEX  = 32
12869371c9d4SSatish Balay } PetscFPTrap;
1287bd2b07b1SBarry 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)
1288014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSetFPTrap(PetscFPTrap);
1289014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPush(PetscFPTrap);
1290014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPop(void);
1291aba4c478SBarry Smith PETSC_EXTERN PetscErrorCode PetscDetermineInitialFPTrap(void);
129254a8ef01SBarry Smith 
12933a40ed3dSBarry Smith /*
12943a40ed3dSBarry Smith       Allows the code to build a stack frame as it runs
12953a40ed3dSBarry Smith */
12963a40ed3dSBarry Smith 
129799cd645aSJed Brown #define PETSCSTACKSIZE 64
12983a40ed3dSBarry Smith typedef struct {
12990e33f6ddSBarry Smith   const char *function[PETSCSTACKSIZE];
13000e33f6ddSBarry Smith   const char *file[PETSCSTACKSIZE];
1301184914b5SBarry Smith   int         line[PETSCSTACKSIZE];
1302362febeeSStefano Zampini   int         petscroutine[PETSCSTACKSIZE]; /* 0 external called from petsc, 1 petsc functions, 2 petsc user functions */
1303184914b5SBarry Smith   int         currentsize;
1304a2f94806SJed Brown   int         hotdepth;
13054be741a6SBarry Smith   PetscBool   check; /* option to check for correct Push/Pop semantics, true for default petscstack but not other stacks */
13063a40ed3dSBarry Smith } PetscStack;
1307dfb7d7afSStefano Zampini #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
130827104ee2SJacob Faibussowitsch PETSC_EXTERN PetscStack petscstack;
130927104ee2SJacob Faibussowitsch #endif
13103a40ed3dSBarry Smith 
13115d12eec7SSatish Balay #if defined(PETSC_SERIALIZE_FUNCTIONS)
13125d12eec7SSatish Balay   #include <petsc/private/petscfptimpl.h>
13135d12eec7SSatish Balay   /*
13145d12eec7SSatish Balay    Registers the current function into the global function pointer to function name table
13155d12eec7SSatish Balay 
13165d12eec7SSatish Balay    Have to fix this to handle errors but cannot return error since used in PETSC_VIEWER_DRAW_() etc
13175d12eec7SSatish Balay */
13189371c9d4SSatish Balay   #define PetscRegister__FUNCT__() \
13199371c9d4SSatish Balay     do { \
13205d12eec7SSatish Balay       static PetscBool __chked = PETSC_FALSE; \
13215d12eec7SSatish Balay       if (!__chked) { \
13229371c9d4SSatish Balay         void *ptr; \
13233ba16761SJacob Faibussowitsch         PetscCallAbort(PETSC_COMM_SELF, PetscDLSym(NULL, PETSC_FUNCTION_NAME, &ptr)); \
13245d12eec7SSatish Balay         __chked = PETSC_TRUE; \
13259371c9d4SSatish Balay       } \
13269371c9d4SSatish Balay     } while (0)
13275d12eec7SSatish Balay #else
13285d12eec7SSatish Balay   #define PetscRegister__FUNCT__()
13295d12eec7SSatish Balay #endif
13305d12eec7SSatish Balay 
1331eae3dc7dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) || defined(__clang_analyzer__)
13326d210af2SJacob Faibussowitsch   #define PetscStackPushNoCheck(funct, petsc_routine, hot)
133337154d10SBarry Smith   #define PetscStackUpdateLine
1334792fecdfSBarry Smith   #define PetscStackPushExternal(funct)
13356d210af2SJacob Faibussowitsch   #define PetscStackPopNoCheck
13366d210af2SJacob Faibussowitsch   #define PetscStackClearTop
13376d210af2SJacob Faibussowitsch   #define PetscFunctionBegin
13386d210af2SJacob Faibussowitsch   #define PetscFunctionBeginUser
13396d210af2SJacob Faibussowitsch   #define PetscFunctionBeginHot
13400a57284eSJacob Faibussowitsch   #define PetscFunctionReturn(...)  return __VA_ARGS__
13416d210af2SJacob Faibussowitsch   #define PetscFunctionReturnVoid() return
13426d210af2SJacob Faibussowitsch   #define PetscStackPop
13436d210af2SJacob Faibussowitsch   #define PetscStackPush(f)
1344dfb7d7afSStefano Zampini #elif defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1345660278c0SBarry Smith 
13469371c9d4SSatish Balay   #define PetscStackPush_Private(stack__, file__, func__, line__, petsc_routine__, hot__) \
13479371c9d4SSatish Balay     do { \
13485a96b57dSJacob Faibussowitsch       if (stack__.currentsize < PETSCSTACKSIZE) { \
13495a96b57dSJacob Faibussowitsch         stack__.function[stack__.currentsize] = func__; \
1350ef1023bdSBarry Smith         if (petsc_routine__) { \
1351ef1023bdSBarry Smith           stack__.file[stack__.currentsize] = file__; \
13525a96b57dSJacob Faibussowitsch           stack__.line[stack__.currentsize] = line__; \
1353ef1023bdSBarry Smith         } else { \
1354648bc8c4SBarry Smith           stack__.file[stack__.currentsize] = PETSC_NULLPTR; \
1355ef1023bdSBarry Smith           stack__.line[stack__.currentsize] = 0; \
1356ef1023bdSBarry Smith         } \
13575a96b57dSJacob Faibussowitsch         stack__.petscroutine[stack__.currentsize] = petsc_routine__; \
13585a96b57dSJacob Faibussowitsch       } \
13595a96b57dSJacob Faibussowitsch       ++stack__.currentsize; \
13605a96b57dSJacob Faibussowitsch       stack__.hotdepth += (hot__ || stack__.hotdepth); \
13615a96b57dSJacob Faibussowitsch     } while (0)
13625a96b57dSJacob Faibussowitsch 
13634be741a6SBarry Smith   /* uses PetscCheckAbort() because may be used in a function that does not return an error code */
13649371c9d4SSatish Balay   #define PetscStackPop_Private(stack__, func__) \
13659371c9d4SSatish Balay     do { \
13664be741a6SBarry 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__); \
13675a96b57dSJacob Faibussowitsch       if (--stack__.currentsize < PETSCSTACKSIZE) { \
13689371c9d4SSatish 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", \
13699371c9d4SSatish Balay                         stack__.function[stack__.currentsize], stack__.file[stack__.currentsize], stack__.line[stack__.currentsize], func__, __FILE__, __LINE__); \
13705a96b57dSJacob Faibussowitsch         stack__.function[stack__.currentsize]     = PETSC_NULLPTR; \
13715a96b57dSJacob Faibussowitsch         stack__.file[stack__.currentsize]         = PETSC_NULLPTR; \
13725a96b57dSJacob Faibussowitsch         stack__.line[stack__.currentsize]         = 0; \
13735a96b57dSJacob Faibussowitsch         stack__.petscroutine[stack__.currentsize] = 0; \
13745a96b57dSJacob Faibussowitsch       } \
13755a96b57dSJacob Faibussowitsch       stack__.hotdepth = PetscMax(stack__.hotdepth - 1, 0); \
13765a96b57dSJacob Faibussowitsch     } while (0)
13775a96b57dSJacob Faibussowitsch 
1378586f9135SBarry Smith   /*MC
1379586f9135SBarry Smith    PetscStackPushNoCheck - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is
1380586f9135SBarry Smith    currently in the source code.
1381586f9135SBarry Smith 
1382586f9135SBarry Smith    Synopsis:
1383586f9135SBarry Smith    #include <petscsys.h>
1384586f9135SBarry Smith    void PetscStackPushNoCheck(char *funct,int petsc_routine,PetscBool hot);
1385586f9135SBarry Smith 
13867cdbe19fSJose E. Roman    Not Collective
13877cdbe19fSJose E. Roman 
1388586f9135SBarry Smith    Input Parameters:
1389586f9135SBarry Smith +  funct - the function name
1390586f9135SBarry Smith .  petsc_routine - 2 user function, 1 PETSc function, 0 some other function
1391586f9135SBarry Smith -  hot - indicates that the function may be called often so expensive error checking should be turned off inside the function
1392586f9135SBarry Smith 
1393586f9135SBarry Smith    Level: developer
1394586f9135SBarry Smith 
1395586f9135SBarry Smith    Notes:
1396586f9135SBarry 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
139787497f52SBarry 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
1398586f9135SBarry Smith    help debug the problem.
1399586f9135SBarry Smith 
1400ef1023bdSBarry Smith    This version does not check the memory corruption (an expensive operation), use `PetscStackPush()` to check the memory.
1401ef1023bdSBarry Smith 
1402792fecdfSBarry Smith    Use `PetscStackPushExternal()` for a function call that is about to be made to a non-PETSc or user function (such as BLAS etc).
1403ef1023bdSBarry Smith 
1404586f9135SBarry Smith    The default stack is a global variable called `petscstack`.
1405586f9135SBarry Smith 
1406586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1407ef1023bdSBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPush()`, `PetscStackPop`,
1408792fecdfSBarry Smith           `PetscStackPushExternal()`
1409586f9135SBarry Smith M*/
14109371c9d4SSatish Balay   #define PetscStackPushNoCheck(funct, petsc_routine, hot) \
14119371c9d4SSatish Balay     do { \
1412e04113cfSBarry Smith       PetscStackSAWsTakeAccess(); \
14135a96b57dSJacob Faibussowitsch       PetscStackPush_Private(petscstack, __FILE__, funct, __LINE__, petsc_routine, hot); \
1414e04113cfSBarry Smith       PetscStackSAWsGrantAccess(); \
1415441dd030SJed Brown     } while (0)
1416441dd030SJed Brown 
1417586f9135SBarry Smith   /*MC
141887497f52SBarry Smith    PetscStackUpdateLine - in a function that has a `PetscFunctionBegin` or `PetscFunctionBeginUser` updates the stack line number to the
141937154d10SBarry Smith    current line number.
142037154d10SBarry Smith 
142137154d10SBarry Smith    Synopsis:
142237154d10SBarry Smith    #include <petscsys.h>
142337154d10SBarry Smith    void PetscStackUpdateLine
142437154d10SBarry Smith 
14257cdbe19fSJose E. Roman    Not Collective
14267cdbe19fSJose E. Roman 
142737154d10SBarry Smith    Level: developer
142837154d10SBarry Smith 
142937154d10SBarry Smith    Notes:
143087497f52SBarry Smith    Using `PetscCall()` and friends automatically handles this process
143187497f52SBarry Smith 
143237154d10SBarry 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
143337154d10SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
143437154d10SBarry Smith    help debug the problem.
143537154d10SBarry Smith 
143695bd0b28SBarry Smith    The default stack is a global variable called `petscstack`.
143737154d10SBarry Smith 
143837154d10SBarry Smith    This is used by `PetscCall()` and is otherwise not like to be needed
143937154d10SBarry Smith 
144037154d10SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`, `PetscCall()`
144137154d10SBarry Smith M*/
144237154d10SBarry Smith   #define PetscStackUpdateLine \
14433ba16761SJacob Faibussowitsch     do { \
1444f1e99121SPierre Jolivet       if (petscstack.currentsize > 0 && petscstack.currentsize < PETSCSTACKSIZE && petscstack.function[petscstack.currentsize - 1] == PETSC_FUNCTION_NAME) { petscstack.line[petscstack.currentsize - 1] = __LINE__; } \
14453ba16761SJacob Faibussowitsch     } while (0)
144637154d10SBarry Smith 
144737154d10SBarry Smith   /*MC
1448792fecdfSBarry Smith    PetscStackPushExternal - Pushes a new function name onto the PETSc default stack that tracks where the running program is
1449ef1023bdSBarry Smith    currently in the source code. Does not include the filename or line number since this is called by the calling routine
1450ef1023bdSBarry Smith    for non-PETSc or user functions.
1451ef1023bdSBarry Smith 
1452ef1023bdSBarry Smith    Synopsis:
1453ef1023bdSBarry Smith    #include <petscsys.h>
1454792fecdfSBarry Smith    void PetscStackPushExternal(char *funct);
1455ef1023bdSBarry Smith 
14567cdbe19fSJose E. Roman    Not Collective
14577cdbe19fSJose E. Roman 
14582fe279fdSBarry Smith    Input Parameter:
1459ef1023bdSBarry Smith .  funct - the function name
1460ef1023bdSBarry Smith 
1461ef1023bdSBarry Smith    Level: developer
1462ef1023bdSBarry Smith 
1463ef1023bdSBarry Smith    Notes:
146487497f52SBarry Smith    Using `PetscCallExternal()` and friends automatically handles this process
146587497f52SBarry Smith 
1466ef1023bdSBarry 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
1467ef1023bdSBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1468ef1023bdSBarry Smith    help debug the problem.
1469ef1023bdSBarry Smith 
1470ef1023bdSBarry Smith    The default stack is a global variable called `petscstack`.
1471ef1023bdSBarry Smith 
1472ef1023bdSBarry Smith    This is to be used when calling an external package function such as a BLAS function.
1473ef1023bdSBarry Smith 
1474ef1023bdSBarry Smith    This also updates the stack line number for the current stack function.
1475ef1023bdSBarry Smith 
1476ef1023bdSBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1477ef1023bdSBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1478ef1023bdSBarry Smith M*/
14799371c9d4SSatish Balay   #define PetscStackPushExternal(funct) \
14809371c9d4SSatish Balay     do { \
14819371c9d4SSatish Balay       PetscStackUpdateLine; \
14829371c9d4SSatish Balay       PetscStackPushNoCheck(funct, 0, PETSC_TRUE); \
1483a8f51744SPierre Jolivet     } while (0)
1484ef1023bdSBarry Smith 
1485ef1023bdSBarry Smith   /*MC
1486586f9135SBarry Smith    PetscStackPopNoCheck - Pops a function name from the PETSc default stack that tracks where the running program is
1487586f9135SBarry Smith    currently in the source code.
1488586f9135SBarry Smith 
1489586f9135SBarry Smith    Synopsis:
1490586f9135SBarry Smith    #include <petscsys.h>
1491586f9135SBarry Smith    void PetscStackPopNoCheck(char *funct);
1492586f9135SBarry Smith 
14937cdbe19fSJose E. Roman    Not Collective
14947cdbe19fSJose E. Roman 
1495586f9135SBarry Smith    Input Parameter:
1496586f9135SBarry Smith .   funct - the function name
1497586f9135SBarry Smith 
1498586f9135SBarry Smith    Level: developer
1499586f9135SBarry Smith 
1500586f9135SBarry Smith    Notes:
150187497f52SBarry Smith    Using `PetscCall()`, `PetscCallExternal()`, `PetscCallBack()` and friends negates the need to call this
150287497f52SBarry Smith 
1503586f9135SBarry 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
1504586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1505586f9135SBarry Smith    help debug the problem.
1506586f9135SBarry Smith 
150795bd0b28SBarry Smith    The default stack is a global variable called `petscstack`.
1508586f9135SBarry Smith 
1509586f9135SBarry Smith    Developer Note:
1510586f9135SBarry Smith    `PetscStackPopNoCheck()` takes a function argument while  `PetscStackPop` does not, this difference is likely just historical.
1511586f9135SBarry Smith 
1512586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1513586f9135SBarry Smith M*/
15149371c9d4SSatish Balay   #define PetscStackPopNoCheck(funct) \
15159371c9d4SSatish Balay     do { \
1516362febeeSStefano Zampini       PetscStackSAWsTakeAccess(); \
15175a96b57dSJacob Faibussowitsch       PetscStackPop_Private(petscstack, funct); \
1518362febeeSStefano Zampini       PetscStackSAWsGrantAccess(); \
1519362febeeSStefano Zampini     } while (0)
1520362febeeSStefano Zampini 
15219371c9d4SSatish Balay   #define PetscStackClearTop \
15229371c9d4SSatish Balay     do { \
1523e04113cfSBarry Smith       PetscStackSAWsTakeAccess(); \
15249371c9d4SSatish Balay       if (petscstack.currentsize > 0 && --petscstack.currentsize < PETSCSTACKSIZE) { \
152527104ee2SJacob Faibussowitsch         petscstack.function[petscstack.currentsize]     = PETSC_NULLPTR; \
152627104ee2SJacob Faibussowitsch         petscstack.file[petscstack.currentsize]         = PETSC_NULLPTR; \
152727104ee2SJacob Faibussowitsch         petscstack.line[petscstack.currentsize]         = 0; \
152827104ee2SJacob Faibussowitsch         petscstack.petscroutine[petscstack.currentsize] = 0; \
1529441dd030SJed Brown       } \
153027104ee2SJacob Faibussowitsch       petscstack.hotdepth = PetscMax(petscstack.hotdepth - 1, 0); \
1531e04113cfSBarry Smith       PetscStackSAWsGrantAccess(); \
1532441dd030SJed Brown     } while (0)
1533441dd030SJed Brown 
153430de9b25SBarry Smith   /*MC
15351957e957SBarry Smith    PetscFunctionBegin - First executable line of each PETSc function,  used for error handling. Final
153687497f52SBarry Smith    line of PETSc functions should be `PetscFunctionReturn`(0);
153730de9b25SBarry Smith 
153830de9b25SBarry Smith    Synopsis:
1539aaa7dc30SBarry Smith    #include <petscsys.h>
154030de9b25SBarry Smith    void PetscFunctionBegin;
154130de9b25SBarry Smith 
154220f4b53cSBarry Smith    Not Collective; No Fortran Support
1543eca87e8dSBarry Smith 
154430de9b25SBarry Smith    Usage:
154530de9b25SBarry Smith .vb
154630de9b25SBarry Smith      int something;
154730de9b25SBarry Smith 
154830de9b25SBarry Smith      PetscFunctionBegin;
154930de9b25SBarry Smith .ve
155030de9b25SBarry Smith 
155130de9b25SBarry Smith    Level: developer
155230de9b25SBarry Smith 
155320f4b53cSBarry Smith    Note:
155420f4b53cSBarry Smith    Use `PetscFunctionBeginUser` for application codes.
155520f4b53cSBarry Smith 
1556586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`
155730de9b25SBarry Smith 
155830de9b25SBarry Smith M*/
15599371c9d4SSatish Balay   #define PetscFunctionBegin \
15609371c9d4SSatish Balay     do { \
1561362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_FALSE); \
1562a2f94806SJed Brown       PetscRegister__FUNCT__(); \
1563a2f94806SJed Brown     } while (0)
1564a2f94806SJed Brown 
1565a2f94806SJed Brown   /*MC
156687497f52SBarry Smith    PetscFunctionBeginHot - Substitute for `PetscFunctionBegin` to be used in functions that are called in
1567a2f94806SJed Brown    performance-critical circumstances.  Use of this function allows for lighter profiling by default.
1568a2f94806SJed Brown 
1569a2f94806SJed Brown    Synopsis:
1570aaa7dc30SBarry Smith    #include <petscsys.h>
1571a2f94806SJed Brown    void PetscFunctionBeginHot;
1572a2f94806SJed Brown 
157320f4b53cSBarry Smith    Not Collective; No Fortran Support
1574a2f94806SJed Brown 
1575a2f94806SJed Brown    Usage:
1576a2f94806SJed Brown .vb
1577a2f94806SJed Brown      int something;
1578a2f94806SJed Brown 
1579a2f94806SJed Brown      PetscFunctionBeginHot;
1580a2f94806SJed Brown .ve
1581a2f94806SJed Brown 
1582a2f94806SJed Brown    Level: developer
1583a2f94806SJed Brown 
1584586f9135SBarry Smith .seealso: `PetscFunctionBegin`, `PetscFunctionReturn()`, `PetscStackPushNoCheck()`
1585a2f94806SJed Brown 
1586a2f94806SJed Brown M*/
15879371c9d4SSatish Balay   #define PetscFunctionBeginHot \
15889371c9d4SSatish Balay     do { \
1589362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_TRUE); \
15902d53ad75SBarry Smith       PetscRegister__FUNCT__(); \
159153c77d0aSJed Brown     } while (0)
159253c77d0aSJed Brown 
1593a8d2bbe5SBarry Smith   /*MC
1594530d85c1SBarry Smith    PetscFunctionBeginUser - First executable line of user provided routines
1595a8d2bbe5SBarry Smith 
1596a8d2bbe5SBarry Smith    Synopsis:
1597aaa7dc30SBarry Smith    #include <petscsys.h>
1598a8d2bbe5SBarry Smith    void PetscFunctionBeginUser;
1599a8d2bbe5SBarry Smith 
160020f4b53cSBarry Smith    Not Collective; No Fortran Support
1601a8d2bbe5SBarry Smith 
1602a8d2bbe5SBarry Smith    Usage:
1603a8d2bbe5SBarry Smith .vb
1604a8d2bbe5SBarry Smith      int something;
1605a8d2bbe5SBarry Smith 
1606ac285190SBarry Smith      PetscFunctionBeginUser;
1607a8d2bbe5SBarry Smith .ve
1608a8d2bbe5SBarry Smith 
160920f4b53cSBarry Smith    Level: intermediate
161020f4b53cSBarry Smith 
1611a8d2bbe5SBarry Smith    Notes:
1612530d85c1SBarry Smith    Functions that incorporate this must call `PetscFunctionReturn()` instead of return except for main().
1613530d85c1SBarry Smith 
1614530d85c1SBarry Smith    May be used before `PetscInitialize()`
16151957e957SBarry Smith 
1616530d85c1SBarry Smith    This is identical to `PetscFunctionBegin` except it labels the routine as a user
1617ac285190SBarry Smith    routine instead of as a PETSc library routine.
1618ac285190SBarry Smith 
1619586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, `PetscFunctionBeginHot`, `PetscStackPushNoCheck()`
1620a8d2bbe5SBarry Smith M*/
16219371c9d4SSatish Balay   #define PetscFunctionBeginUser \
16229371c9d4SSatish Balay     do { \
1623362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 2, PETSC_FALSE); \
1624a8d2bbe5SBarry Smith       PetscRegister__FUNCT__(); \
1625a8d2bbe5SBarry Smith     } while (0)
1626a8d2bbe5SBarry Smith 
1627586f9135SBarry Smith   /*MC
1628586f9135SBarry Smith    PetscStackPush - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is
1629586f9135SBarry Smith    currently in the source code and verifies the memory is not corrupted.
1630586f9135SBarry Smith 
1631586f9135SBarry Smith    Synopsis:
1632586f9135SBarry Smith    #include <petscsys.h>
1633586f9135SBarry Smith    void PetscStackPush(char *funct)
1634586f9135SBarry Smith 
16357cdbe19fSJose E. Roman    Not Collective
16367cdbe19fSJose E. Roman 
1637586f9135SBarry Smith    Input Parameter:
1638586f9135SBarry Smith .  funct - the function name
1639586f9135SBarry Smith 
1640586f9135SBarry Smith    Level: developer
1641586f9135SBarry Smith 
1642586f9135SBarry Smith    Notes:
1643586f9135SBarry 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
1644586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1645586f9135SBarry Smith    help debug the problem.
1646586f9135SBarry Smith 
164795bd0b28SBarry Smith    The default stack is a global variable called `petscstack`.
1648586f9135SBarry Smith 
1649586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1650586f9135SBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1651586f9135SBarry Smith M*/
16529371c9d4SSatish Balay   #define PetscStackPush(n) \
16539371c9d4SSatish Balay     do { \
1654362febeeSStefano Zampini       PetscStackPushNoCheck(n, 0, PETSC_FALSE); \
165515681b3cSBarry Smith       CHKMEMQ; \
165615681b3cSBarry Smith     } while (0)
16573a40ed3dSBarry Smith 
1658586f9135SBarry Smith   /*MC
1659586f9135SBarry Smith    PetscStackPop - Pops a function name from the PETSc default stack that tracks where the running program is
1660586f9135SBarry Smith    currently in the source code and verifies the memory is not corrupted.
1661586f9135SBarry Smith 
1662586f9135SBarry Smith    Synopsis:
1663586f9135SBarry Smith    #include <petscsys.h>
1664586f9135SBarry Smith    void PetscStackPop
1665586f9135SBarry Smith 
16667cdbe19fSJose E. Roman    Not Collective
16677cdbe19fSJose E. Roman 
1668586f9135SBarry Smith    Level: developer
1669586f9135SBarry Smith 
1670586f9135SBarry Smith    Notes:
1671586f9135SBarry 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
1672586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1673586f9135SBarry Smith    help debug the problem.
1674586f9135SBarry Smith 
167595bd0b28SBarry Smith    The default stack is a global variable called `petscstack`.
1676586f9135SBarry Smith 
1677586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPopNoCheck()`, `PetscStackPush()`
1678586f9135SBarry Smith M*/
16799371c9d4SSatish Balay   #define PetscStackPop \
16809371c9d4SSatish Balay     do { \
1681441dd030SJed Brown       CHKMEMQ; \
1682362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
168315681b3cSBarry Smith     } while (0)
1684d64ed03dSBarry Smith 
168530de9b25SBarry Smith   /*MC
16860a57284eSJacob Faibussowitsch    PetscFunctionReturn - Last executable line of each PETSc function used for error
16870a57284eSJacob Faibussowitsch    handling. Replaces `return()`.
168830de9b25SBarry Smith 
168930de9b25SBarry Smith    Synopsis:
16900a57284eSJacob Faibussowitsch    #include <petscerror.h>
16910a57284eSJacob Faibussowitsch    void PetscFunctionReturn(...)
169230de9b25SBarry Smith 
1693667f096bSBarry Smith    Not Collective; No Fortran Support
1694eca87e8dSBarry Smith 
16955687f061SJacob Faibussowitsch    Level: beginner
16965687f061SJacob Faibussowitsch 
16970a57284eSJacob Faibussowitsch    Notes:
16980a57284eSJacob Faibussowitsch    This routine is a macro, so while it does not "return" anything itself, it does return from
16990a57284eSJacob Faibussowitsch    the function in the literal sense.
17000a57284eSJacob Faibussowitsch 
17010a57284eSJacob Faibussowitsch    Usually the return value is the integer literal `0` (for example in any function returning
17020a57284eSJacob Faibussowitsch    `PetscErrorCode`), however it is possible to return any arbitrary type. The arguments of
17030a57284eSJacob Faibussowitsch    this macro are placed before the `return` statement as-is.
17040a57284eSJacob Faibussowitsch 
17050a57284eSJacob Faibussowitsch    Any routine which returns via `PetscFunctionReturn()` must begin with a corresponding
17060a57284eSJacob Faibussowitsch    `PetscFunctionBegin`.
17070a57284eSJacob Faibussowitsch 
17080a57284eSJacob Faibussowitsch    For routines which return `void` use `PetscFunctionReturnVoid()` instead.
17090a57284eSJacob Faibussowitsch 
17100a57284eSJacob Faibussowitsch    Example Usage:
171130de9b25SBarry Smith .vb
17120a57284eSJacob Faibussowitsch    PetscErrorCode foo(int *x)
17130a57284eSJacob Faibussowitsch    {
17140a57284eSJacob Faibussowitsch      PetscFunctionBegin; // don't forget the begin!
17150a57284eSJacob Faibussowitsch      *x = 10;
17163ba16761SJacob Faibussowitsch      PetscFunctionReturn(PETSC_SUCCESS);
171730de9b25SBarry Smith    }
171830de9b25SBarry Smith .ve
171930de9b25SBarry Smith 
17200a57284eSJacob Faibussowitsch    May return any arbitrary type\:
17210a57284eSJacob Faibussowitsch .vb
17220a57284eSJacob Faibussowitsch   struct Foo
17230a57284eSJacob Faibussowitsch   {
17240a57284eSJacob Faibussowitsch     int x;
17250a57284eSJacob Faibussowitsch   };
17260a57284eSJacob Faibussowitsch 
17270a57284eSJacob Faibussowitsch   struct Foo make_foo(int value)
17280a57284eSJacob Faibussowitsch   {
17290a57284eSJacob Faibussowitsch     struct Foo f;
17300a57284eSJacob Faibussowitsch 
17310a57284eSJacob Faibussowitsch     PetscFunctionBegin;
17320a57284eSJacob Faibussowitsch     f.x = value;
17330a57284eSJacob Faibussowitsch     PetscFunctionReturn(f);
17340a57284eSJacob Faibussowitsch   }
17350a57284eSJacob Faibussowitsch .ve
17360a57284eSJacob Faibussowitsch 
17370a57284eSJacob Faibussowitsch .seealso: `PetscFunctionBegin`, `PetscFunctionBeginUser`, `PetscFunctionReturnVoid()`,
17380a57284eSJacob Faibussowitsch           `PetscStackPopNoCheck()`
173930de9b25SBarry Smith M*/
17405687f061SJacob Faibussowitsch   #define PetscFunctionReturn(...) \
17419371c9d4SSatish Balay     do { \
1742362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
17435687f061SJacob Faibussowitsch       return __VA_ARGS__; \
174427104ee2SJacob Faibussowitsch     } while (0)
1745d64ed03dSBarry Smith 
17460a57284eSJacob Faibussowitsch   /*MC
17470a57284eSJacob Faibussowitsch   PetscFunctionReturnVoid - Like `PetscFunctionReturn()` but returns `void`
17480a57284eSJacob Faibussowitsch 
17490a57284eSJacob Faibussowitsch   Synopsis:
17500a57284eSJacob Faibussowitsch   #include <petscerror.h>
17510a57284eSJacob Faibussowitsch   void PetscFunctionReturnVoid()
17520a57284eSJacob Faibussowitsch 
17530a57284eSJacob Faibussowitsch   Not Collective
17540a57284eSJacob Faibussowitsch 
17555687f061SJacob Faibussowitsch   Level: beginner
17565687f061SJacob Faibussowitsch 
17575687f061SJacob Faibussowitsch   Note:
17580a57284eSJacob Faibussowitsch   Behaves identically to `PetscFunctionReturn()` except that it returns `void`. That is, this
17595687f061SJacob Faibussowitsch   macro culminates with `return`.
17600a57284eSJacob Faibussowitsch 
17610a57284eSJacob Faibussowitsch   Example Usage:
17620a57284eSJacob Faibussowitsch .vb
17630a57284eSJacob Faibussowitsch   void foo()
17640a57284eSJacob Faibussowitsch   {
17650a57284eSJacob Faibussowitsch     PetscFunctionBegin; // must start with PetscFunctionBegin!
17660a57284eSJacob Faibussowitsch     bar();
17670a57284eSJacob Faibussowitsch     baz();
17680a57284eSJacob Faibussowitsch     PetscFunctionReturnVoid();
17690a57284eSJacob Faibussowitsch   }
17700a57284eSJacob Faibussowitsch .ve
17710a57284eSJacob Faibussowitsch 
17720a57284eSJacob Faibussowitsch .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, PetscFunctionBeginUser`
17730a57284eSJacob Faibussowitsch M*/
17749371c9d4SSatish Balay   #define PetscFunctionReturnVoid() \
17759371c9d4SSatish Balay     do { \
1776362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
177727104ee2SJacob Faibussowitsch       return; \
177827104ee2SJacob Faibussowitsch     } while (0)
177927104ee2SJacob Faibussowitsch #else /* PETSC_USE_DEBUG */
178027104ee2SJacob Faibussowitsch   #define PetscStackPushNoCheck(funct, petsc_routine, hot)
178137154d10SBarry Smith   #define PetscStackUpdateLine
1782792fecdfSBarry Smith   #define PetscStackPushExternal(funct)
17833ba16761SJacob Faibussowitsch   #define PetscStackPopNoCheck(...)
178427104ee2SJacob Faibussowitsch   #define PetscStackClearTop
17853a40ed3dSBarry Smith   #define PetscFunctionBegin
17860bdf7c52SPeter Brune   #define PetscFunctionBeginUser
1787a2f94806SJed Brown   #define PetscFunctionBeginHot
17880a57284eSJacob Faibussowitsch   #define PetscFunctionReturn(...)  return __VA_ARGS__
17895665465eSBarry Smith   #define PetscFunctionReturnVoid() return
1790812af9f3SBarry Smith   #define PetscStackPop             CHKMEMQ
1791812af9f3SBarry Smith   #define PetscStackPush(f)         CHKMEMQ
179227104ee2SJacob Faibussowitsch #endif /* PETSC_USE_DEBUG */
17933a40ed3dSBarry Smith 
17946d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
17953ba16761SJacob Faibussowitsch   #define PetscStackCallExternalVoid(...)
17963ba16761SJacob Faibussowitsch template <typename F, typename... Args>
17973ba16761SJacob Faibussowitsch void PetscCallExternal(F, Args...);
17986d210af2SJacob Faibussowitsch #else
1799586f9135SBarry Smith   /*MC
1800e77caa6dSBarry Smith     PetscStackCallExternalVoid - Calls an external library routine or user function after pushing the name of the routine on the stack.
1801eb6b5d47SBarry Smith 
1802eb6b5d47SBarry Smith    Input Parameters:
1803eb6b5d47SBarry Smith +   name    - string that gives the name of the function being called
1804586f9135SBarry Smith -   routine - actual call to the routine, for example, functionname(a,b)
1805fd3f9acdSBarry Smith 
1806586f9135SBarry Smith    Level: developer
1807eb6b5d47SBarry Smith 
180895bd0b28SBarry Smith    Notes:
1809792fecdfSBarry Smith    Often one should use `PetscCallExternal()` instead. This routine is intended for external library routines that DO NOT return error codes
1810eb6b5d47SBarry Smith 
1811586f9135SBarry Smith    In debug mode this also checks the memory for corruption at the end of the function call.
1812586f9135SBarry Smith 
181395bd0b28SBarry Smith    Certain external packages, such as BLAS/LAPACK may have their own macros, `PetscCallBLAS()` for managing the call, error checking, etc.
1814586f9135SBarry Smith 
1815586f9135SBarry Smith    Developer Note:
1816586f9135SBarry 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.
1817586f9135SBarry Smith 
1818792fecdfSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscCallExternal()`, `PetscCallBLAS()`
1819586f9135SBarry Smith @*/
18203ba16761SJacob Faibussowitsch   #define PetscStackCallExternalVoid(name, ...) \
18219371c9d4SSatish Balay     do { \
182228770333SStefano Zampini       PetscStackPushExternal(name); \
18233ba16761SJacob Faibussowitsch       __VA_ARGS__; \
18249371c9d4SSatish Balay       PetscStackPop; \
18259371c9d4SSatish Balay     } while (0)
1826eb6b5d47SBarry Smith 
1827586f9135SBarry Smith   /*MC
1828792fecdfSBarry Smith     PetscCallExternal - Calls an external library routine that returns an error code after pushing the name of the routine on the stack.
1829fd3f9acdSBarry Smith 
1830fd3f9acdSBarry Smith    Input Parameters:
1831fd3f9acdSBarry Smith +  func - name of the routine
1832586f9135SBarry Smith -  args - arguments to the routine
1833586f9135SBarry Smith 
1834586f9135SBarry Smith    Level: developer
1835fd3f9acdSBarry Smith 
183695452b02SPatrick Sanan    Notes:
1837e77caa6dSBarry Smith    This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not.
1838dbf62e16SBarry Smith 
1839586f9135SBarry Smith    In debug mode this also checks the memory for corruption at the end of the function call.
1840fd3f9acdSBarry Smith 
184187497f52SBarry Smith    Assumes the error return code of the function is an integer and that a value of 0 indicates success
184287497f52SBarry Smith 
1843586f9135SBarry Smith    Developer Note:
1844d5b43468SJose E. Roman    This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc.
1845586f9135SBarry Smith 
1846e77caa6dSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()`
1847586f9135SBarry Smith M*/
18489371c9d4SSatish Balay   #define PetscCallExternal(func, ...) \
18499371c9d4SSatish Balay     do { \
1850a74df02fSJacob Faibussowitsch       PetscStackPush(PetscStringize(func)); \
18513ba16761SJacob Faibussowitsch       int ierr_petsc_call_external_ = func(__VA_ARGS__); \
18521d4906efSStefano Zampini       PetscStackPop; \
18533ba16761SJacob 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_); \
1854fd3f9acdSBarry Smith     } while (0)
18556d210af2SJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */
1856