xref: /petsc/include/petscerror.h (revision ede9db9363e1fdaaa09befd664c8164883ccce80)
154a8ef01SBarry Smith /*
2f621e05eSBarry Smith     Contains all error handling interfaces for PETSc.
354a8ef01SBarry Smith */
4a4963045SJacob Faibussowitsch #pragma once
57233ce55SJed Brown // IWYU pragma: private, include "petscsys.h"
66c7e564aSBarry Smith 
7e1bf4ed2SJacob Faibussowitsch #include <petscmacros.h>
8e1bf4ed2SJacob Faibussowitsch #include <petscsystypes.h>
9e1bf4ed2SJacob Faibussowitsch 
10af75f258SJacob Faibussowitsch #if defined(__cplusplus)
11af75f258SJacob Faibussowitsch   #include <exception> // std::exception
12af75f258SJacob Faibussowitsch #endif
13af75f258SJacob Faibussowitsch 
14ac09b921SBarry Smith /* SUBMANSEC = Sys */
15ac09b921SBarry Smith 
16edd03b47SJacob Faibussowitsch #define SETERRQ1(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
17edd03b47SJacob Faibussowitsch #define SETERRQ2(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
18edd03b47SJacob Faibussowitsch #define SETERRQ3(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
19edd03b47SJacob Faibussowitsch #define SETERRQ4(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
20edd03b47SJacob Faibussowitsch #define SETERRQ5(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
21edd03b47SJacob Faibussowitsch #define SETERRQ6(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
22edd03b47SJacob Faibussowitsch #define SETERRQ7(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
23edd03b47SJacob Faibussowitsch #define SETERRQ8(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
24edd03b47SJacob Faibussowitsch #define SETERRQ9(...) PETSC_DEPRECATED_MACRO(3, 17, 0, "SETERRQ", ) SETERRQ(__VA_ARGS__)
2598921bdaSJacob Faibussowitsch 
2630de9b25SBarry Smith /*MC
271957e957SBarry Smith    SETERRQ - Macro to be called when an error has been detected,
2830de9b25SBarry Smith 
2930de9b25SBarry Smith    Synopsis:
30aaa7dc30SBarry Smith    #include <petscsys.h>
3198921bdaSJacob Faibussowitsch    PetscErrorCode SETERRQ(MPI_Comm comm, PetscErrorCode ierr, char *message, ...)
3230de9b25SBarry Smith 
33d083f849SBarry Smith    Collective
3430de9b25SBarry Smith 
3530de9b25SBarry Smith    Input Parameters:
3695bd0b28SBarry Smith +  comm    - An MPI communicator, use `PETSC_COMM_SELF` unless you know all ranks of another communicator will detect the error
373af045c5SBarry Smith .  ierr    - nonzero error code, see the list of standard error codes in include/petscerror.h
3830de9b25SBarry Smith -  message - error message
3930de9b25SBarry Smith 
4030de9b25SBarry Smith   Level: beginner
4130de9b25SBarry Smith 
4230de9b25SBarry Smith    Notes:
4387497f52SBarry Smith    This is rarely needed, one should use `PetscCheck()` and `PetscCall()` and friends to automatically handle error conditions.
4430de9b25SBarry Smith    Once the error handler is called the calling function is then returned from with the given error code.
4530de9b25SBarry Smith 
4687497f52SBarry Smith    Experienced users can set the error handler with `PetscPushErrorHandler()`.
4730de9b25SBarry Smith 
4816a05f60SBarry Smith    Fortran Note:
49989c49a2SBarry Smith    `SETERRQ()` may be called from Fortran subroutines but `SETERRA()` must be called from the
503b1008b8SBarry Smith    Fortran main program.
513b1008b8SBarry Smith 
52db781477SPatrick Sanan .seealso: `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
53af27ebaaSBarry Smith           `PetscError()`, `PetscCall()`, `CHKMEMQ`, `CHKERRA()`, `PetscCallMPI()`, `PetscErrorCode`
5430de9b25SBarry Smith M*/
55e809461dSJacob Faibussowitsch #define SETERRQ(comm, ierr, ...) \
56e809461dSJacob Faibussowitsch   do { \
57e809461dSJacob Faibussowitsch     PetscErrorCode ierr_seterrq_petsc_ = PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \
58e809461dSJacob Faibussowitsch     return ierr_seterrq_petsc_ ? ierr_seterrq_petsc_ : PETSC_ERR_RETURN; \
59e809461dSJacob Faibussowitsch   } while (0)
60986eef2eSBarry Smith 
61648c30bcSBarry Smith #define SETERRQNULL(comm, ierr, ...) \
62648c30bcSBarry Smith   do { \
63648c30bcSBarry Smith     (void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \
64648c30bcSBarry Smith     return NULL; \
65648c30bcSBarry Smith   } while (0)
66648c30bcSBarry Smith 
67ffc4695bSBarry Smith /*
68ffc4695bSBarry Smith     Returned from PETSc functions that are called from MPI, such as related to attributes
69ffc4695bSBarry Smith       Do not confuse PETSC_MPI_ERROR_CODE and PETSC_ERR_MPI, the first is registered with MPI and returned to MPI as
70888a1fb3SBarry 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.
71ffc4695bSBarry Smith */
72ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CLASS;
73ffc4695bSBarry Smith PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CODE;
74ffc4695bSBarry Smith 
75986eef2eSBarry Smith /*MC
76986eef2eSBarry Smith    SETERRMPI - Macro to be called when an error has been detected within an MPI callback function
77986eef2eSBarry Smith 
78989c49a2SBarry Smith    No Fortran Support
79989c49a2SBarry Smith 
80986eef2eSBarry Smith    Synopsis:
81986eef2eSBarry Smith    #include <petscsys.h>
8298921bdaSJacob Faibussowitsch    PetscErrorCode SETERRMPI(MPI_Comm comm, PetscErrorCode ierr, char *message, ...)
83986eef2eSBarry Smith 
84d083f849SBarry Smith    Collective
85986eef2eSBarry Smith 
86986eef2eSBarry Smith    Input Parameters:
8795bd0b28SBarry Smith +  comm    - An MPI communicator, use `PETSC_COMM_SELF` unless you know all ranks of another communicator will detect the error
88986eef2eSBarry Smith .  ierr    - nonzero error code, see the list of standard error codes in include/petscerror.h
89986eef2eSBarry Smith -  message - error message
90986eef2eSBarry Smith 
91986eef2eSBarry Smith   Level: developer
92986eef2eSBarry Smith 
9395bd0b28SBarry Smith    Note:
9487497f52SBarry 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`
9587497f52SBarry Smith   which is registered with `MPI_Add_error_code()` when PETSc is initialized.
96986eef2eSBarry Smith 
97af27ebaaSBarry Smith .seealso: `SETERRQ()`, `PetscCall()`, `PetscCallMPI()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `PetscErrorCode`
98986eef2eSBarry Smith M*/
993ba16761SJacob Faibussowitsch #define SETERRMPI(comm, ierr, ...) return ((void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__), PETSC_MPI_ERROR_CODE)
100ee8199e6SMatthew G. Knepley 
101ee8199e6SMatthew G. Knepley /*MC
102f388eb8bSPatrick Sanan    SETERRA - Fortran-only macro that can be called when an error has been detected from the main program
103f388eb8bSPatrick Sanan 
104f388eb8bSPatrick Sanan    Synopsis:
105f388eb8bSPatrick Sanan    #include <petscsys.h>
106f388eb8bSPatrick Sanan    PetscErrorCode SETERRA(MPI_Comm comm, PetscErrorCode ierr, char *message)
107f388eb8bSPatrick Sanan 
108f388eb8bSPatrick Sanan    Collective
109f388eb8bSPatrick Sanan 
110f388eb8bSPatrick Sanan    Input Parameters:
11195bd0b28SBarry Smith +  comm    - An MPI communicator, so that the error can be collective
112f388eb8bSPatrick Sanan .  ierr    - nonzero error code, see the list of standard error codes in include/petscerror.h
113e55b9c91SBarry Smith -  message - error message in the `printf()` format
114f388eb8bSPatrick Sanan 
115f388eb8bSPatrick Sanan    Level: beginner
116f388eb8bSPatrick Sanan 
117f388eb8bSPatrick Sanan    Notes:
11887497f52SBarry Smith    This should only be used with Fortran. With C/C++, use `SETERRQ()`.
119f388eb8bSPatrick Sanan 
12087497f52SBarry Smith    `SETERRQ()` may be called from Fortran subroutines but `SETERRA()` must be called from the
121f388eb8bSPatrick Sanan     Fortran main program.
122f388eb8bSPatrick Sanan 
123af27ebaaSBarry Smith .seealso: `SETERRQ()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()`, `PetscErrorCode`
124f388eb8bSPatrick Sanan M*/
125f388eb8bSPatrick Sanan 
126f388eb8bSPatrick Sanan /*MC
127fa190f98SMatthew G. Knepley    SETERRABORT - Macro that can be called when an error has been detected,
128fa190f98SMatthew G. Knepley 
129fa190f98SMatthew G. Knepley    Synopsis:
130fa190f98SMatthew G. Knepley    #include <petscsys.h>
13198921bdaSJacob Faibussowitsch    PetscErrorCode SETERRABORT(MPI_Comm comm, PetscErrorCode ierr, char *message, ...)
132fa190f98SMatthew G. Knepley 
133d083f849SBarry Smith    Collective
134fa190f98SMatthew G. Knepley 
135fa190f98SMatthew G. Knepley    Input Parameters:
13695bd0b28SBarry Smith +  comm    - An MPI communicator, so that the error can be collective
1373af045c5SBarry Smith .  ierr    - nonzero error code, see the list of standard error codes in include/petscerror.h
138e55b9c91SBarry Smith -  message - error message in the `printf()` format
139fa190f98SMatthew G. Knepley 
140fa190f98SMatthew G. Knepley    Level: beginner
141fa190f98SMatthew G. Knepley 
142fa190f98SMatthew G. Knepley    Notes:
14387497f52SBarry Smith    This function just calls `MPI_Abort()`.
14487497f52SBarry Smith 
14587497f52SBarry Smith    This should only be called in routines that cannot return an error code, such as in C++ constructors.
146fa190f98SMatthew G. Knepley 
147989c49a2SBarry Smith    Fortran Note:
148989c49a2SBarry Smith    Use `SETERRA()` in Fortran main program and `SETERRQ()` in Fortran subroutines
149989c49a2SBarry Smith 
150989c49a2SBarry Smith    Developer Note:
151989c49a2SBarry Smith    In Fortran `SETERRA()` could be called `SETERRABORT()` since they serve the same purpose
152989c49a2SBarry Smith 
153af27ebaaSBarry Smith .seealso: `SETERRQ()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `PetscCall()`, `CHKMEMQ`, `PetscErrorCode`
154fa190f98SMatthew G. Knepley M*/
1559371c9d4SSatish Balay #define SETERRABORT(comm, ierr, ...) \
1569371c9d4SSatish Balay   do { \
1573ba16761SJacob Faibussowitsch     (void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr, PETSC_ERROR_INITIAL, __VA_ARGS__); \
158dd460d27SBarry Smith     (void)MPI_Abort(comm, ierr); \
15998921bdaSJacob Faibussowitsch   } while (0)
1609a00fa46SSatish Balay 
16130de9b25SBarry Smith /*MC
162139c8099SSatish Balay   PetscCheck - Checks that a particular condition is true; if not true, then returns the provided error code
1632c71b3e2SJacob Faibussowitsch 
1642c71b3e2SJacob Faibussowitsch   Synopsis:
1657233ce55SJed Brown   #include <petscsys.h>
1662c71b3e2SJacob Faibussowitsch   void PetscCheck(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
1672c71b3e2SJacob Faibussowitsch 
1687cdbe19fSJose E. Roman   Collective; No Fortran Support
1692c71b3e2SJacob Faibussowitsch 
1702c71b3e2SJacob Faibussowitsch   Input Parameters:
1712c71b3e2SJacob Faibussowitsch + cond    - The boolean condition
1722c71b3e2SJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
1732c71b3e2SJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
174e55b9c91SBarry Smith - message - Error message in the `printf()` format
1752c71b3e2SJacob Faibussowitsch 
176667f096bSBarry Smith   Level: beginner
177667f096bSBarry Smith 
1782c71b3e2SJacob Faibussowitsch   Notes:
1792c71b3e2SJacob Faibussowitsch   Enabled in both optimized and debug builds.
1802c71b3e2SJacob Faibussowitsch 
181a2454668SBarry Smith   As a general rule, `PetscCheck()` is used to check "usage error" (for example, passing an incorrect value as a function argument),
182a2454668SBarry Smith   `PetscAssert()` is used to "check for bugs in PETSc" (for example, is a value in a PETSc data structure nonsensical).
183a2454668SBarry Smith   However, for functions that are called in a "hot spot", for example, thousands of times in a loop, `PetscAssert()` should be used instead
184a2454668SBarry Smith   of `PetscCheck()` since the former is compiled out in PETSc's optimization code.
185a2454668SBarry Smith 
18687497f52SBarry Smith   Calls `SETERRQ()` if the assertion fails, so can only be called from functions returning a
18787497f52SBarry Smith   `PetscErrorCode` (or equivalent type after conversion).
1882c71b3e2SJacob Faibussowitsch 
1897c5b2466SBarry Smith .seealso: `PetscAssert()`, `PetscCheckReturnMPI()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheckAbort()`, `PetscErrorCode`
1909566063dSJacob Faibussowitsch M*/
1919371c9d4SSatish Balay #define PetscCheck(cond, comm, ierr, ...) \
1929371c9d4SSatish Balay   do { \
1939371c9d4SSatish Balay     if (PetscUnlikely(!(cond))) SETERRQ(comm, ierr, __VA_ARGS__); \
1949371c9d4SSatish Balay   } while (0)
1952c71b3e2SJacob Faibussowitsch 
1962c71b3e2SJacob Faibussowitsch /*MC
1977c5b2466SBarry Smith   PetscCheckReturnMPI - Checks that a particular condition is true; if not true, then returns an MPI error code.
1987c5b2466SBarry Smith   To check for errors in PETSc-provided MPI callbacks.
1997c5b2466SBarry Smith 
2007c5b2466SBarry Smith   Synopsis:
2017233ce55SJed Brown   #include <petscsys.h>
2027c5b2466SBarry Smith   void PetscCheckReturnMPI(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
2037c5b2466SBarry Smith 
2047c5b2466SBarry Smith   Collective; No Fortran Support
2057c5b2466SBarry Smith 
2067c5b2466SBarry Smith   Input Parameters:
2077c5b2466SBarry Smith + cond    - The boolean condition
2087c5b2466SBarry Smith . comm    - The communicator on which the check can be collective on
2097c5b2466SBarry Smith . ierr    - A nonzero error code, see include/petscerror.h for the complete list
2107c5b2466SBarry Smith - message - Error message in the `printf()` format
2117c5b2466SBarry Smith 
2127c5b2466SBarry Smith   Level: beginner
2137c5b2466SBarry Smith 
2147c5b2466SBarry Smith   Note:
2157c5b2466SBarry Smith   Enabled in both optimized and debug builds.
2167c5b2466SBarry Smith 
2177c5b2466SBarry Smith .seealso: `PetscCheck()`, `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheckAbort()`, `PetscErrorCode`
2187c5b2466SBarry Smith M*/
2197c5b2466SBarry Smith #define PetscCheckReturnMPI(cond, comm, ierr, ...) \
2207c5b2466SBarry Smith   do { \
2217c5b2466SBarry Smith     if (PetscUnlikely(!(cond))) SETERRMPI(comm, ierr, __VA_ARGS__); \
2227c5b2466SBarry Smith   } while (0)
2237c5b2466SBarry Smith 
2247c5b2466SBarry Smith /*MC
2254be741a6SBarry Smith   PetscCheckAbort - Check that a particular condition is true, otherwise prints error and aborts
2264be741a6SBarry Smith 
2274be741a6SBarry Smith   Synopsis:
2287233ce55SJed Brown   #include <petscsys.h>
2294be741a6SBarry Smith   void PetscCheckAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
2304be741a6SBarry Smith 
2317cdbe19fSJose E. Roman   Collective; No Fortran Support
2324be741a6SBarry Smith 
2334be741a6SBarry Smith   Input Parameters:
2344be741a6SBarry Smith + cond    - The boolean condition
2354be741a6SBarry Smith . comm    - The communicator on which the check can be collective on
2364be741a6SBarry Smith . ierr    - A nonzero error code, see include/petscerror.h for the complete list
237e55b9c91SBarry Smith - message - Error message in the `printf()` format
2384be741a6SBarry Smith 
239667f096bSBarry Smith   Level: developer
240667f096bSBarry Smith 
2414be741a6SBarry Smith   Notes:
2424be741a6SBarry Smith   Enabled in both optimized and debug builds.
2434be741a6SBarry Smith 
24487497f52SBarry Smith   Calls `SETERRABORT()` if the assertion fails, can be called from a function that does not return an
24587497f52SBarry Smith   error code, such as a C++ constructor. usually `PetscCheck()` should be used.
2464be741a6SBarry Smith 
247af27ebaaSBarry Smith .seealso: `PetscAssertAbort()`, `PetscAssert()`, `SETERRQ()`, `PetscError()`, `PetscCall()`, `PetscCheck()`, `SETERRABORT()`, `PetscErrorCode`
2484be741a6SBarry Smith M*/
2499371c9d4SSatish Balay #define PetscCheckAbort(cond, comm, ierr, ...) \
2500e6b6b59SJacob Faibussowitsch   do { \
2510e6b6b59SJacob Faibussowitsch     if (PetscUnlikely(!(cond))) SETERRABORT(comm, ierr, __VA_ARGS__); \
252c7481402SJacob Faibussowitsch   } while (0)
2534be741a6SBarry Smith 
2544be741a6SBarry Smith /*MC
2559ace16cdSJacob Faibussowitsch   PetscAssert - Assert that a particular condition is true
2569ace16cdSJacob Faibussowitsch 
2579ace16cdSJacob Faibussowitsch   Synopsis:
2587233ce55SJed Brown   #include <petscsys.h>
2599ace16cdSJacob Faibussowitsch   void PetscAssert(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
2609ace16cdSJacob Faibussowitsch 
2617cdbe19fSJose E. Roman   Collective; No Fortran Support
2629ace16cdSJacob Faibussowitsch 
2639ace16cdSJacob Faibussowitsch   Input Parameters:
2649ace16cdSJacob Faibussowitsch + cond    - The boolean condition
2659ace16cdSJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
2669ace16cdSJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
267af27ebaaSBarry Smith - message - Error message in `printf()` format
2689ace16cdSJacob Faibussowitsch 
269667f096bSBarry Smith   Level: beginner
270667f096bSBarry Smith 
2719ace16cdSJacob Faibussowitsch   Notes:
272c7481402SJacob Faibussowitsch   Equivalent to `PetscCheck()` if debugging is enabled, and `PetscAssume(cond)` otherwise.
2732c71b3e2SJacob Faibussowitsch 
27487497f52SBarry Smith   See `PetscCheck()` for usage and behaviour.
27587497f52SBarry Smith 
27687497f52SBarry Smith   This is needed instead of simply using `assert()` because this correctly handles the collective nature of errors under MPI
2779ace16cdSJacob Faibussowitsch 
278af27ebaaSBarry Smith .seealso: `PetscCheck()`, `SETERRQ()`, `PetscError()`, `PetscAssertAbort()`, `PetscErrorCode`
2799566063dSJacob Faibussowitsch M*/
280c7481402SJacob Faibussowitsch #if PetscDefined(USE_DEBUG)
281c7481402SJacob Faibussowitsch   #define PetscAssert(cond, comm, ierr, ...) PetscCheck(cond, comm, ierr, __VA_ARGS__)
282c7481402SJacob Faibussowitsch #else
283c7481402SJacob Faibussowitsch   #define PetscAssert(cond, ...) PetscAssume(cond)
284c7481402SJacob Faibussowitsch #endif
2859ace16cdSJacob Faibussowitsch 
2869ace16cdSJacob Faibussowitsch /*MC
2870e6b6b59SJacob Faibussowitsch   PetscAssertAbort - Assert that a particular condition is true, otherwise prints error and aborts
2880e6b6b59SJacob Faibussowitsch 
2890e6b6b59SJacob Faibussowitsch   Synopsis:
2907233ce55SJed Brown   #include <petscsys.h>
2910e6b6b59SJacob Faibussowitsch   void PetscAssertAbort(bool cond, MPI_Comm comm, PetscErrorCode ierr, const char *message, ...)
2920e6b6b59SJacob Faibussowitsch 
2937cdbe19fSJose E. Roman   Collective; No Fortran Support
2940e6b6b59SJacob Faibussowitsch 
2950e6b6b59SJacob Faibussowitsch   Input Parameters:
2960e6b6b59SJacob Faibussowitsch + cond    - The boolean condition
2970e6b6b59SJacob Faibussowitsch . comm    - The communicator on which the check can be collective on
2980e6b6b59SJacob Faibussowitsch . ierr    - A nonzero error code, see include/petscerror.h for the complete list
299e55b9c91SBarry Smith - message - Error message in the `printf()` format
3000e6b6b59SJacob Faibussowitsch 
301667f096bSBarry Smith   Level: beginner
302667f096bSBarry Smith 
30395bd0b28SBarry Smith   Note:
3040e6b6b59SJacob Faibussowitsch   Enabled only in debug builds. See `PetscCheckAbort()` for usage.
3050e6b6b59SJacob Faibussowitsch 
3060e6b6b59SJacob Faibussowitsch .seealso: `PetscCheckAbort()`, `PetscAssert()`, `PetscCheck()`, `SETERRABORT()`, `PetscError()`
3070e6b6b59SJacob Faibussowitsch M*/
3083ba16761SJacob Faibussowitsch #if PetscDefined(USE_DEBUG)
3093ba16761SJacob Faibussowitsch   #define PetscAssertAbort(cond, comm, ierr, ...) PetscCheckAbort(cond, comm, ierr, __VA_ARGS__)
3103ba16761SJacob Faibussowitsch #else
3113ba16761SJacob Faibussowitsch   #define PetscAssertAbort(cond, comm, ierr, ...) PetscAssume(cond)
3123ba16761SJacob Faibussowitsch #endif
3130e6b6b59SJacob Faibussowitsch 
3140e6b6b59SJacob Faibussowitsch /*MC
3153ba16761SJacob Faibussowitsch   PetscCall - Calls a PETSc function and then checks the resulting error code, if it is
3163ba16761SJacob Faibussowitsch   non-zero it calls the error handler and returns from the current function with the error
3173ba16761SJacob Faibussowitsch   code.
3189566063dSJacob Faibussowitsch 
3199566063dSJacob Faibussowitsch   Synopsis:
3207233ce55SJed Brown   #include <petscsys.h>
32149c86fc7SBarry Smith   void PetscCall(PetscFunction(args))
3229566063dSJacob Faibussowitsch 
3239566063dSJacob Faibussowitsch   Not Collective
3249566063dSJacob Faibussowitsch 
3259566063dSJacob Faibussowitsch   Input Parameter:
32649c86fc7SBarry Smith . PetscFunction - any PETSc function that returns an error code
3279566063dSJacob Faibussowitsch 
328667f096bSBarry Smith   Level: beginner
329667f096bSBarry Smith 
3309566063dSJacob Faibussowitsch   Notes:
3319566063dSJacob Faibussowitsch   Once the error handler is called the calling function is then returned from with the given
33287497f52SBarry Smith   error code. Experienced users can set the error handler with `PetscPushErrorHandler()`.
3339566063dSJacob Faibussowitsch 
33487497f52SBarry Smith   `PetscCall()` cannot be used in functions returning a datatype not convertible to
3358af9f190SBarry Smith   `PetscErrorCode`. For example, `PetscCall()` may not be used in functions returning `void`, use
3363ba16761SJacob Faibussowitsch   `PetscCallAbort()` or `PetscCallVoid()` in this case.
3379566063dSJacob Faibussowitsch 
3389566063dSJacob Faibussowitsch   Example Usage:
3399566063dSJacob Faibussowitsch .vb
3409566063dSJacob Faibussowitsch   PetscCall(PetscInitiailize(...)); // OK to call even when PETSc is not yet initialized!
3419566063dSJacob Faibussowitsch 
3429566063dSJacob Faibussowitsch   struct my_struct
3439566063dSJacob Faibussowitsch   {
3449566063dSJacob Faibussowitsch     void *data;
3459566063dSJacob Faibussowitsch   } my_complex_type;
3469566063dSJacob Faibussowitsch 
3479566063dSJacob Faibussowitsch   struct my_struct bar(void)
3489566063dSJacob Faibussowitsch   {
3496aad120cSJose E. Roman     PetscCall(foo(15)); // ERROR PetscErrorCode not convertible to struct my_struct!
3509566063dSJacob Faibussowitsch   }
3519566063dSJacob Faibussowitsch 
3526aad120cSJose E. Roman   PetscCall(bar()) // ERROR input not convertible to PetscErrorCode
3539566063dSJacob Faibussowitsch .ve
3549566063dSJacob Faibussowitsch 
35587497f52SBarry Smith   It is also possible to call this directly on a `PetscErrorCode` variable
35649c86fc7SBarry Smith .vb
35749c86fc7SBarry Smith   PetscCall(ierr);  // check if ierr is nonzero
35849c86fc7SBarry Smith .ve
35949c86fc7SBarry Smith 
360792fecdfSBarry Smith   Should not be used to call callback functions provided by users, `PetscCallBack()` should be used in that situation.
361ef1023bdSBarry Smith 
3626a8be23eSBarry Smith   `PetscUseTypeMethod()` or `PetscTryTypeMethod()` should be used when calling functions pointers contained in a PETSc object's `ops` array
3636a8be23eSBarry Smith 
36449c86fc7SBarry Smith   Fortran Notes:
36598d50953SPierre Jolivet   The Fortran function in which this is used must declare a `PetscErrorCode` variable necessarily named `ierr`, and `ierr` must be
36687497f52SBarry Smith   the final argument to the PETSc function being called.
36749c86fc7SBarry Smith 
36898d50953SPierre Jolivet   In the main program and in Fortran subroutines that do not have `ierr` as the final return parameter, one
36987497f52SBarry Smith   should use `PetscCallA()`
37049c86fc7SBarry Smith 
37149c86fc7SBarry Smith   Example Fortran Usage:
37249c86fc7SBarry Smith .vb
37349c86fc7SBarry Smith   PetscErrorCode ierr
37449c86fc7SBarry Smith   Vec v
37549c86fc7SBarry Smith 
37649c86fc7SBarry Smith   ...
37749c86fc7SBarry Smith   PetscCall(VecShift(v, 1.0, ierr))
37849c86fc7SBarry Smith   PetscCallA(VecShift(v, 1.0, ierr))
37949c86fc7SBarry Smith .ve
38049c86fc7SBarry Smith 
381667f096bSBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`,
382667f096bSBarry Smith           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`,
383648c30bcSBarry Smith           `CHKERRMPI()`, `PetscCallBack()`, `PetscCallAbort()`, `PetscCallVoid()`, `PetscCallNull()`
384648c30bcSBarry Smith M*/
385648c30bcSBarry Smith 
386648c30bcSBarry Smith /*MC
387648c30bcSBarry Smith   PetscCallNull - Calls a PETSc function and then checks the resulting error code, if it is
388648c30bcSBarry Smith   non-zero it calls the error handler and returns a `NULL`
389648c30bcSBarry Smith 
390648c30bcSBarry Smith   Synopsis:
3917233ce55SJed Brown   #include <petscsys.h>
392648c30bcSBarry Smith   void PetscCallNull(PetscFunction(args))
393648c30bcSBarry Smith 
394648c30bcSBarry Smith   Not Collective; No Fortran Support
395648c30bcSBarry Smith 
396648c30bcSBarry Smith   Input Parameter:
397648c30bcSBarry Smith . PetscFunction - any PETSc function that returns something that can be returned as a `NULL`
398648c30bcSBarry Smith 
399648c30bcSBarry Smith   Level: developer
400648c30bcSBarry Smith 
401648c30bcSBarry Smith .seealso: `PetscCall()`, `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`,
402648c30bcSBarry Smith           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`,
403648c30bcSBarry Smith           `CHKERRMPI()`, `PetscCallBack()`, `PetscCallAbort()`, `PetscCallVoid()`, `PetscCall()`
4049566063dSJacob Faibussowitsch M*/
405ef1023bdSBarry Smith 
406ef1023bdSBarry Smith /*MC
40798d50953SPierre 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
40898d50953SPierre Jolivet    `PetscCall()` which should be used in other Fortran subroutines
409989c49a2SBarry Smith 
410989c49a2SBarry Smith    Synopsis:
411989c49a2SBarry Smith    #include <petscsys.h>
412989c49a2SBarry Smith    PetscErrorCode PetscCallA(PetscFunction(arguments, ierr))
413989c49a2SBarry Smith 
414989c49a2SBarry Smith    Collective
415989c49a2SBarry Smith 
416989c49a2SBarry Smith    Input Parameter:
417989c49a2SBarry Smith .  PetscFunction(arguments,ierr) - the call to the function
418989c49a2SBarry Smith 
419989c49a2SBarry Smith   Level: beginner
420989c49a2SBarry Smith 
421989c49a2SBarry Smith    Notes:
422989c49a2SBarry Smith    This should only be used with Fortran. With C/C++, use `PetscCall()` always.
423989c49a2SBarry Smith 
42498d50953SPierre Jolivet    The Fortran function in which this is used must declare a `PetscErrorCode` variable necessarily named `ierr`
425989c49a2SBarry Smith    Use `SETERRA()` to set an error in a Fortran main program and `SETERRQ()` in Fortran subroutines
426989c49a2SBarry Smith 
427989c49a2SBarry Smith .seealso: `SETERRQ()`, `SETERRA()`, `SETERRABORT()`, `PetscCall()`, `CHKERRA()`, `PetscCallAbort()`
428989c49a2SBarry Smith M*/
429989c49a2SBarry Smith 
430989c49a2SBarry Smith /*MC
431792fecdfSBarry 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
432ef1023bdSBarry Smith   handler and returns from the current function with the error code.
433ef1023bdSBarry Smith 
434ef1023bdSBarry Smith   Synopsis:
4357233ce55SJed Brown   #include <petscsys.h>
436792fecdfSBarry Smith   void PetscCallBack(const char *functionname, PetscFunction(args))
437ef1023bdSBarry Smith 
4387cdbe19fSJose E. Roman   Not Collective; No Fortran Support
439ef1023bdSBarry Smith 
440ef1023bdSBarry Smith   Input Parameters:
441ef1023bdSBarry Smith + functionname - the name of the function being called, this can be a string with spaces that describes the meaning of the callback
44287497f52SBarry Smith - PetscFunction - user provided callback function that returns an error code
443ef1023bdSBarry Smith 
444ef1023bdSBarry Smith   Example Usage:
445ef1023bdSBarry Smith .vb
446792fecdfSBarry Smith   PetscCallBack("XXX callback to do something", a->callback(...));
447ef1023bdSBarry Smith .ve
448ef1023bdSBarry Smith 
449ef1023bdSBarry Smith   Level: developer
450ef1023bdSBarry Smith 
451667f096bSBarry Smith   Notes:
4527e6f8dd6SBarry Smith   `PetscUseTypeMethod()` and ` PetscTryTypeMethod()` are the preferred API for this functionality. But when the callback functions are associated with a
4537e6f8dd6SBarry Smith   `DMSNES` or `DMTS` this API must be used.
4547e6f8dd6SBarry Smith 
455667f096bSBarry Smith   Once the error handler is called the calling function is then returned from with the given
456667f096bSBarry Smith   error code. Experienced users can set the error handler with `PetscPushErrorHandler()`.
457667f096bSBarry Smith 
458667f096bSBarry Smith   `PetscCallBack()` should only be called in PETSc when a call is being made to a user provided call-back routine.
459667f096bSBarry Smith 
4607e6f8dd6SBarry Smith   Developer Note:
4617e6f8dd6SBarry 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
4627e6f8dd6SBarry Smith 
46387497f52SBarry Smith .seealso: `SETERRQ()`, `PetscCheck()`, `PetscCall()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`
4647e6f8dd6SBarry Smith           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`, `CHKERRMPI()`, `PetscCall()`,  `PetscUseTypeMethod()`, `PetscTryTypeMethod()`
465ef1023bdSBarry Smith M*/
466ef1023bdSBarry Smith 
4673ba16761SJacob Faibussowitsch /*MC
4688af9f190SBarry Smith   PetscCallVoid - Like `PetscCall()` but for use in functions that return `void`
4693ba16761SJacob Faibussowitsch 
4703ba16761SJacob Faibussowitsch   Synopsis:
4717233ce55SJed Brown   #include <petscsys.h>
4728af9f190SBarry Smith   void PetscCallVoid(PetscFunction(args))
4733ba16761SJacob Faibussowitsch 
4747cdbe19fSJose E. Roman   Not Collective; No Fortran Support
4753ba16761SJacob Faibussowitsch 
4763ba16761SJacob Faibussowitsch   Input Parameter:
4773ba16761SJacob Faibussowitsch . PetscFunction - any PETSc function that returns an error code
4783ba16761SJacob Faibussowitsch 
4793ba16761SJacob Faibussowitsch   Example Usage:
4803ba16761SJacob Faibussowitsch .vb
4813ba16761SJacob Faibussowitsch   void foo()
4823ba16761SJacob Faibussowitsch   {
4833ba16761SJacob Faibussowitsch     KSP ksp;
4843ba16761SJacob Faibussowitsch 
4853ba16761SJacob Faibussowitsch     PetscFunctionBeginUser;
4863ba16761SJacob Faibussowitsch     // OK, properly handles PETSc error codes
4873ba16761SJacob Faibussowitsch     PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp));
4888af9f190SBarry Smith     PetscFunctionReturnVoid();
4893ba16761SJacob Faibussowitsch   }
4903ba16761SJacob Faibussowitsch 
4913ba16761SJacob Faibussowitsch   PetscErrorCode bar()
4923ba16761SJacob Faibussowitsch   {
4933ba16761SJacob Faibussowitsch     KSP ksp;
4943ba16761SJacob Faibussowitsch 
4953ba16761SJacob Faibussowitsch     PetscFunctionBeginUser;
4963ba16761SJacob Faibussowitsch     // ERROR, Non-void function 'bar' should return a value
4973ba16761SJacob Faibussowitsch     PetscCallVoid(KSPCreate(PETSC_COMM_WORLD, &ksp));
4983ba16761SJacob Faibussowitsch     // OK, returning PetscErrorCode
4993ba16761SJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
5003ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
5013ba16761SJacob Faibussowitsch   }
502667f096bSBarry Smith .ve
5033ba16761SJacob Faibussowitsch 
5043ba16761SJacob Faibussowitsch   Level: beginner
5053ba16761SJacob Faibussowitsch 
506667f096bSBarry Smith   Notes:
507667f096bSBarry Smith   Has identical usage to `PetscCall()`, except that it returns `void` on error instead of a
508667f096bSBarry Smith   `PetscErrorCode`. See `PetscCall()` for more detailed discussion.
509667f096bSBarry Smith 
510667f096bSBarry Smith   Note that users should prefer `PetscCallAbort()` to this routine. While this routine does
511667f096bSBarry Smith   "handle" errors by returning from the enclosing function, it effectively gobbles the
512667f096bSBarry Smith   error. Since the enclosing function itself returns `void`, its callers have no way of knowing
513667f096bSBarry Smith   that the routine returned early due to an error. `PetscCallAbort()` at least ensures that the
514667f096bSBarry Smith   program crashes gracefully.
515667f096bSBarry Smith 
516648c30bcSBarry Smith .seealso: `PetscCall()`, `PetscErrorCode`, `PetscCallAbort()`, `PetscCallNull()`
5173ba16761SJacob Faibussowitsch M*/
5187c5b2466SBarry Smith 
5197c5b2466SBarry Smith /*MC
5207c5b2466SBarry Smith   PetscCallReturnMPI - Calls a PETSc function and then checks the resulting error code, if it is
5217c5b2466SBarry Smith   non-zero it calls the error handler and returns from the current function with an MPI error code.
5227c5b2466SBarry Smith   To check for errors in PETSc provided MPI callbacks.
5237c5b2466SBarry Smith 
5247c5b2466SBarry Smith   Synopsis:
5257233ce55SJed Brown   #include <petscsys.h>
5267c5b2466SBarry Smith   void PetscCallReturnMPI(PetscFunction(args))
5277c5b2466SBarry Smith 
5287c5b2466SBarry Smith   Not Collective
5297c5b2466SBarry Smith 
5307c5b2466SBarry Smith   Input Parameter:
5317c5b2466SBarry Smith . PetscFunction - any PETSc function that returns an error code
5327c5b2466SBarry Smith 
5337c5b2466SBarry Smith   Level: advanced
5347c5b2466SBarry Smith 
5357c5b2466SBarry Smith   Notes:
5367c5b2466SBarry Smith   Note to be confused with `PetscCallMPI()`.
5377c5b2466SBarry Smith 
5387c5b2466SBarry Smith   This is be used in a PETSc-provided MPI callback function, such as `MPI_Comm_delete_attr_function function()`.
5397c5b2466SBarry Smith 
5407c5b2466SBarry Smith   Currently, it always returns `MPI_ERR_OTHER` on failure
5417c5b2466SBarry Smith 
5427c5b2466SBarry Smith .seealso: `PetscCall()`, `PetscCallMPI()`, `SETERRQ()`, `PetscCheck()`, `PetscAssert()`, `PetscTraceBackErrorHandler()`, `PetscCallMPI()`,
5437c5b2466SBarry Smith           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`, `CHKERRA()`,
5447c5b2466SBarry Smith           `CHKERRMPI()`, `PetscCallBack()`, `PetscCallAbort()`, `PetscCallVoid()`, `PetscCallNull()`
5457c5b2466SBarry Smith M*/
5467c5b2466SBarry Smith 
5479566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
5489566063dSJacob Faibussowitsch void PetscCall(PetscErrorCode);
549792fecdfSBarry Smith void PetscCallBack(const char *, PetscErrorCode);
5509566063dSJacob Faibussowitsch void PetscCallVoid(PetscErrorCode);
551648c30bcSBarry Smith void PetscCallNull(PetscErrorCode);
5527c5b2466SBarry Smith void PetscCallReturnMPI(PetscErrorCode);
5539566063dSJacob Faibussowitsch #else
5549371c9d4SSatish Balay   #define PetscCall(...) \
5559371c9d4SSatish Balay     do { \
5563ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_q_; \
55737154d10SBarry Smith       PetscStackUpdateLine; \
5583ba16761SJacob Faibussowitsch       ierr_petsc_call_q_ = __VA_ARGS__; \
5593ba16761SJacob 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, " "); \
5609566063dSJacob Faibussowitsch     } while (0)
561648c30bcSBarry Smith   #define PetscCallNull(...) \
562648c30bcSBarry Smith     do { \
563648c30bcSBarry Smith       PetscErrorCode ierr_petsc_call_q_; \
564648c30bcSBarry Smith       PetscStackUpdateLine; \
565648c30bcSBarry Smith       ierr_petsc_call_q_ = __VA_ARGS__; \
566648c30bcSBarry Smith       if (PetscUnlikely(ierr_petsc_call_q_ != PETSC_SUCCESS)) { \
567648c30bcSBarry Smith         (void)PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); \
568648c30bcSBarry Smith         PetscFunctionReturn(NULL); \
569648c30bcSBarry Smith       } \
570648c30bcSBarry Smith     } while (0)
5719371c9d4SSatish Balay   #define PetscCallBack(function, ...) \
5729371c9d4SSatish Balay     do { \
5733ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_q_; \
574e33ced7fSLisandro Dalcin       PetscStackUpdateLine; \
575e33ced7fSLisandro Dalcin       PetscStackPushExternal(function); \
5763ba16761SJacob Faibussowitsch       ierr_petsc_call_q_ = __VA_ARGS__; \
577ef1023bdSBarry Smith       PetscStackPop; \
5783ba16761SJacob 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, " "); \
579ef1023bdSBarry Smith     } while (0)
5809371c9d4SSatish Balay   #define PetscCallVoid(...) \
5819371c9d4SSatish Balay     do { \
5823ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_void_; \
583e33ced7fSLisandro Dalcin       PetscStackUpdateLine; \
5843ba16761SJacob Faibussowitsch       ierr_petsc_call_void_ = __VA_ARGS__; \
5853ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_void_ != PETSC_SUCCESS)) { \
586dd460d27SBarry Smith         (void)PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_void_, PETSC_ERROR_REPEAT, " "); \
5879371c9d4SSatish Balay         return; \
5889371c9d4SSatish Balay       } \
5899566063dSJacob Faibussowitsch     } while (0)
5907c5b2466SBarry Smith   #define PetscCallReturnMPI(...) \
5917c5b2466SBarry Smith     do { \
5927c5b2466SBarry Smith       PetscErrorCode ierr_petsc_call_q_; \
5937c5b2466SBarry Smith       PetscStackUpdateLine; \
5947c5b2466SBarry Smith       ierr_petsc_call_q_ = __VA_ARGS__; \
5957c5b2466SBarry Smith       if (PetscUnlikely(ierr_petsc_call_q_ != PETSC_SUCCESS)) { \
5967c5b2466SBarry Smith         (void)PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_q_, PETSC_ERROR_REPEAT, " "); \
5977c5b2466SBarry Smith         return MPI_ERR_OTHER; \
5987c5b2466SBarry Smith       } \
5997c5b2466SBarry Smith     } while (0)
6009566063dSJacob Faibussowitsch #endif
6019566063dSJacob Faibussowitsch 
6029566063dSJacob Faibussowitsch /*MC
6039566063dSJacob Faibussowitsch   CHKERRQ - Checks error code returned from PETSc function
60430de9b25SBarry Smith 
60530de9b25SBarry Smith   Synopsis:
606aaa7dc30SBarry Smith   #include <petscsys.h>
6079566063dSJacob Faibussowitsch   void CHKERRQ(PetscErrorCode ierr)
6089566063dSJacob Faibussowitsch 
6099566063dSJacob Faibussowitsch   Not Collective
6109566063dSJacob Faibussowitsch 
6112fe279fdSBarry Smith   Input Parameter:
6129566063dSJacob Faibussowitsch . ierr - nonzero error code
6139566063dSJacob Faibussowitsch 
6149566063dSJacob Faibussowitsch   Level: deprecated
6159566063dSJacob Faibussowitsch 
616667f096bSBarry Smith   Note:
617667f096bSBarry Smith   Deprecated in favor of `PetscCall()`. This routine behaves identically to it.
618667f096bSBarry Smith 
619db781477SPatrick Sanan .seealso: `PetscCall()`
6209566063dSJacob Faibussowitsch M*/
6219566063dSJacob Faibussowitsch #define CHKERRQ(...) PetscCall(__VA_ARGS__)
6229566063dSJacob Faibussowitsch #define CHKERRV(...) PetscCallVoid(__VA_ARGS__)
6239566063dSJacob Faibussowitsch 
624a1fd7ae3SBarry Smith PETSC_EXTERN void PetscMPIErrorString(PetscMPIInt, size_t, char *);
625db9cea48SBarry Smith 
6269566063dSJacob Faibussowitsch /*MC
6279566063dSJacob Faibussowitsch   PetscCallMPI - Checks error code returned from MPI calls, if non-zero it calls the error
628648c30bcSBarry Smith   handler and then returns a `PetscErrorCode`
6299566063dSJacob Faibussowitsch 
6309566063dSJacob Faibussowitsch   Synopsis:
6317233ce55SJed Brown   #include <petscsys.h>
63249c86fc7SBarry Smith   void PetscCallMPI(MPI_Function(args))
63330de9b25SBarry Smith 
634eca87e8dSBarry Smith   Not Collective
63530de9b25SBarry Smith 
6362fe279fdSBarry Smith   Input Parameter:
63749c86fc7SBarry Smith . MPI_Function - an MPI function that returns an MPI error code
63830de9b25SBarry Smith 
639667f096bSBarry Smith   Level: beginner
640667f096bSBarry Smith 
6419566063dSJacob Faibussowitsch   Notes:
64287497f52SBarry Smith   Always returns the error code `PETSC_ERR_MPI`; the MPI error code and string are embedded in
6439566063dSJacob Faibussowitsch   the string error message. Do not use this to call any other routines (for example PETSc
6443ba16761SJacob Faibussowitsch   routines), it should only be used for direct MPI calls. The user may configure PETSc with the
6453ba16761SJacob Faibussowitsch   `--with-strict-petscerrorcode` option to check this at compile-time, otherwise they must
6469566063dSJacob Faibussowitsch   check this themselves.
6479566063dSJacob Faibussowitsch 
648aaa8cc7dSPierre Jolivet   This routine can only be used in functions returning `PetscErrorCode` themselves. If the
6493ba16761SJacob Faibussowitsch   calling function returns a different type, use `PetscCallMPIAbort()` instead.
6503ba16761SJacob Faibussowitsch 
6519566063dSJacob Faibussowitsch   Example Usage:
6529566063dSJacob Faibussowitsch .vb
6539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(...)); // OK, calling MPI function
6549566063dSJacob Faibussowitsch 
6559566063dSJacob Faibussowitsch   PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead!
6569566063dSJacob Faibussowitsch .ve
6579566063dSJacob Faibussowitsch 
65849c86fc7SBarry Smith   Fortran Notes:
65987497f52SBarry Smith   The Fortran function from which this is used must declare a variable `PetscErrorCode` ierr and ierr must be
66049c86fc7SBarry Smith   the final argument to the MPI function being called.
66149c86fc7SBarry Smith 
66249c86fc7SBarry Smith   In the main program and in Fortran subroutines that do not have ierr as the final return parameter one
66387497f52SBarry Smith   should use `PetscCallMPIA()`
66449c86fc7SBarry Smith 
66549c86fc7SBarry Smith   Fortran Usage:
66649c86fc7SBarry Smith .vb
66749c86fc7SBarry Smith   PetscErrorCode ierr or integer ierr
66849c86fc7SBarry Smith   ...
66949c86fc7SBarry Smith   PetscCallMPI(MPI_Comm_size(...,ierr))
67049c86fc7SBarry Smith   PetscCallMPIA(MPI_Comm_size(...,ierr)) ! Will abort after calling error handler
67149c86fc7SBarry Smith 
67249c86fc7SBarry Smith   PetscCallMPI(MPI_Comm_size(...,eflag)) ! ERROR, final argument must be ierr
67349c86fc7SBarry Smith .ve
67449c86fc7SBarry Smith 
6757c5b2466SBarry Smith .seealso: `SETERRMPI()`, `PetscCall()`, `PetscCallReturnMPI()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`,
6767c5b2466SBarry Smith           `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
6777c5b2466SBarry Smith           `PetscError()`, `CHKMEMQ`, `PetscCallMPINull()`
6787c5b2466SBarry Smith M*/
6797c5b2466SBarry Smith 
6807c5b2466SBarry Smith /*MC
6817c5b2466SBarry Smith   PetscCallMPIReturnMPI - Checks error code returned from MPI calls, if non-zero it calls the error
6827c5b2466SBarry Smith   handler and then returns an MPI error code. To check for errors in PETSc-provided MPI callbacks.
6837c5b2466SBarry Smith 
6847c5b2466SBarry Smith   Synopsis:
6857233ce55SJed Brown   #include <petscsys.h>
6867c5b2466SBarry Smith   void PetscCallMPIReturnMPI(MPI_Function(args))
6877c5b2466SBarry Smith 
6887c5b2466SBarry Smith   Not Collective
6897c5b2466SBarry Smith 
6907c5b2466SBarry Smith   Input Parameter:
6917c5b2466SBarry Smith . MPI_Function - an MPI function that returns an MPI error code
6927c5b2466SBarry Smith 
6937c5b2466SBarry Smith   Level: advanced
6947c5b2466SBarry Smith 
6957c5b2466SBarry Smith .seealso: `SETERRMPI()`, `PetscCall()`, `PetscCallMPI()`, `PetscCallReturnMPI()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`,
6963ba16761SJacob Faibussowitsch           `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
697648c30bcSBarry Smith           `PetscError()`, `CHKMEMQ`, `PetscCallMPINull()`
698648c30bcSBarry Smith M*/
699648c30bcSBarry Smith 
700648c30bcSBarry Smith /*MC
701648c30bcSBarry Smith   PetscCallMPINull - Checks error code returned from MPI calls, if non-zero it calls the error
702648c30bcSBarry Smith   handler and then returns a `NULL`
703648c30bcSBarry Smith 
704648c30bcSBarry Smith   Synopsis:
7057233ce55SJed Brown   #include <petscsys.h>
706648c30bcSBarry Smith   void PetscCallMPINull(MPI_Function(args))
707648c30bcSBarry Smith 
708648c30bcSBarry Smith   Not Collective; No Fortran Support
709648c30bcSBarry Smith 
710648c30bcSBarry Smith   Input Parameter:
711648c30bcSBarry Smith . MPI_Function - an MPI function that returns an MPI error code
712648c30bcSBarry Smith 
713648c30bcSBarry Smith   Level: beginner
714648c30bcSBarry Smith 
715648c30bcSBarry Smith   Notes:
716648c30bcSBarry Smith   Always passes the error code `PETSC_ERR_MPI` to the error handler `PetscError()`; the MPI error code and string are embedded in
717648c30bcSBarry Smith   the string error message. Do not use this to call any other routines (for example PETSc
718648c30bcSBarry Smith   routines), it should only be used for direct MPI calls.
719648c30bcSBarry Smith 
720648c30bcSBarry Smith   This routine can only be used in functions returning anything that can be returned as a `NULL` themselves. If the
721648c30bcSBarry Smith   calling function returns a different type, use `PetscCallMPIAbort()` instead.
722648c30bcSBarry Smith 
723648c30bcSBarry Smith   Example Usage:
724648c30bcSBarry Smith .vb
725648c30bcSBarry Smith   PetscCallMPINull(MPI_Comm_size(...)); // OK, calling MPI function
726648c30bcSBarry Smith 
727648c30bcSBarry Smith   PetscCallMPI(PetscFunction(...)); // ERROR, use PetscCall() instead!
728648c30bcSBarry Smith .ve
729648c30bcSBarry Smith 
730648c30bcSBarry Smith .seealso: `SETERRMPI()`, `PetscCall()`, `SETERRQ()`, `SETERRABORT()`, `PetscCallAbort()`,
731648c30bcSBarry Smith           `PetscCallMPIAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
732648c30bcSBarry Smith           `PetscError()`, `CHKMEMQ`, `PetscCallMPI()`
7333ba16761SJacob Faibussowitsch M*/
7343ba16761SJacob Faibussowitsch 
7353ba16761SJacob Faibussowitsch /*MC
7363ba16761SJacob Faibussowitsch   PetscCallMPIAbort - Like `PetscCallMPI()` but calls `MPI_Abort()` on error
7373ba16761SJacob Faibussowitsch 
7383ba16761SJacob Faibussowitsch   Synopsis:
7397233ce55SJed Brown   #include <petscsys.h>
7403ba16761SJacob Faibussowitsch   void PetscCallMPIAbort(MPI_Comm comm, MPI_Function(args))
7413ba16761SJacob Faibussowitsch 
7423ba16761SJacob Faibussowitsch   Not Collective
7433ba16761SJacob Faibussowitsch 
7443ba16761SJacob Faibussowitsch   Input Parameters:
7453ba16761SJacob Faibussowitsch + comm         - the MPI communicator to abort on
7463ba16761SJacob Faibussowitsch - MPI_Function - an MPI function that returns an MPI error code
7473ba16761SJacob Faibussowitsch 
748667f096bSBarry Smith   Level: beginner
749667f096bSBarry Smith 
7503ba16761SJacob Faibussowitsch   Notes:
7513ba16761SJacob Faibussowitsch   Usage is identical to `PetscCallMPI()`. See `PetscCallMPI()` for detailed discussion.
7523ba16761SJacob Faibussowitsch 
7533ba16761SJacob Faibussowitsch   This routine may be used in functions returning `void` or other non-`PetscErrorCode` types.
7543ba16761SJacob Faibussowitsch 
755989c49a2SBarry Smith   Fortran Note:
756989c49a2SBarry Smith   In Fortran this is called `PetscCallMPIA()` and is intended to be used in the main program while `PetscCallMPI()` is
757989c49a2SBarry Smith   used in Fortran subroutines.
758989c49a2SBarry Smith 
759989c49a2SBarry Smith   Developer Note:
760989c49a2SBarry Smith   This should have the same name in Fortran.
761989c49a2SBarry Smith 
7623ba16761SJacob Faibussowitsch .seealso: `PetscCallMPI()`, `PetscCallAbort()`, `SETERRABORT()`
76330de9b25SBarry Smith M*/
7643fcd9f07SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
76563a3b9bcSJacob Faibussowitsch void PetscCallMPI(PetscMPIInt);
7663ba16761SJacob Faibussowitsch void PetscCallMPIAbort(MPI_Comm, PetscMPIInt);
767648c30bcSBarry Smith void PetscCallMPINull(PetscMPIInt);
768064a246eSJacob Faibussowitsch #else
7693ba16761SJacob Faibussowitsch   #define PetscCallMPI_Private(__PETSC_STACK_POP_FUNC__, __SETERR_FUNC__, __COMM__, ...) \
7709371c9d4SSatish Balay     do { \
7713ba16761SJacob Faibussowitsch       PetscMPIInt ierr_petsc_call_mpi_; \
772ef1023bdSBarry Smith       PetscStackUpdateLine; \
773792fecdfSBarry Smith       PetscStackPushExternal("MPI function"); \
774d71ae5a4SJacob Faibussowitsch       { \
7753ba16761SJacob Faibussowitsch         ierr_petsc_call_mpi_ = __VA_ARGS__; \
776d71ae5a4SJacob Faibussowitsch       } \
7773ba16761SJacob Faibussowitsch       __PETSC_STACK_POP_FUNC__; \
7783ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_mpi_ != MPI_SUCCESS)) { \
7793ba16761SJacob Faibussowitsch         char petsc_mpi_7_errorstring[2 * MPI_MAX_ERROR_STRING]; \
780a1fd7ae3SBarry Smith         PetscMPIErrorString(ierr_petsc_call_mpi_, 2 * MPI_MAX_ERROR_STRING, (char *)petsc_mpi_7_errorstring); \
781a1fd7ae3SBarry Smith         __SETERR_FUNC__(__COMM__, PETSC_ERR_MPI, "MPI error %d %s", ierr_petsc_call_mpi_, petsc_mpi_7_errorstring); \
7823fcd9f07SJacob Faibussowitsch       } \
7833fcd9f07SJacob Faibussowitsch     } while (0)
7843ba16761SJacob Faibussowitsch 
7853ba16761SJacob Faibussowitsch   #define PetscCallMPI(...)            PetscCallMPI_Private(PetscStackPop, SETERRQ, PETSC_COMM_SELF, __VA_ARGS__)
7867c5b2466SBarry Smith   #define PetscCallMPIReturnMPI(...)   PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRMPI, PETSC_COMM_SELF, __VA_ARGS__)
7873ba16761SJacob Faibussowitsch   #define PetscCallMPIAbort(comm, ...) PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRABORT, comm, __VA_ARGS__)
788648c30bcSBarry Smith   #define PetscCallMPINull(...)        PetscCallMPI_Private(PetscStackPopNoCheck(PETSC_FUNCTION_NAME), SETERRQNULL, PETSC_COMM_SELF, __VA_ARGS__)
7899566063dSJacob Faibussowitsch #endif
790064a246eSJacob Faibussowitsch 
7917037b0edSPatrick Sanan /*MC
7922acb36afSStefano Zampini   CHKERRMPI - Checks the error code returned from MPI calls, if different from `MPI_SUCCESS` it calls the error handler and then returns
7939566063dSJacob Faibussowitsch 
7949566063dSJacob Faibussowitsch   Synopsis:
7957233ce55SJed Brown   #include <petscsys.h>
7962acb36afSStefano Zampini   void CHKERRMPI(PetscMPIInt ierr)
7979566063dSJacob Faibussowitsch 
7989566063dSJacob Faibussowitsch   Not Collective
7999566063dSJacob Faibussowitsch 
8009566063dSJacob Faibussowitsch   Input Parameter:
8019566063dSJacob Faibussowitsch . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
8029566063dSJacob Faibussowitsch 
8039566063dSJacob Faibussowitsch   Level: deprecated
8049566063dSJacob Faibussowitsch 
805667f096bSBarry Smith   Note:
806667f096bSBarry Smith   Deprecated in favor of `PetscCallMPI()`. This routine behaves identically to it.
807667f096bSBarry Smith 
808db781477SPatrick Sanan .seealso: `PetscCallMPI()`
8099566063dSJacob Faibussowitsch M*/
8109566063dSJacob Faibussowitsch #define CHKERRMPI(...) PetscCallMPI(__VA_ARGS__)
8119566063dSJacob Faibussowitsch 
8129566063dSJacob Faibussowitsch /*MC
813989c49a2SBarry Smith   PetscCallAbort - Checks error code returned from PETSc function, if non-zero it aborts immediately by calling `MPI_Abort()`
8149566063dSJacob Faibussowitsch 
8159566063dSJacob Faibussowitsch   Synopsis:
8167233ce55SJed Brown   #include <petscsys.h>
8179566063dSJacob Faibussowitsch   void PetscCallAbort(MPI_Comm comm, PetscErrorCode ierr)
8189566063dSJacob Faibussowitsch 
819c3339decSBarry Smith   Collective
8209566063dSJacob Faibussowitsch 
8219566063dSJacob Faibussowitsch   Input Parameters:
8229566063dSJacob Faibussowitsch + comm - the MPI communicator on which to abort
8239566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
8249566063dSJacob Faibussowitsch 
825667f096bSBarry Smith   Level: intermediate
826667f096bSBarry Smith 
8279566063dSJacob Faibussowitsch   Notes:
82887497f52SBarry Smith   This macro has identical type and usage semantics to `PetscCall()` with the important caveat
8299566063dSJacob Faibussowitsch   that this macro does not return. Instead, if ierr is nonzero it calls the PETSc error handler
83087497f52SBarry Smith   and then immediately calls `MPI_Abort()`. It can therefore be used anywhere.
8319566063dSJacob Faibussowitsch 
83287497f52SBarry Smith   As per `MPI_Abort()` semantics the communicator passed must be valid, although there is currently
83387497f52SBarry Smith   no attempt made at handling any potential errors from `MPI_Abort()`. Note that while
83487497f52SBarry Smith   `MPI_Abort()` is required to terminate only those processes which reside on comm, it is often
83587497f52SBarry Smith   the case that `MPI_Abort()` terminates *all* processes.
8369566063dSJacob Faibussowitsch 
8379566063dSJacob Faibussowitsch   Example Usage:
8389566063dSJacob Faibussowitsch .vb
8399566063dSJacob Faibussowitsch   PetscErrorCode boom(void) { return PETSC_ERR_MEM; }
8409566063dSJacob Faibussowitsch 
8419566063dSJacob Faibussowitsch   void foo(void)
8429566063dSJacob Faibussowitsch   {
8439566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type
8449566063dSJacob Faibussowitsch   }
8459566063dSJacob Faibussowitsch 
8469566063dSJacob Faibussowitsch   double bar(void)
8479566063dSJacob Faibussowitsch   {
8489566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_WORLD,boom()); // OK, does not return a type
8499566063dSJacob Faibussowitsch   }
8509566063dSJacob Faibussowitsch 
8519566063dSJacob Faibussowitsch   PetscCallAbort(MPI_COMM_NULL,boom()); // ERROR, communicator should be valid
8529566063dSJacob Faibussowitsch 
8539566063dSJacob Faibussowitsch   struct baz
8549566063dSJacob Faibussowitsch   {
8559566063dSJacob Faibussowitsch     baz()
8569566063dSJacob Faibussowitsch     {
8579566063dSJacob Faibussowitsch       PetscCallAbort(PETSC_COMM_SELF,boom()); // OK
8589566063dSJacob Faibussowitsch     }
8599566063dSJacob Faibussowitsch 
8609566063dSJacob Faibussowitsch     ~baz()
8619566063dSJacob Faibussowitsch     {
8629566063dSJacob Faibussowitsch       PetscCallAbort(PETSC_COMM_SELF,boom()); // OK (in fact the only way to handle PETSc errors)
8639566063dSJacob Faibussowitsch     }
8649566063dSJacob Faibussowitsch   };
8659566063dSJacob Faibussowitsch .ve
8669566063dSJacob Faibussowitsch 
867989c49a2SBarry Smith   Fortran Note:
868989c49a2SBarry Smith   Use `PetscCallA()`.
869989c49a2SBarry Smith 
870989c49a2SBarry Smith   Developer Note:
871989c49a2SBarry Smith   This should have the same name in Fortran as in C.
872989c49a2SBarry Smith 
873db781477SPatrick Sanan .seealso: `SETERRABORT()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`,
8745687f061SJacob Faibussowitsch           `SETERRQ()`, `CHKMEMQ`, `PetscCallMPI()`, `PetscCallCXXAbort()`
8759566063dSJacob Faibussowitsch M*/
8769566063dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
8779566063dSJacob Faibussowitsch void PetscCallAbort(MPI_Comm, PetscErrorCode);
8789566063dSJacob Faibussowitsch void PetscCallContinue(PetscErrorCode);
8799566063dSJacob Faibussowitsch #else
8809371c9d4SSatish Balay   #define PetscCallAbort(comm, ...) \
8819371c9d4SSatish Balay     do { \
8827da6b3ccSLisandro Dalcin       PetscErrorCode ierr_petsc_call_abort_; \
8833ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
8847da6b3ccSLisandro Dalcin       ierr_petsc_call_abort_ = __VA_ARGS__; \
8853ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_abort_ != PETSC_SUCCESS)) { \
8863ba16761SJacob Faibussowitsch         ierr_petsc_call_abort_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_abort_, PETSC_ERROR_REPEAT, " "); \
8873ba16761SJacob Faibussowitsch         (void)MPI_Abort(comm, (PetscMPIInt)ierr_petsc_call_abort_); \
8889566063dSJacob Faibussowitsch       } \
8899566063dSJacob Faibussowitsch     } while (0)
8909371c9d4SSatish Balay   #define PetscCallContinue(...) \
8919371c9d4SSatish Balay     do { \
8927da6b3ccSLisandro Dalcin       PetscErrorCode ierr_petsc_call_continue_; \
8933ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
8947da6b3ccSLisandro Dalcin       ierr_petsc_call_continue_ = __VA_ARGS__; \
8953ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_call_continue_ != PETSC_SUCCESS)) { \
8963ba16761SJacob Faibussowitsch         ierr_petsc_call_continue_ = PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_call_continue_, PETSC_ERROR_REPEAT, " "); \
8973ba16761SJacob Faibussowitsch         (void)ierr_petsc_call_continue_; \
8983ba16761SJacob Faibussowitsch       } \
8999566063dSJacob Faibussowitsch     } while (0)
9009566063dSJacob Faibussowitsch #endif
9019566063dSJacob Faibussowitsch 
9029566063dSJacob Faibussowitsch /*MC
9039566063dSJacob Faibussowitsch   CHKERRABORT - Checks error code returned from PETSc function. If non-zero it aborts immediately.
9049566063dSJacob Faibussowitsch 
9059566063dSJacob Faibussowitsch   Synopsis:
9067233ce55SJed Brown   #include <petscsys.h>
9079566063dSJacob Faibussowitsch   void CHKERRABORT(MPI_Comm comm, PetscErrorCode ierr)
9089566063dSJacob Faibussowitsch 
9099566063dSJacob Faibussowitsch   Not Collective
9109566063dSJacob Faibussowitsch 
9119566063dSJacob Faibussowitsch   Input Parameters:
9129566063dSJacob Faibussowitsch + comm - the MPI communicator
9139566063dSJacob Faibussowitsch - ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
9149566063dSJacob Faibussowitsch 
9159566063dSJacob Faibussowitsch   Level: deprecated
9169566063dSJacob Faibussowitsch 
917667f096bSBarry Smith   Note:
918667f096bSBarry Smith   Deprecated in favor of `PetscCallAbort()`. This routine behaves identically to it.
919667f096bSBarry Smith 
920af27ebaaSBarry Smith .seealso: `PetscCallAbort()`, `PetscErrorCode`
9219566063dSJacob Faibussowitsch M*/
9229566063dSJacob Faibussowitsch #define CHKERRABORT(comm, ...) PetscCallAbort(comm, __VA_ARGS__)
9239566063dSJacob Faibussowitsch #define CHKERRCONTINUE(...)    PetscCallContinue(__VA_ARGS__)
9249566063dSJacob Faibussowitsch 
9259566063dSJacob Faibussowitsch /*MC
92687497f52SBarry Smith    CHKERRA - Fortran-only replacement for use of `CHKERRQ()` in the main program, which aborts immediately
927f388eb8bSPatrick Sanan 
928f388eb8bSPatrick Sanan    Synopsis:
929f388eb8bSPatrick Sanan    #include <petscsys.h>
930f388eb8bSPatrick Sanan    PetscErrorCode CHKERRA(PetscErrorCode ierr)
931f388eb8bSPatrick Sanan 
932f388eb8bSPatrick Sanan    Not Collective
933f388eb8bSPatrick Sanan 
9342fe279fdSBarry Smith    Input Parameter:
935f388eb8bSPatrick Sanan .  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
936f388eb8bSPatrick Sanan 
93787497f52SBarry Smith    Level: deprecated
938f388eb8bSPatrick Sanan 
93987497f52SBarry Smith    Note:
94087497f52SBarry Smith    This macro is rarely needed, normal usage is `PetscCallA()` in the main Fortran program.
941f388eb8bSPatrick Sanan 
942989c49a2SBarry Smith    Developer Note:
943989c49a2SBarry Smith    Why isn't this named `CHKERRABORT()` in Fortran?
944989c49a2SBarry Smith 
94587497f52SBarry Smith .seealso: `PetscCall()`, `PetscCallA()`, `PetscCallAbort()`, `CHKERRQ()`, `SETERRA()`, `SETERRQ()`, `SETERRABORT()`
946f388eb8bSPatrick Sanan M*/
947f388eb8bSPatrick Sanan 
9481e4f893aSSatish Balay PETSC_EXTERN PetscBool petscwaitonerrorflg;
9491e4f893aSSatish Balay PETSC_EXTERN PetscBool petscindebugger;
95035f00c14SToby Isaac PETSC_EXTERN PetscBool petscabortmpifinalize;
9517c66cc67SJunchao Zhang 
9521fc6d13cSAlexander #if defined(PETSC_CLANG_STATIC_ANALYZER)
9531fc6d13cSAlexander void PETSCABORTWITHERR_Private(MPI_Comm, PetscErrorCode);
9541fc6d13cSAlexander #else
9551fc6d13cSAlexander   #define PETSCABORTWITHIERR_Private(comm, ierr) \
9561fc6d13cSAlexander     do { \
9571fc6d13cSAlexander       PetscMPIInt size_; \
958dd460d27SBarry Smith       (void)MPI_Comm_size(comm, &size_); \
9591fc6d13cSAlexander       if (PetscCIEnabledPortableErrorOutput && (size_ == PetscGlobalSize || petscabortmpifinalize) && ierr != PETSC_ERR_SIG) { \
960dd460d27SBarry Smith         (void)MPI_Finalize(); \
9611fc6d13cSAlexander         exit(0); \
9621fc6d13cSAlexander       } else if (PetscCIEnabledPortableErrorOutput && PetscGlobalSize == 1) { \
9631fc6d13cSAlexander         exit(0); \
9641fc6d13cSAlexander       } else { \
965dd460d27SBarry Smith         (void)MPI_Abort(comm, (PetscMPIInt)ierr); \
9661fc6d13cSAlexander       } \
9671fc6d13cSAlexander     } while (0)
9681fc6d13cSAlexander #endif
9691fc6d13cSAlexander 
9707c66cc67SJunchao Zhang /*MC
971667f096bSBarry Smith    PETSCABORT - Call `MPI_Abort()` with an informative error code
9727c66cc67SJunchao Zhang 
9737c66cc67SJunchao Zhang    Synopsis:
9747c66cc67SJunchao Zhang    #include <petscsys.h>
9757c66cc67SJunchao Zhang    PETSCABORT(MPI_Comm comm, PetscErrorCode ierr)
9767c66cc67SJunchao Zhang 
9777cdbe19fSJose E. Roman    Collective; No Fortran Support
9787c66cc67SJunchao Zhang 
9797c66cc67SJunchao Zhang    Input Parameters:
98095bd0b28SBarry Smith +  comm - An MPI communicator, so that the error can be collective
9817c66cc67SJunchao Zhang -  ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
9827c66cc67SJunchao Zhang 
983bf4d2887SBarry Smith    Level: advanced
9847c66cc67SJunchao Zhang 
985bf4d2887SBarry Smith    Notes:
986989c49a2SBarry Smith    If the option `-start_in_debugger` was used then this calls `abort()` to stop the program in the debugger.
987bf4d2887SBarry Smith 
988989c49a2SBarry Smith    if `PetscCIEnabledPortableErrorOutput` is set, which means the code is running in the PETSc test harness (make test),
989989c49a2SBarry Smith    and `comm` is `MPI_COMM_WORLD` it strives to exit cleanly without calling `MPI_Abort()` and instead calling `MPI_Finalize()`.
990660278c0SBarry Smith 
991989c49a2SBarry Smith    This is currently only used when an error propagates up to the C `main()` program and is detected by a `PetscCall()`, `PetscCallMPI()`,
992989c49a2SBarry Smith    or is set in `main()` with `SETERRQ()`. Abort calls such as `SETERRABORT()`,
993989c49a2SBarry Smith    `PetscCheckAbort()`, `PetscCallMPIAbort()`, and `PetscCallAbort()` always call `MPI_Abort()` and do not have any special
994989c49a2SBarry Smith    handling for the test harness.
995989c49a2SBarry Smith 
996989c49a2SBarry Smith    Developer Note:
997989c49a2SBarry Smith    Should the other abort calls also pass through this call instead of calling `MPI_Abort()` directly?
998989c49a2SBarry Smith 
999989c49a2SBarry Smith .seealso: `PetscError()`, `PetscCall()`, `SETERRABORT()`, `PetscCheckAbort()`, `PetscCallMPIAbort()`, `PetscCall()`, `PetscCallMPI()`,
1000af27ebaaSBarry Smith           `PetscCallAbort()`, `MPI_Abort()`, `PetscErrorCode`
10017c66cc67SJunchao Zhang M*/
100229f88feaSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
100329f88feaSJacob Faibussowitsch void PETSCABORT(MPI_Comm, PetscErrorCode);
100429f88feaSJacob Faibussowitsch #else
10059371c9d4SSatish Balay   #define PETSCABORT(comm, ...) \
10069371c9d4SSatish Balay     do { \
10073ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_abort_; \
10080b36fe0aSPierre Jolivet       if (petscwaitonerrorflg) ierr_petsc_abort_ = PetscSleep(1000); \
10090b36fe0aSPierre Jolivet       if (petscindebugger) abort(); \
10100b36fe0aSPierre Jolivet       else { \
10113ba16761SJacob Faibussowitsch         ierr_petsc_abort_ = __VA_ARGS__; \
10121fc6d13cSAlexander         PETSCABORTWITHIERR_Private(comm, ierr_petsc_abort_); \
10133fcd9f07SJacob Faibussowitsch       } \
10147c66cc67SJunchao Zhang     } while (0)
101529f88feaSJacob Faibussowitsch #endif
1016986eef2eSBarry Smith 
10179566063dSJacob Faibussowitsch #ifdef PETSC_CLANGUAGE_CXX
1018986eef2eSBarry Smith   /*MC
10199566063dSJacob Faibussowitsch   PetscCallThrow - Checks error code, if non-zero it calls the C++ error handler which throws
10209566063dSJacob Faibussowitsch   an exception
1021986eef2eSBarry Smith 
1022986eef2eSBarry Smith   Synopsis:
10237233ce55SJed Brown   #include <petscsys.h>
10249566063dSJacob Faibussowitsch   void PetscCallThrow(PetscErrorCode ierr)
1025986eef2eSBarry Smith 
1026986eef2eSBarry Smith   Not Collective
1027986eef2eSBarry Smith 
10289566063dSJacob Faibussowitsch   Input Parameter:
1029986eef2eSBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
1030986eef2eSBarry Smith 
1031667f096bSBarry Smith   Level: beginner
1032667f096bSBarry Smith 
1033986eef2eSBarry Smith   Notes:
1034af27ebaaSBarry Smith   Requires PETSc to be configured with clanguage of c++. Throws a std::runtime_error() on error.
1035986eef2eSBarry Smith 
103687497f52SBarry Smith   Once the error handler throws the exception you can use `PetscCallVoid()` which returns without
103787497f52SBarry Smith   an error code (bad idea since the error is ignored) or `PetscCallAbort()` to have `MPI_Abort()`
10389566063dSJacob Faibussowitsch   called immediately.
10399566063dSJacob Faibussowitsch 
1040db781477SPatrick Sanan .seealso: `SETERRQ()`, `PetscCall()`, `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`,
1041db781477SPatrick Sanan           `PetscPushErrorHandler()`, `PetscError()`, `CHKMEMQ`
1042986eef2eSBarry Smith M*/
10439371c9d4SSatish Balay   #define PetscCallThrow(...) \
10449371c9d4SSatish Balay     do { \
10453ba16761SJacob Faibussowitsch       PetscStackUpdateLine; \
10463ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_call_throw_ = __VA_ARGS__; \
10473ba16761SJacob 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); \
1048ffc4695bSBarry Smith     } while (0)
104985614651SBarry Smith 
1050cc26af49SMatthew Knepley   /*MC
1051cc26af49SMatthew Knepley   CHKERRXX - Checks error code, if non-zero it calls the C++ error handler which throws an exception
1052cc26af49SMatthew Knepley 
1053cc26af49SMatthew Knepley   Synopsis:
10547233ce55SJed Brown   #include <petscsys.h>
10553af045c5SBarry Smith   void CHKERRXX(PetscErrorCode ierr)
1056cc26af49SMatthew Knepley 
1057eca87e8dSBarry Smith   Not Collective
1058cc26af49SMatthew Knepley 
10599566063dSJacob Faibussowitsch   Input Parameter:
10603af045c5SBarry Smith . ierr - nonzero error code, see the list of standard error codes in include/petscerror.h
1061cc26af49SMatthew Knepley 
10629566063dSJacob Faibussowitsch   Level: deprecated
1063cc26af49SMatthew Knepley 
1064667f096bSBarry Smith   Note:
1065667f096bSBarry Smith   Deprecated in favor of `PetscCallThrow()`. This routine behaves identically to it.
1066667f096bSBarry Smith 
1067db781477SPatrick Sanan .seealso: `PetscCallThrow()`
1068cc26af49SMatthew Knepley M*/
10699566063dSJacob Faibussowitsch   #define CHKERRXX(...) PetscCallThrow(__VA_ARGS__)
1070fd705b32SMatthew Knepley #endif
1071fd705b32SMatthew Knepley 
10723ba16761SJacob Faibussowitsch #define PetscCallCXX_Private(__SETERR_FUNC__, __COMM__, ...) \
10733ba16761SJacob Faibussowitsch   do { \
10743ba16761SJacob Faibussowitsch     PetscStackUpdateLine; \
10753ba16761SJacob Faibussowitsch     try { \
10763ba16761SJacob Faibussowitsch       __VA_ARGS__; \
10773ba16761SJacob Faibussowitsch     } catch (const std::exception &e) { \
10783ba16761SJacob Faibussowitsch       __SETERR_FUNC__(__COMM__, PETSC_ERR_LIB, "%s", e.what()); \
10793ba16761SJacob Faibussowitsch     } \
10803ba16761SJacob Faibussowitsch   } while (0)
10813ba16761SJacob Faibussowitsch 
10823f520e80SJunchao Zhang /*MC
10839566063dSJacob Faibussowitsch   PetscCallCXX - Checks C++ function calls and if they throw an exception, catch it and then
10849566063dSJacob Faibussowitsch   return a PETSc error code
10853f520e80SJunchao Zhang 
10863f520e80SJunchao Zhang   Synopsis:
10877233ce55SJed Brown   #include <petscsys.h>
10885687f061SJacob Faibussowitsch   void PetscCallCXX(...) noexcept;
10893f520e80SJunchao Zhang 
10903f520e80SJunchao Zhang   Not Collective
10913f520e80SJunchao Zhang 
10929566063dSJacob Faibussowitsch   Input Parameter:
10935687f061SJacob Faibussowitsch . __VA_ARGS__ - An arbitrary expression
10945687f061SJacob Faibussowitsch 
10955687f061SJacob Faibussowitsch   Level: beginner
10963f520e80SJunchao Zhang 
10973f520e80SJunchao Zhang   Notes:
10985687f061SJacob Faibussowitsch   `PetscCallCXX(...)` is a macro replacement for
10999566063dSJacob Faibussowitsch .vb
11009566063dSJacob Faibussowitsch   try {
11015687f061SJacob Faibussowitsch     __VA_ARGS__;
11029566063dSJacob Faibussowitsch   } catch (const std::exception& e) {
11039566063dSJacob Faibussowitsch     return ConvertToPetscErrorCode(e);
11049566063dSJacob Faibussowitsch   }
11059566063dSJacob Faibussowitsch .ve
11069566063dSJacob Faibussowitsch   Due to the fact that it catches any (reasonable) exception, it is essentially noexcept.
11073f520e80SJunchao Zhang 
11085687f061SJacob Faibussowitsch   If you cannot return a `PetscErrorCode` use `PetscCallCXXAbort()` instead.
11095687f061SJacob Faibussowitsch 
11109566063dSJacob Faibussowitsch   Example Usage:
11119566063dSJacob Faibussowitsch .vb
11129566063dSJacob Faibussowitsch   void foo(void) { throw std::runtime_error("error"); }
11133f520e80SJunchao Zhang 
11149566063dSJacob Faibussowitsch   void bar()
11159566063dSJacob Faibussowitsch   {
1116e8952933SJacob Faibussowitsch     PetscCallCXX(foo()); // ERROR bar() does not return PetscErrorCode
11179566063dSJacob Faibussowitsch   }
11189566063dSJacob Faibussowitsch 
11199566063dSJacob Faibussowitsch   PetscErrorCode baz()
11209566063dSJacob Faibussowitsch   {
11219566063dSJacob Faibussowitsch     PetscCallCXX(foo()); // OK
11229566063dSJacob Faibussowitsch 
11239566063dSJacob Faibussowitsch     PetscCallCXX(
11249566063dSJacob Faibussowitsch       bar();
1125d5b43468SJose E. Roman       foo(); // OK multiple statements allowed
11269566063dSJacob Faibussowitsch     );
11279566063dSJacob Faibussowitsch   }
11289566063dSJacob Faibussowitsch 
11299566063dSJacob Faibussowitsch   struct bop
11309566063dSJacob Faibussowitsch   {
11319566063dSJacob Faibussowitsch     bop()
11329566063dSJacob Faibussowitsch     {
1133e8952933SJacob Faibussowitsch       PetscCallCXX(foo()); // ERROR returns PetscErrorCode, cannot be used in constructors
11349566063dSJacob Faibussowitsch     }
11359566063dSJacob Faibussowitsch   };
11369566063dSJacob Faibussowitsch 
1137e8952933SJacob Faibussowitsch   // ERROR contains do-while, cannot be used as function-try block
11389566063dSJacob Faibussowitsch   PetscErrorCode qux() PetscCallCXX(
11399566063dSJacob Faibussowitsch     bar();
11409566063dSJacob Faibussowitsch     baz();
11419566063dSJacob Faibussowitsch     foo();
11429566063dSJacob Faibussowitsch     return 0;
11439566063dSJacob Faibussowitsch   )
11449566063dSJacob Faibussowitsch .ve
11459566063dSJacob Faibussowitsch 
11465687f061SJacob Faibussowitsch .seealso: `PetscCallCXXAbort()`, `PetscCallThrow()`, `SETERRQ()`, `PetscCall()`,
11475687f061SJacob Faibussowitsch           `SETERRABORT()`, `PetscCallAbort()`, `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`,
11485687f061SJacob Faibussowitsch           `PetscError()`, `CHKMEMQ`
11493f520e80SJunchao Zhang M*/
11503ba16761SJacob Faibussowitsch #define PetscCallCXX(...) PetscCallCXX_Private(SETERRQ, PETSC_COMM_SELF, __VA_ARGS__)
11515687f061SJacob Faibussowitsch 
11525687f061SJacob Faibussowitsch /*MC
11535687f061SJacob Faibussowitsch   PetscCallCXXAbort - Like `PetscCallCXX()` but calls `MPI_Abort()` instead of returning an
11545687f061SJacob Faibussowitsch   error-code
11555687f061SJacob Faibussowitsch 
11565687f061SJacob Faibussowitsch   Synopsis:
11577233ce55SJed Brown   #include <petscsys.h>
11585687f061SJacob Faibussowitsch   void PetscCallCXXAbort(MPI_Comm comm, ...) noexcept;
11595687f061SJacob Faibussowitsch 
116020f4b53cSBarry Smith   Collective; No Fortran Support
11615687f061SJacob Faibussowitsch 
11625687f061SJacob Faibussowitsch   Input Parameters:
11635687f061SJacob Faibussowitsch + comm        - The MPI communicator to abort on
11645687f061SJacob Faibussowitsch - __VA_ARGS__ - An arbitrary expression
11655687f061SJacob Faibussowitsch 
11665687f061SJacob Faibussowitsch   Level: beginner
11675687f061SJacob Faibussowitsch 
11685687f061SJacob Faibussowitsch   Notes:
11695687f061SJacob Faibussowitsch   This macro may be used to check C++ expressions for exceptions in cases where you cannot
11705687f061SJacob Faibussowitsch   return an error code. This includes constructors, destructors, copy/move assignment functions
11715687f061SJacob Faibussowitsch   or constructors among others.
11725687f061SJacob Faibussowitsch 
11735687f061SJacob Faibussowitsch   If an exception is caught, the macro calls `SETERRABORT()` on `comm`. The exception must
11745687f061SJacob Faibussowitsch   derive from `std::exception` in order to be caught.
11755687f061SJacob Faibussowitsch 
11765687f061SJacob Faibussowitsch   If the routine _can_ return an error-code it is highly advised to use `PetscCallCXX()`
11775687f061SJacob Faibussowitsch   instead.
11785687f061SJacob Faibussowitsch 
11795687f061SJacob Faibussowitsch   See `PetscCallCXX()` for additional discussion.
11805687f061SJacob Faibussowitsch 
11815687f061SJacob Faibussowitsch   Example Usage:
11825687f061SJacob Faibussowitsch .vb
11835687f061SJacob Faibussowitsch   class Foo
11845687f061SJacob Faibussowitsch   {
11855687f061SJacob Faibussowitsch     std::vector<int> data_;
11865687f061SJacob Faibussowitsch 
11875687f061SJacob Faibussowitsch   public:
11885687f061SJacob Faibussowitsch     // normally std::vector::reserve() may raise an exception, but since we handle it with
11895687f061SJacob Faibussowitsch     // PetscCallCXXAbort() we may mark this routine as noexcept!
11905687f061SJacob Faibussowitsch     Foo() noexcept
11915687f061SJacob Faibussowitsch     {
11925687f061SJacob Faibussowitsch       PetscCallCXXAbort(PETSC_COMM_SELF, data_.reserve(10));
11935687f061SJacob Faibussowitsch     }
11945687f061SJacob Faibussowitsch   };
11955687f061SJacob Faibussowitsch 
11965687f061SJacob Faibussowitsch   std::vector<int> bar()
11975687f061SJacob Faibussowitsch   {
11985687f061SJacob Faibussowitsch     std::vector<int> v;
11995687f061SJacob Faibussowitsch 
12005687f061SJacob Faibussowitsch     PetscFunctionBegin;
12015687f061SJacob Faibussowitsch     // OK!
12025687f061SJacob Faibussowitsch     PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1));
12035687f061SJacob Faibussowitsch     PetscFunctionReturn(v);
12045687f061SJacob Faibussowitsch   }
12055687f061SJacob Faibussowitsch 
12065687f061SJacob Faibussowitsch   PetscErrorCode baz()
12075687f061SJacob Faibussowitsch   {
12085687f061SJacob Faibussowitsch     std::vector<int> v;
12095687f061SJacob Faibussowitsch 
12105687f061SJacob Faibussowitsch     PetscFunctionBegin;
12115687f061SJacob Faibussowitsch     // WRONG! baz() returns a PetscErrorCode, prefer PetscCallCXX() instead
12125687f061SJacob Faibussowitsch     PetscCallCXXAbort(PETSC_COMM_SELF, v.emplace_back(1));
12133ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
12145687f061SJacob Faibussowitsch   }
12155687f061SJacob Faibussowitsch .ve
12165687f061SJacob Faibussowitsch 
12175687f061SJacob Faibussowitsch .seealso: `PetscCallCXX()`, `SETERRABORT()`, `PetscCallAbort()`
12185687f061SJacob Faibussowitsch M*/
12193ba16761SJacob Faibussowitsch #define PetscCallCXXAbort(comm, ...) PetscCallCXX_Private(SETERRABORT, comm, __VA_ARGS__)
12203f520e80SJunchao Zhang 
122130de9b25SBarry Smith /*MC
12229566063dSJacob Faibussowitsch   CHKERRCXX - Checks C++ function calls and if they throw an exception, catch it and then
12239566063dSJacob Faibussowitsch   return a PETSc error code
12249566063dSJacob Faibussowitsch 
12259566063dSJacob Faibussowitsch   Synopsis:
12267233ce55SJed Brown   #include <petscsys.h>
12279566063dSJacob Faibussowitsch   void CHKERRCXX(func) noexcept;
12289566063dSJacob Faibussowitsch 
12299566063dSJacob Faibussowitsch   Not Collective
12309566063dSJacob Faibussowitsch 
12319566063dSJacob Faibussowitsch   Input Parameter:
12329566063dSJacob Faibussowitsch . func - C++ function calls
12339566063dSJacob Faibussowitsch 
12349566063dSJacob Faibussowitsch   Level: deprecated
12359566063dSJacob Faibussowitsch 
1236667f096bSBarry Smith   Note:
1237667f096bSBarry Smith   Deprecated in favor of `PetscCallCXX()`. This routine behaves identically to it.
1238667f096bSBarry Smith 
1239db781477SPatrick Sanan .seealso: `PetscCallCXX()`
12409566063dSJacob Faibussowitsch M*/
12419566063dSJacob Faibussowitsch #define CHKERRCXX(...) PetscCallCXX(__VA_ARGS__)
12429566063dSJacob Faibussowitsch 
12439566063dSJacob Faibussowitsch /*MC
124430de9b25SBarry Smith    CHKMEMQ - Checks the memory for corruption, calls error handler if any is detected
124530de9b25SBarry Smith 
124630de9b25SBarry Smith    Synopsis:
1247aaa7dc30SBarry Smith    #include <petscsys.h>
124891d3bdf4SKris Buschelman    CHKMEMQ;
124930de9b25SBarry Smith 
1250eca87e8dSBarry Smith    Not Collective
1251eca87e8dSBarry Smith 
125230de9b25SBarry Smith    Level: beginner
125330de9b25SBarry Smith 
125430de9b25SBarry Smith    Notes:
1255edf220e6SMartin Diehl    We recommend using Valgrind <https://petsc.org/release/faq/#valgrind> or Compute Sanitizer
1256edf220e6SMartin Diehl    <https://developer.nvidia.com/compute-sanitizer> on NVIDIA CUDA systems for finding memory problems. The ``CHKMEMQ``
1257*74df5e01SJunchao Zhang    macro is useful on systems that do not have `valgrind`, but is not as good as `valgrind` or `compute-sanitizer`.
12581957e957SBarry Smith 
1259667f096bSBarry Smith    Must run with the option `-malloc_debug` (`-malloc_test` in debug mode; or if `PetscMallocSetDebug()` called) to enable this option
126030de9b25SBarry Smith 
126130de9b25SBarry Smith    Once the error handler is called the calling function is then returned from with the given error code.
126230de9b25SBarry Smith 
126330de9b25SBarry Smith    By defaults prints location where memory that is corrupted was allocated.
126430de9b25SBarry Smith 
12658af9f190SBarry Smith    Use `CHKMEMA` for functions that return `void`
1266f621e05eSBarry Smith 
1267db781477SPatrick Sanan .seealso: `PetscTraceBackErrorHandler()`, `PetscPushErrorHandler()`, `PetscError()`, `SETERRQ()`, `PetscMallocValidate()`
126830de9b25SBarry Smith M*/
12696d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
1270064a246eSJacob Faibussowitsch   #define CHKMEMQ
1271064a246eSJacob Faibussowitsch   #define CHKMEMA
12726d210af2SJacob Faibussowitsch #else
12739371c9d4SSatish Balay   #define CHKMEMQ \
12749371c9d4SSatish Balay     do { \
12753ba16761SJacob Faibussowitsch       PetscErrorCode ierr_petsc_memq_ = PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__); \
12763ba16761SJacob Faibussowitsch       if (PetscUnlikely(ierr_petsc_memq_ != PETSC_SUCCESS)) return PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ierr_petsc_memq_, PETSC_ERROR_REPEAT, " "); \
127786d09637SLisandro Dalcin     } while (0)
12786d210af2SJacob Faibussowitsch   #define CHKMEMA PetscMallocValidate(__LINE__, PETSC_FUNCTION_NAME, __FILE__)
1279064a246eSJacob Faibussowitsch #endif
12809566063dSJacob Faibussowitsch 
1281668f157eSBarry Smith /*E
1282668f157eSBarry Smith   PetscErrorType - passed to the PETSc error handling routines indicating if this is the first or a later call to the error handlers
1283668f157eSBarry Smith 
1284668f157eSBarry Smith   Level: advanced
1285668f157eSBarry Smith 
1286667f096bSBarry Smith   Note:
128787497f52SBarry Smith   `PETSC_ERROR_IN_CXX` indicates the error was detected in C++ and an exception should be generated
1288d736bfebSBarry Smith 
1289667f096bSBarry Smith   Developer Note:
1290667f096bSBarry Smith   This is currently used to decide when to print the detailed information about the run in `PetscTraceBackErrorHandler()`
1291668f157eSBarry Smith 
129287497f52SBarry Smith .seealso: `PetscError()`, `SETERRQ()`
1293668f157eSBarry Smith E*/
12949371c9d4SSatish Balay typedef enum {
12959371c9d4SSatish Balay   PETSC_ERROR_INITIAL = 0,
12969371c9d4SSatish Balay   PETSC_ERROR_REPEAT  = 1,
12979371c9d4SSatish Balay   PETSC_ERROR_IN_CXX  = 2
12989371c9d4SSatish Balay } PetscErrorType;
12994b209cf6SBarry Smith 
1300eb9e708aSLisandro Dalcin #if defined(__clang_analyzer__)
1301eb9e708aSLisandro Dalcin __attribute__((analyzer_noreturn))
1302eb9e708aSLisandro Dalcin #endif
1303b859f8c0SPierre Jolivet PETSC_EXTERN PetscErrorCode PetscError(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, ...) PETSC_ATTRIBUTE_COLD PETSC_ATTRIBUTE_FORMAT(7, 8);
1304eb9e708aSLisandro Dalcin 
1305014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscErrorPrintfInitialize(void);
13063ba16761SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscErrorMessage(PetscErrorCode, const char *[], char **);
1307d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscTraceBackErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1308d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscIgnoreErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1309d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscEmacsClientErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1310d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscMPIAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1311d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAbortErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1312d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscAttachDebuggerErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1313d8e4614bSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscReturnErrorHandler(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *) PETSC_ATTRIBUTE_COLD;
1314efca3c55SSatish Balay PETSC_EXTERN PetscErrorCode PetscPushErrorHandler(PetscErrorCode (*handler)(MPI_Comm, int, const char *, const char *, PetscErrorCode, PetscErrorType, const char *, void *), void *);
1315014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopErrorHandler(void);
13168d359177SBarry Smith PETSC_EXTERN PetscErrorCode PetscSignalHandlerDefault(int, void *);
1317014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPushSignalHandler(PetscErrorCode (*)(int, void *), void *);
1318014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscPopSignalHandler(void);
131928559dc8SJed Brown PETSC_EXTERN PetscErrorCode PetscCheckPointerSetIntensity(PetscInt);
1320c2a741eeSJunchao Zhang PETSC_EXTERN void           PetscSignalSegvCheckPointerOrMpi(void);
PetscSignalSegvCheckPointer(void)1321edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 13, 0, "PetscSignalSegvCheckPointerOrMpi()", ) static inline void PetscSignalSegvCheckPointer(void)
1322d71ae5a4SJacob Faibussowitsch {
13239371c9d4SSatish Balay   PetscSignalSegvCheckPointerOrMpi();
13249371c9d4SSatish Balay }
1325329f5518SBarry Smith 
1326639ff905SBarry Smith /*MC
1327639ff905SBarry Smith     PetscErrorPrintf - Prints error messages.
1328639ff905SBarry Smith 
1329639ff905SBarry Smith    Synopsis:
1330aaa7dc30SBarry Smith     #include <petscsys.h>
1331639ff905SBarry Smith      PetscErrorCode (*PetscErrorPrintf)(const char format[], ...);
1332639ff905SBarry Smith 
13337cdbe19fSJose E. Roman     Not Collective; No Fortran Support
13347cdbe19fSJose E. Roman 
1335f899ff85SJose E. Roman     Input Parameter:
1336cd05f99aSJacob Faibussowitsch .   format - the usual `printf()` format string
1337639ff905SBarry Smith 
1338639ff905SBarry Smith    Options Database Keys:
13391957e957SBarry Smith +  -error_output_stdout - cause error messages to be printed to stdout instead of the (default) stderr
1340e1bc860dSBarry Smith -  -error_output_none   - to turn off all printing of error messages (does not change the way the error is handled.)
1341639ff905SBarry Smith 
1342cf53795eSBarry Smith    Level: developer
1343cf53795eSBarry Smith 
134495452b02SPatrick Sanan    Notes:
134595452b02SPatrick Sanan    Use
1346667f096bSBarry Smith .vb
134795bd0b28SBarry Smith      PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the error is handled) and
1348667f096bSBarry Smith      PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on or you can use your own function
1349667f096bSBarry Smith .ve
1350639ff905SBarry Smith    Use
1351667f096bSBarry Smith .vb
135287497f52SBarry Smith      `PETSC_STDERR` = FILE* obtained from a file open etc. to have stderr printed to the file.
135387497f52SBarry Smith      `PETSC_STDOUT` = FILE* obtained from a file open etc. to have stdout printed to the file.
1354667f096bSBarry Smith .ve
1355639ff905SBarry Smith    Use
1356af27ebaaSBarry Smith .vb
135787497f52SBarry Smith       `PetscPushErrorHandler()` to provide your own error handler that determines what kind of messages to print
1358af27ebaaSBarry Smith .ve
1359639ff905SBarry Smith 
1360db781477SPatrick Sanan .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscHelpPrintf()`, `PetscPrintf()`, `PetscPushErrorHandler()`, `PetscVFPrintf()`, `PetscHelpPrintf()`
1361639ff905SBarry Smith M*/
13623ca90d2dSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode (*PetscErrorPrintf)(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1363639ff905SBarry Smith 
1364cf0818bdSBarry Smith /*E
1365cf0818bdSBarry Smith      PetscFPTrap - types of floating point exceptions that may be trapped
1366cf0818bdSBarry Smith 
136787497f52SBarry Smith      Currently only `PETSC_FP_TRAP_OFF` and `PETSC_FP_TRAP_ON` are handled. All others are treated as `PETSC_FP_TRAP_ON`.
1368cf0818bdSBarry Smith 
1369cf0818bdSBarry Smith      Level: intermediate
1370cf0818bdSBarry Smith 
13717de69702SBarry Smith .seealso: `PetscSetFPTrap()`, `PetscFPTrapPush()`
1372cf0818bdSBarry Smith  E*/
13739371c9d4SSatish Balay typedef enum {
13749371c9d4SSatish Balay   PETSC_FP_TRAP_OFF      = 0,
13759371c9d4SSatish Balay   PETSC_FP_TRAP_INDIV    = 1,
13769371c9d4SSatish Balay   PETSC_FP_TRAP_FLTOPERR = 2,
13779371c9d4SSatish Balay   PETSC_FP_TRAP_FLTOVF   = 4,
13789371c9d4SSatish Balay   PETSC_FP_TRAP_FLTUND   = 8,
13799371c9d4SSatish Balay   PETSC_FP_TRAP_FLTDIV   = 16,
1380dd460d27SBarry Smith   PETSC_FP_TRAP_FLTINEX  = 32,
1381ce78bad3SBarry Smith   PETSC_FP_TRAP_ON       = 63
13829371c9d4SSatish Balay } PetscFPTrap;
1383dd460d27SBarry Smith 
1384014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSetFPTrap(PetscFPTrap);
1385014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPush(PetscFPTrap);
1386014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscFPTrapPop(void);
1387aba4c478SBarry Smith PETSC_EXTERN PetscErrorCode PetscDetermineInitialFPTrap(void);
138854a8ef01SBarry Smith 
13893a40ed3dSBarry Smith /*
13903a40ed3dSBarry Smith       Allows the code to build a stack frame as it runs
13913a40ed3dSBarry Smith */
13923a40ed3dSBarry Smith 
139399cd645aSJed Brown #define PETSCSTACKSIZE 64
13943a40ed3dSBarry Smith typedef struct {
13950e33f6ddSBarry Smith   const char *function[PETSCSTACKSIZE];
13960e33f6ddSBarry Smith   const char *file[PETSCSTACKSIZE];
1397184914b5SBarry Smith   int         line[PETSCSTACKSIZE];
1398f0b74427SPierre Jolivet   int         petscroutine[PETSCSTACKSIZE]; /* 0 external called from PETSc, 1 PETSc functions, 2 PETSc user functions */
1399184914b5SBarry Smith   int         currentsize;
1400a2f94806SJed Brown   int         hotdepth;
14014be741a6SBarry Smith   PetscBool   check; /* option to check for correct Push/Pop semantics, true for default petscstack but not other stacks */
14023a40ed3dSBarry Smith } PetscStack;
1403dfb7d7afSStefano Zampini #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
140427104ee2SJacob Faibussowitsch PETSC_EXTERN PetscStack petscstack;
140527104ee2SJacob Faibussowitsch #endif
14063a40ed3dSBarry Smith 
14075d12eec7SSatish Balay #if defined(PETSC_SERIALIZE_FUNCTIONS)
14085d12eec7SSatish Balay   #include <petsc/private/petscfptimpl.h>
14095d12eec7SSatish Balay   /*
14105d12eec7SSatish Balay    Registers the current function into the global function pointer to function name table
14115d12eec7SSatish Balay 
14125d12eec7SSatish Balay    Have to fix this to handle errors but cannot return error since used in PETSC_VIEWER_DRAW_() etc
14135d12eec7SSatish Balay */
14149371c9d4SSatish Balay   #define PetscRegister__FUNCT__() \
14159371c9d4SSatish Balay     do { \
14165d12eec7SSatish Balay       static PetscBool __chked = PETSC_FALSE; \
14175d12eec7SSatish Balay       if (!__chked) { \
14189371c9d4SSatish Balay         void *ptr; \
14193ba16761SJacob Faibussowitsch         PetscCallAbort(PETSC_COMM_SELF, PetscDLSym(NULL, PETSC_FUNCTION_NAME, &ptr)); \
14205d12eec7SSatish Balay         __chked = PETSC_TRUE; \
14219371c9d4SSatish Balay       } \
14229371c9d4SSatish Balay     } while (0)
14235d12eec7SSatish Balay #else
14245d12eec7SSatish Balay   #define PetscRegister__FUNCT__()
14255d12eec7SSatish Balay #endif
14265d12eec7SSatish Balay 
1427eae3dc7dSJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) || defined(__clang_analyzer__)
14286d210af2SJacob Faibussowitsch   #define PetscStackPushNoCheck(funct, petsc_routine, hot)
142937154d10SBarry Smith   #define PetscStackUpdateLine
1430792fecdfSBarry Smith   #define PetscStackPushExternal(funct)
1431dd460d27SBarry Smith   #define PetscStackPopNoCheck(funct)
14326d210af2SJacob Faibussowitsch   #define PetscStackClearTop
14336d210af2SJacob Faibussowitsch   #define PetscFunctionBegin
14346d210af2SJacob Faibussowitsch   #define PetscFunctionBeginUser
14356d210af2SJacob Faibussowitsch   #define PetscFunctionBeginHot
14360a57284eSJacob Faibussowitsch   #define PetscFunctionReturn(...)  return __VA_ARGS__
14376d210af2SJacob Faibussowitsch   #define PetscFunctionReturnVoid() return
14386d210af2SJacob Faibussowitsch   #define PetscStackPop
14396d210af2SJacob Faibussowitsch   #define PetscStackPush(f)
144063bfac88SBarry Smith   #define PetscStackPush_Private(stack__, file__, func__, line__, petsc_routine__, hot__) \
144163bfac88SBarry Smith     (void)file__; \
144263bfac88SBarry Smith     (void)func__; \
144363bfac88SBarry Smith     (void)line__
144463bfac88SBarry Smith   #define PetscStackPop_Private(stack__, func__) (void)func__
1445dfb7d7afSStefano Zampini #elif defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_THREADSAFETY)
1446660278c0SBarry Smith 
14479371c9d4SSatish Balay   #define PetscStackPush_Private(stack__, file__, func__, line__, petsc_routine__, hot__) \
14489371c9d4SSatish Balay     do { \
14495a96b57dSJacob Faibussowitsch       if (stack__.currentsize < PETSCSTACKSIZE) { \
14505a96b57dSJacob Faibussowitsch         stack__.function[stack__.currentsize] = func__; \
1451ef1023bdSBarry Smith         if (petsc_routine__) { \
1452ef1023bdSBarry Smith           stack__.file[stack__.currentsize] = file__; \
14535a96b57dSJacob Faibussowitsch           stack__.line[stack__.currentsize] = line__; \
1454ef1023bdSBarry Smith         } else { \
1455648bc8c4SBarry Smith           stack__.file[stack__.currentsize] = PETSC_NULLPTR; \
1456ef1023bdSBarry Smith           stack__.line[stack__.currentsize] = 0; \
1457ef1023bdSBarry Smith         } \
14585a96b57dSJacob Faibussowitsch         stack__.petscroutine[stack__.currentsize] = petsc_routine__; \
14595a96b57dSJacob Faibussowitsch       } \
14605a96b57dSJacob Faibussowitsch       ++stack__.currentsize; \
14615a96b57dSJacob Faibussowitsch       stack__.hotdepth += (hot__ || stack__.hotdepth); \
14625a96b57dSJacob Faibussowitsch     } while (0)
14635a96b57dSJacob Faibussowitsch 
14644be741a6SBarry Smith   /* uses PetscCheckAbort() because may be used in a function that does not return an error code */
14659371c9d4SSatish Balay   #define PetscStackPop_Private(stack__, func__) \
14669371c9d4SSatish Balay     do { \
14674be741a6SBarry 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__); \
14685a96b57dSJacob Faibussowitsch       if (--stack__.currentsize < PETSCSTACKSIZE) { \
14691564b4b1SPierre Jolivet         PetscCheckAbort(!stack__.check || stack__.petscroutine[stack__.currentsize] != 1 || stack__.function[stack__.currentsize] == func__, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid stack: push from %s %s:%d. Pop from %s %s:%d.\n", \
14709371c9d4SSatish Balay                         stack__.function[stack__.currentsize], stack__.file[stack__.currentsize], stack__.line[stack__.currentsize], func__, __FILE__, __LINE__); \
14715a96b57dSJacob Faibussowitsch         stack__.function[stack__.currentsize]     = PETSC_NULLPTR; \
14725a96b57dSJacob Faibussowitsch         stack__.file[stack__.currentsize]         = PETSC_NULLPTR; \
14735a96b57dSJacob Faibussowitsch         stack__.line[stack__.currentsize]         = 0; \
14745a96b57dSJacob Faibussowitsch         stack__.petscroutine[stack__.currentsize] = 0; \
14755a96b57dSJacob Faibussowitsch       } \
14765a96b57dSJacob Faibussowitsch       stack__.hotdepth = PetscMax(stack__.hotdepth - 1, 0); \
14775a96b57dSJacob Faibussowitsch     } while (0)
14785a96b57dSJacob Faibussowitsch 
1479586f9135SBarry Smith   /*MC
1480586f9135SBarry Smith    PetscStackPushNoCheck - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is
1481586f9135SBarry Smith    currently in the source code.
1482586f9135SBarry Smith 
1483586f9135SBarry Smith    Synopsis:
1484586f9135SBarry Smith    #include <petscsys.h>
1485586f9135SBarry Smith    void PetscStackPushNoCheck(char *funct,int petsc_routine,PetscBool hot);
1486586f9135SBarry Smith 
14877cdbe19fSJose E. Roman    Not Collective
14887cdbe19fSJose E. Roman 
1489586f9135SBarry Smith    Input Parameters:
1490586f9135SBarry Smith +  funct - the function name
1491586f9135SBarry Smith .  petsc_routine - 2 user function, 1 PETSc function, 0 some other function
1492586f9135SBarry Smith -  hot - indicates that the function may be called often so expensive error checking should be turned off inside the function
1493586f9135SBarry Smith 
1494586f9135SBarry Smith    Level: developer
1495586f9135SBarry Smith 
1496586f9135SBarry Smith    Notes:
1497586f9135SBarry 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
149887497f52SBarry 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
1499586f9135SBarry Smith    help debug the problem.
1500586f9135SBarry Smith 
1501ef1023bdSBarry Smith    This version does not check the memory corruption (an expensive operation), use `PetscStackPush()` to check the memory.
1502ef1023bdSBarry Smith 
1503792fecdfSBarry Smith    Use `PetscStackPushExternal()` for a function call that is about to be made to a non-PETSc or user function (such as BLAS etc).
1504ef1023bdSBarry Smith 
1505586f9135SBarry Smith    The default stack is a global variable called `petscstack`.
1506586f9135SBarry Smith 
1507586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1508ef1023bdSBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPush()`, `PetscStackPop`,
1509792fecdfSBarry Smith           `PetscStackPushExternal()`
1510586f9135SBarry Smith M*/
15119371c9d4SSatish Balay   #define PetscStackPushNoCheck(funct, petsc_routine, hot) \
15129371c9d4SSatish Balay     do { \
1513e04113cfSBarry Smith       PetscStackSAWsTakeAccess(); \
15145a96b57dSJacob Faibussowitsch       PetscStackPush_Private(petscstack, __FILE__, funct, __LINE__, petsc_routine, hot); \
1515e04113cfSBarry Smith       PetscStackSAWsGrantAccess(); \
1516441dd030SJed Brown     } while (0)
1517441dd030SJed Brown 
1518586f9135SBarry Smith   /*MC
151987497f52SBarry Smith    PetscStackUpdateLine - in a function that has a `PetscFunctionBegin` or `PetscFunctionBeginUser` updates the stack line number to the
152037154d10SBarry Smith    current line number.
152137154d10SBarry Smith 
152237154d10SBarry Smith    Synopsis:
152337154d10SBarry Smith    #include <petscsys.h>
152437154d10SBarry Smith    void PetscStackUpdateLine
152537154d10SBarry Smith 
15267cdbe19fSJose E. Roman    Not Collective
15277cdbe19fSJose E. Roman 
152837154d10SBarry Smith    Level: developer
152937154d10SBarry Smith 
153037154d10SBarry Smith    Notes:
153187497f52SBarry Smith    Using `PetscCall()` and friends automatically handles this process
153287497f52SBarry Smith 
153337154d10SBarry 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
153437154d10SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
153537154d10SBarry Smith    help debug the problem.
153637154d10SBarry Smith 
153795bd0b28SBarry Smith    The default stack is a global variable called `petscstack`.
153837154d10SBarry Smith 
153937154d10SBarry Smith    This is used by `PetscCall()` and is otherwise not like to be needed
154037154d10SBarry Smith 
154137154d10SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`, `PetscCall()`
154237154d10SBarry Smith M*/
154337154d10SBarry Smith   #define PetscStackUpdateLine \
15443ba16761SJacob Faibussowitsch     do { \
15450b36fe0aSPierre Jolivet       if (petscstack.currentsize > 0 && petscstack.currentsize < PETSCSTACKSIZE && petscstack.function[petscstack.currentsize - 1] == PETSC_FUNCTION_NAME) petscstack.line[petscstack.currentsize - 1] = __LINE__; \
15463ba16761SJacob Faibussowitsch     } while (0)
154737154d10SBarry Smith 
154837154d10SBarry Smith   /*MC
1549792fecdfSBarry Smith    PetscStackPushExternal - Pushes a new function name onto the PETSc default stack that tracks where the running program is
1550ef1023bdSBarry Smith    currently in the source code. Does not include the filename or line number since this is called by the calling routine
1551ef1023bdSBarry Smith    for non-PETSc or user functions.
1552ef1023bdSBarry Smith 
1553ef1023bdSBarry Smith    Synopsis:
1554ef1023bdSBarry Smith    #include <petscsys.h>
1555792fecdfSBarry Smith    void PetscStackPushExternal(char *funct);
1556ef1023bdSBarry Smith 
15577cdbe19fSJose E. Roman    Not Collective
15587cdbe19fSJose E. Roman 
15592fe279fdSBarry Smith    Input Parameter:
1560ef1023bdSBarry Smith .  funct - the function name
1561ef1023bdSBarry Smith 
1562ef1023bdSBarry Smith    Level: developer
1563ef1023bdSBarry Smith 
1564ef1023bdSBarry Smith    Notes:
156587497f52SBarry Smith    Using `PetscCallExternal()` and friends automatically handles this process
156687497f52SBarry Smith 
1567ef1023bdSBarry 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
1568ef1023bdSBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1569ef1023bdSBarry Smith    help debug the problem.
1570ef1023bdSBarry Smith 
1571ef1023bdSBarry Smith    The default stack is a global variable called `petscstack`.
1572ef1023bdSBarry Smith 
1573ef1023bdSBarry Smith    This is to be used when calling an external package function such as a BLAS function.
1574ef1023bdSBarry Smith 
1575ef1023bdSBarry Smith    This also updates the stack line number for the current stack function.
1576ef1023bdSBarry Smith 
1577ef1023bdSBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1578ef1023bdSBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1579ef1023bdSBarry Smith M*/
15809371c9d4SSatish Balay   #define PetscStackPushExternal(funct) \
15819371c9d4SSatish Balay     do { \
15829371c9d4SSatish Balay       PetscStackUpdateLine; \
15839371c9d4SSatish Balay       PetscStackPushNoCheck(funct, 0, PETSC_TRUE); \
1584a8f51744SPierre Jolivet     } while (0)
1585ef1023bdSBarry Smith 
1586ef1023bdSBarry Smith   /*MC
1587586f9135SBarry Smith    PetscStackPopNoCheck - Pops a function name from the PETSc default stack that tracks where the running program is
1588586f9135SBarry Smith    currently in the source code.
1589586f9135SBarry Smith 
1590586f9135SBarry Smith    Synopsis:
1591586f9135SBarry Smith    #include <petscsys.h>
1592586f9135SBarry Smith    void PetscStackPopNoCheck(char *funct);
1593586f9135SBarry Smith 
15947cdbe19fSJose E. Roman    Not Collective
15957cdbe19fSJose E. Roman 
1596586f9135SBarry Smith    Input Parameter:
1597586f9135SBarry Smith .   funct - the function name
1598586f9135SBarry Smith 
1599586f9135SBarry Smith    Level: developer
1600586f9135SBarry Smith 
1601586f9135SBarry Smith    Notes:
160287497f52SBarry Smith    Using `PetscCall()`, `PetscCallExternal()`, `PetscCallBack()` and friends negates the need to call this
160387497f52SBarry Smith 
1604586f9135SBarry 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
1605586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1606586f9135SBarry Smith    help debug the problem.
1607586f9135SBarry Smith 
160895bd0b28SBarry Smith    The default stack is a global variable called `petscstack`.
1609586f9135SBarry Smith 
1610586f9135SBarry Smith    Developer Note:
1611586f9135SBarry Smith    `PetscStackPopNoCheck()` takes a function argument while  `PetscStackPop` does not, this difference is likely just historical.
1612586f9135SBarry Smith 
1613586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1614586f9135SBarry Smith M*/
16159371c9d4SSatish Balay   #define PetscStackPopNoCheck(funct) \
16169371c9d4SSatish Balay     do { \
1617362febeeSStefano Zampini       PetscStackSAWsTakeAccess(); \
16185a96b57dSJacob Faibussowitsch       PetscStackPop_Private(petscstack, funct); \
1619362febeeSStefano Zampini       PetscStackSAWsGrantAccess(); \
1620362febeeSStefano Zampini     } while (0)
1621362febeeSStefano Zampini 
16229371c9d4SSatish Balay   #define PetscStackClearTop \
16239371c9d4SSatish Balay     do { \
1624e04113cfSBarry Smith       PetscStackSAWsTakeAccess(); \
16259371c9d4SSatish Balay       if (petscstack.currentsize > 0 && --petscstack.currentsize < PETSCSTACKSIZE) { \
162627104ee2SJacob Faibussowitsch         petscstack.function[petscstack.currentsize]     = PETSC_NULLPTR; \
162727104ee2SJacob Faibussowitsch         petscstack.file[petscstack.currentsize]         = PETSC_NULLPTR; \
162827104ee2SJacob Faibussowitsch         petscstack.line[petscstack.currentsize]         = 0; \
162927104ee2SJacob Faibussowitsch         petscstack.petscroutine[petscstack.currentsize] = 0; \
1630441dd030SJed Brown       } \
163127104ee2SJacob Faibussowitsch       petscstack.hotdepth = PetscMax(petscstack.hotdepth - 1, 0); \
1632e04113cfSBarry Smith       PetscStackSAWsGrantAccess(); \
1633441dd030SJed Brown     } while (0)
1634441dd030SJed Brown 
163530de9b25SBarry Smith   /*MC
16361957e957SBarry Smith    PetscFunctionBegin - First executable line of each PETSc function,  used for error handling. Final
163787497f52SBarry Smith    line of PETSc functions should be `PetscFunctionReturn`(0);
163830de9b25SBarry Smith 
163930de9b25SBarry Smith    Synopsis:
1640aaa7dc30SBarry Smith    #include <petscsys.h>
164130de9b25SBarry Smith    void PetscFunctionBegin;
164230de9b25SBarry Smith 
164320f4b53cSBarry Smith    Not Collective; No Fortran Support
1644eca87e8dSBarry Smith 
164530de9b25SBarry Smith    Usage:
164630de9b25SBarry Smith .vb
164730de9b25SBarry Smith      int something;
164830de9b25SBarry Smith 
164930de9b25SBarry Smith      PetscFunctionBegin;
165030de9b25SBarry Smith .ve
165130de9b25SBarry Smith 
165230de9b25SBarry Smith    Level: developer
165330de9b25SBarry Smith 
165420f4b53cSBarry Smith    Note:
165520f4b53cSBarry Smith    Use `PetscFunctionBeginUser` for application codes.
165620f4b53cSBarry Smith 
1657586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`
165830de9b25SBarry Smith 
165930de9b25SBarry Smith M*/
16609371c9d4SSatish Balay   #define PetscFunctionBegin \
16619371c9d4SSatish Balay     do { \
1662362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_FALSE); \
1663a2f94806SJed Brown       PetscRegister__FUNCT__(); \
1664a2f94806SJed Brown     } while (0)
1665a2f94806SJed Brown 
1666a2f94806SJed Brown   /*MC
166787497f52SBarry Smith    PetscFunctionBeginHot - Substitute for `PetscFunctionBegin` to be used in functions that are called in
1668a2f94806SJed Brown    performance-critical circumstances.  Use of this function allows for lighter profiling by default.
1669a2f94806SJed Brown 
1670a2f94806SJed Brown    Synopsis:
1671aaa7dc30SBarry Smith    #include <petscsys.h>
1672a2f94806SJed Brown    void PetscFunctionBeginHot;
1673a2f94806SJed Brown 
167420f4b53cSBarry Smith    Not Collective; No Fortran Support
1675a2f94806SJed Brown 
1676a2f94806SJed Brown    Usage:
1677a2f94806SJed Brown .vb
1678a2f94806SJed Brown      int something;
1679a2f94806SJed Brown 
1680a2f94806SJed Brown      PetscFunctionBeginHot;
1681a2f94806SJed Brown .ve
1682a2f94806SJed Brown 
1683a2f94806SJed Brown    Level: developer
1684a2f94806SJed Brown 
1685586f9135SBarry Smith .seealso: `PetscFunctionBegin`, `PetscFunctionReturn()`, `PetscStackPushNoCheck()`
1686a2f94806SJed Brown 
1687a2f94806SJed Brown M*/
16889371c9d4SSatish Balay   #define PetscFunctionBeginHot \
16899371c9d4SSatish Balay     do { \
1690362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 1, PETSC_TRUE); \
16912d53ad75SBarry Smith       PetscRegister__FUNCT__(); \
169253c77d0aSJed Brown     } while (0)
169353c77d0aSJed Brown 
1694a8d2bbe5SBarry Smith   /*MC
1695530d85c1SBarry Smith    PetscFunctionBeginUser - First executable line of user provided routines
1696a8d2bbe5SBarry Smith 
1697a8d2bbe5SBarry Smith    Synopsis:
1698aaa7dc30SBarry Smith    #include <petscsys.h>
1699a8d2bbe5SBarry Smith    void PetscFunctionBeginUser;
1700a8d2bbe5SBarry Smith 
170120f4b53cSBarry Smith    Not Collective; No Fortran Support
1702a8d2bbe5SBarry Smith 
1703a8d2bbe5SBarry Smith    Usage:
1704a8d2bbe5SBarry Smith .vb
1705a8d2bbe5SBarry Smith      int something;
1706a8d2bbe5SBarry Smith 
1707ac285190SBarry Smith      PetscFunctionBeginUser;
1708a8d2bbe5SBarry Smith .ve
1709a8d2bbe5SBarry Smith 
171020f4b53cSBarry Smith    Level: intermediate
171120f4b53cSBarry Smith 
1712a8d2bbe5SBarry Smith    Notes:
1713530d85c1SBarry Smith    Functions that incorporate this must call `PetscFunctionReturn()` instead of return except for main().
1714530d85c1SBarry Smith 
1715530d85c1SBarry Smith    May be used before `PetscInitialize()`
17161957e957SBarry Smith 
1717530d85c1SBarry Smith    This is identical to `PetscFunctionBegin` except it labels the routine as a user
1718ac285190SBarry Smith    routine instead of as a PETSc library routine.
1719ac285190SBarry Smith 
1720586f9135SBarry Smith .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, `PetscFunctionBeginHot`, `PetscStackPushNoCheck()`
1721a8d2bbe5SBarry Smith M*/
17229371c9d4SSatish Balay   #define PetscFunctionBeginUser \
17239371c9d4SSatish Balay     do { \
1724362febeeSStefano Zampini       PetscStackPushNoCheck(PETSC_FUNCTION_NAME, 2, PETSC_FALSE); \
1725a8d2bbe5SBarry Smith       PetscRegister__FUNCT__(); \
1726a8d2bbe5SBarry Smith     } while (0)
1727a8d2bbe5SBarry Smith 
1728586f9135SBarry Smith   /*MC
1729586f9135SBarry Smith    PetscStackPush - Pushes a new function name and line number onto the PETSc default stack that tracks where the running program is
1730586f9135SBarry Smith    currently in the source code and verifies the memory is not corrupted.
1731586f9135SBarry Smith 
1732586f9135SBarry Smith    Synopsis:
1733586f9135SBarry Smith    #include <petscsys.h>
1734586f9135SBarry Smith    void PetscStackPush(char *funct)
1735586f9135SBarry Smith 
17367cdbe19fSJose E. Roman    Not Collective
17377cdbe19fSJose E. Roman 
1738586f9135SBarry Smith    Input Parameter:
1739586f9135SBarry Smith .  funct - the function name
1740586f9135SBarry Smith 
1741586f9135SBarry Smith    Level: developer
1742586f9135SBarry Smith 
1743586f9135SBarry Smith    Notes:
1744586f9135SBarry 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
1745586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1746586f9135SBarry Smith    help debug the problem.
1747586f9135SBarry Smith 
174895bd0b28SBarry Smith    The default stack is a global variable called `petscstack`.
1749586f9135SBarry Smith 
1750586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPopNoCheck()`, `PetscCall()`, `PetscFunctionBegin()`,
1751586f9135SBarry Smith           `PetscFunctionReturn()`, `PetscFunctionBeginHot()`, `PetscFunctionBeginUser()`, `PetscStackPushNoCheck()`, `PetscStackPop`
1752586f9135SBarry Smith M*/
17539371c9d4SSatish Balay   #define PetscStackPush(n) \
17549371c9d4SSatish Balay     do { \
1755362febeeSStefano Zampini       PetscStackPushNoCheck(n, 0, PETSC_FALSE); \
175615681b3cSBarry Smith       CHKMEMQ; \
175715681b3cSBarry Smith     } while (0)
17583a40ed3dSBarry Smith 
1759586f9135SBarry Smith   /*MC
1760586f9135SBarry Smith    PetscStackPop - Pops a function name from the PETSc default stack that tracks where the running program is
1761586f9135SBarry Smith    currently in the source code and verifies the memory is not corrupted.
1762586f9135SBarry Smith 
1763586f9135SBarry Smith    Synopsis:
1764586f9135SBarry Smith    #include <petscsys.h>
1765586f9135SBarry Smith    void PetscStackPop
1766586f9135SBarry Smith 
17677cdbe19fSJose E. Roman    Not Collective
17687cdbe19fSJose E. Roman 
1769586f9135SBarry Smith    Level: developer
1770586f9135SBarry Smith 
1771586f9135SBarry Smith    Notes:
1772586f9135SBarry 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
1773586f9135SBarry Smith    occurred, for example, when a signal is received. It is recommended to use the debugger if extensive information is needed to
1774586f9135SBarry Smith    help debug the problem.
1775586f9135SBarry Smith 
177695bd0b28SBarry Smith    The default stack is a global variable called `petscstack`.
1777586f9135SBarry Smith 
1778586f9135SBarry Smith .seealso: `PetscAttachDebugger()`, `PetscStackCopy()`, `PetscStackView()`, `PetscStackPushNoCheck()`, `PetscStackPopNoCheck()`, `PetscStackPush()`
1779586f9135SBarry Smith M*/
17809371c9d4SSatish Balay   #define PetscStackPop \
17819371c9d4SSatish Balay     do { \
1782441dd030SJed Brown       CHKMEMQ; \
1783362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
178415681b3cSBarry Smith     } while (0)
1785d64ed03dSBarry Smith 
178630de9b25SBarry Smith   /*MC
17870a57284eSJacob Faibussowitsch    PetscFunctionReturn - Last executable line of each PETSc function used for error
17880a57284eSJacob Faibussowitsch    handling. Replaces `return()`.
178930de9b25SBarry Smith 
179030de9b25SBarry Smith    Synopsis:
17917233ce55SJed Brown    #include <petscsys.h>
17920a57284eSJacob Faibussowitsch    void PetscFunctionReturn(...)
179330de9b25SBarry Smith 
1794667f096bSBarry Smith    Not Collective; No Fortran Support
1795eca87e8dSBarry Smith 
17965687f061SJacob Faibussowitsch    Level: beginner
17975687f061SJacob Faibussowitsch 
17980a57284eSJacob Faibussowitsch    Notes:
17990a57284eSJacob Faibussowitsch    This routine is a macro, so while it does not "return" anything itself, it does return from
18000a57284eSJacob Faibussowitsch    the function in the literal sense.
18010a57284eSJacob Faibussowitsch 
18020a57284eSJacob Faibussowitsch    Usually the return value is the integer literal `0` (for example in any function returning
18030a57284eSJacob Faibussowitsch    `PetscErrorCode`), however it is possible to return any arbitrary type. The arguments of
18040a57284eSJacob Faibussowitsch    this macro are placed before the `return` statement as-is.
18050a57284eSJacob Faibussowitsch 
18060a57284eSJacob Faibussowitsch    Any routine which returns via `PetscFunctionReturn()` must begin with a corresponding
18070a57284eSJacob Faibussowitsch    `PetscFunctionBegin`.
18080a57284eSJacob Faibussowitsch 
18090a57284eSJacob Faibussowitsch    For routines which return `void` use `PetscFunctionReturnVoid()` instead.
18100a57284eSJacob Faibussowitsch 
18110a57284eSJacob Faibussowitsch    Example Usage:
181230de9b25SBarry Smith .vb
18130a57284eSJacob Faibussowitsch    PetscErrorCode foo(int *x)
18140a57284eSJacob Faibussowitsch    {
18150a57284eSJacob Faibussowitsch      PetscFunctionBegin; // don't forget the begin!
18160a57284eSJacob Faibussowitsch      *x = 10;
18173ba16761SJacob Faibussowitsch      PetscFunctionReturn(PETSC_SUCCESS);
181830de9b25SBarry Smith    }
181930de9b25SBarry Smith .ve
182030de9b25SBarry Smith 
18210a57284eSJacob Faibussowitsch    May return any arbitrary type\:
18220a57284eSJacob Faibussowitsch .vb
18230a57284eSJacob Faibussowitsch   struct Foo
18240a57284eSJacob Faibussowitsch   {
18250a57284eSJacob Faibussowitsch     int x;
18260a57284eSJacob Faibussowitsch   };
18270a57284eSJacob Faibussowitsch 
18280a57284eSJacob Faibussowitsch   struct Foo make_foo(int value)
18290a57284eSJacob Faibussowitsch   {
18300a57284eSJacob Faibussowitsch     struct Foo f;
18310a57284eSJacob Faibussowitsch 
18320a57284eSJacob Faibussowitsch     PetscFunctionBegin;
18330a57284eSJacob Faibussowitsch     f.x = value;
18340a57284eSJacob Faibussowitsch     PetscFunctionReturn(f);
18350a57284eSJacob Faibussowitsch   }
18360a57284eSJacob Faibussowitsch .ve
18370a57284eSJacob Faibussowitsch 
18380a57284eSJacob Faibussowitsch .seealso: `PetscFunctionBegin`, `PetscFunctionBeginUser`, `PetscFunctionReturnVoid()`,
18390a57284eSJacob Faibussowitsch           `PetscStackPopNoCheck()`
184030de9b25SBarry Smith M*/
18415687f061SJacob Faibussowitsch   #define PetscFunctionReturn(...) \
18429371c9d4SSatish Balay     do { \
1843362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
18445687f061SJacob Faibussowitsch       return __VA_ARGS__; \
184527104ee2SJacob Faibussowitsch     } while (0)
1846d64ed03dSBarry Smith 
18470a57284eSJacob Faibussowitsch   /*MC
18480a57284eSJacob Faibussowitsch   PetscFunctionReturnVoid - Like `PetscFunctionReturn()` but returns `void`
18490a57284eSJacob Faibussowitsch 
18500a57284eSJacob Faibussowitsch   Synopsis:
18517233ce55SJed Brown   #include <petscsys.h>
18520a57284eSJacob Faibussowitsch   void PetscFunctionReturnVoid()
18530a57284eSJacob Faibussowitsch 
18540a57284eSJacob Faibussowitsch   Not Collective
18550a57284eSJacob Faibussowitsch 
18565687f061SJacob Faibussowitsch   Level: beginner
18575687f061SJacob Faibussowitsch 
18585687f061SJacob Faibussowitsch   Note:
18590a57284eSJacob Faibussowitsch   Behaves identically to `PetscFunctionReturn()` except that it returns `void`. That is, this
18605687f061SJacob Faibussowitsch   macro culminates with `return`.
18610a57284eSJacob Faibussowitsch 
18620a57284eSJacob Faibussowitsch   Example Usage:
18630a57284eSJacob Faibussowitsch .vb
18640a57284eSJacob Faibussowitsch   void foo()
18650a57284eSJacob Faibussowitsch   {
18660a57284eSJacob Faibussowitsch     PetscFunctionBegin; // must start with PetscFunctionBegin!
18670a57284eSJacob Faibussowitsch     bar();
18680a57284eSJacob Faibussowitsch     baz();
18690a57284eSJacob Faibussowitsch     PetscFunctionReturnVoid();
18700a57284eSJacob Faibussowitsch   }
18710a57284eSJacob Faibussowitsch .ve
18720a57284eSJacob Faibussowitsch 
18730a57284eSJacob Faibussowitsch .seealso: `PetscFunctionReturn()`, `PetscFunctionBegin`, PetscFunctionBeginUser`
18740a57284eSJacob Faibussowitsch M*/
18759371c9d4SSatish Balay   #define PetscFunctionReturnVoid() \
18769371c9d4SSatish Balay     do { \
1877362febeeSStefano Zampini       PetscStackPopNoCheck(PETSC_FUNCTION_NAME); \
187827104ee2SJacob Faibussowitsch       return; \
187927104ee2SJacob Faibussowitsch     } while (0)
188027104ee2SJacob Faibussowitsch #else /* PETSC_USE_DEBUG */
188127104ee2SJacob Faibussowitsch   #define PetscStackPushNoCheck(funct, petsc_routine, hot)
188237154d10SBarry Smith   #define PetscStackUpdateLine
1883792fecdfSBarry Smith   #define PetscStackPushExternal(funct)
18843ba16761SJacob Faibussowitsch   #define PetscStackPopNoCheck(...)
188527104ee2SJacob Faibussowitsch   #define PetscStackClearTop
18863a40ed3dSBarry Smith   #define PetscFunctionBegin
18870bdf7c52SPeter Brune   #define PetscFunctionBeginUser
1888a2f94806SJed Brown   #define PetscFunctionBeginHot
18890a57284eSJacob Faibussowitsch   #define PetscFunctionReturn(...)  return __VA_ARGS__
18905665465eSBarry Smith   #define PetscFunctionReturnVoid() return
1891812af9f3SBarry Smith   #define PetscStackPop             CHKMEMQ
1892812af9f3SBarry Smith   #define PetscStackPush(f)         CHKMEMQ
189327104ee2SJacob Faibussowitsch #endif /* PETSC_USE_DEBUG */
18943a40ed3dSBarry Smith 
18956d210af2SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
18963ba16761SJacob Faibussowitsch   #define PetscStackCallExternalVoid(...)
18973ba16761SJacob Faibussowitsch template <typename F, typename... Args>
18983ba16761SJacob Faibussowitsch void PetscCallExternal(F, Args...);
18991fc6d13cSAlexander template <typename F, typename... Args>
19001fc6d13cSAlexander void PetscCallExternalAbort(F, Args...);
19016d210af2SJacob Faibussowitsch #else
1902586f9135SBarry Smith   /*MC
1903e77caa6dSBarry Smith     PetscStackCallExternalVoid - Calls an external library routine or user function after pushing the name of the routine on the stack.
1904eb6b5d47SBarry Smith 
1905eb6b5d47SBarry Smith    Input Parameters:
1906eb6b5d47SBarry Smith +   name    - string that gives the name of the function being called
1907586f9135SBarry Smith -   routine - actual call to the routine, for example, functionname(a,b)
1908fd3f9acdSBarry Smith 
1909586f9135SBarry Smith    Level: developer
1910eb6b5d47SBarry Smith 
191195bd0b28SBarry Smith    Notes:
1912792fecdfSBarry Smith    Often one should use `PetscCallExternal()` instead. This routine is intended for external library routines that DO NOT return error codes
1913eb6b5d47SBarry Smith 
1914586f9135SBarry Smith    In debug mode this also checks the memory for corruption at the end of the function call.
1915586f9135SBarry Smith 
191695bd0b28SBarry Smith    Certain external packages, such as BLAS/LAPACK may have their own macros, `PetscCallBLAS()` for managing the call, error checking, etc.
1917586f9135SBarry Smith 
1918586f9135SBarry Smith    Developer Note:
1919586f9135SBarry 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.
1920586f9135SBarry Smith 
1921792fecdfSBarry Smith .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscCallExternal()`, `PetscCallBLAS()`
1922586f9135SBarry Smith @*/
19233ba16761SJacob Faibussowitsch   #define PetscStackCallExternalVoid(name, ...) \
19249371c9d4SSatish Balay     do { \
192528770333SStefano Zampini       PetscStackPushExternal(name); \
19263ba16761SJacob Faibussowitsch       __VA_ARGS__; \
19279371c9d4SSatish Balay       PetscStackPop; \
19289371c9d4SSatish Balay     } while (0)
1929eb6b5d47SBarry Smith 
1930586f9135SBarry Smith   /*MC
1931792fecdfSBarry Smith     PetscCallExternal - Calls an external library routine that returns an error code after pushing the name of the routine on the stack.
1932fd3f9acdSBarry Smith 
1933fd3f9acdSBarry Smith    Input Parameters:
1934fd3f9acdSBarry Smith +  func - name of the routine
1935586f9135SBarry Smith -  args - arguments to the routine
1936586f9135SBarry Smith 
1937586f9135SBarry Smith    Level: developer
1938fd3f9acdSBarry Smith 
193995452b02SPatrick Sanan    Notes:
1940e77caa6dSBarry Smith    This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not.
1941dbf62e16SBarry Smith 
1942586f9135SBarry Smith    In debug mode this also checks the memory for corruption at the end of the function call.
1943fd3f9acdSBarry Smith 
194487497f52SBarry Smith    Assumes the error return code of the function is an integer and that a value of 0 indicates success
194587497f52SBarry Smith 
1946586f9135SBarry Smith    Developer Note:
1947d5b43468SJose E. Roman    This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc.
1948586f9135SBarry Smith 
19491fc6d13cSAlexander .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()`, `PetscCallExternalAbort()`
1950586f9135SBarry Smith M*/
19519371c9d4SSatish Balay   #define PetscCallExternal(func, ...) \
19529371c9d4SSatish Balay     do { \
1953a74df02fSJacob Faibussowitsch       PetscStackPush(PetscStringize(func)); \
19549ad2cedaSPierre Jolivet       int ierr_petsc_call_external_ = (int)func(__VA_ARGS__); \
19551d4906efSStefano Zampini       PetscStackPop; \
19563ba16761SJacob 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_); \
1957fd3f9acdSBarry Smith     } while (0)
19581fc6d13cSAlexander 
19591fc6d13cSAlexander   /*MC
19601fc6d13cSAlexander     PetscCallExternalAbort - Calls an external library routine that returns an error code after pushing the name of the routine on the stack. If the external library function return code indicates an error, this prints the error and aborts
19611fc6d13cSAlexander 
19621fc6d13cSAlexander    Input Parameters:
19631fc6d13cSAlexander +  func - name of the routine
19641fc6d13cSAlexander -  args - arguments to the routine
19651fc6d13cSAlexander 
19661fc6d13cSAlexander    Level: developer
19671fc6d13cSAlexander 
19681fc6d13cSAlexander    Notes:
19691fc6d13cSAlexander    This is intended for external package routines that return error codes. Use `PetscStackCallExternalVoid()` for those that do not.
19701fc6d13cSAlexander 
19711fc6d13cSAlexander    In debug mode this also checks the memory for corruption at the end of the function call.
19721fc6d13cSAlexander 
19731fc6d13cSAlexander    Assumes the error return code of the function is an integer and that a value of 0 indicates success
19741fc6d13cSAlexander 
19751fc6d13cSAlexander    Developer Note:
19761fc6d13cSAlexander    This is so that when an external package routine results in a crash or corrupts memory, they get blamed instead of PETSc.
19771fc6d13cSAlexander 
19781fc6d13cSAlexander .seealso: `PetscCall()`, `PetscStackPushNoCheck()`, `PetscStackPush()`, `PetscStackCallExternalVoid()`, `PetscCallExternal()`
19791fc6d13cSAlexander M*/
19801fc6d13cSAlexander   #define PetscCallExternalAbort(func, ...) \
19811fc6d13cSAlexander     do { \
19821fc6d13cSAlexander       PetscStackUpdateLine; \
19831fc6d13cSAlexander       int ierr_petsc_call_external_ = func(__VA_ARGS__); \
19841fc6d13cSAlexander       if (PetscUnlikely(ierr_petsc_call_external_ != 0)) { \
19851fc6d13cSAlexander         (void)PetscError(PETSC_COMM_SELF, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_LIB, PETSC_ERROR_INITIAL, "Error in %s(): error code %d", PetscStringize(func), ierr_petsc_call_external_); \
19861fc6d13cSAlexander         PETSCABORTWITHIERR_Private(PETSC_COMM_SELF, PETSC_ERR_LIB); \
19871fc6d13cSAlexander       } \
19881fc6d13cSAlexander     } while (0)
19896d210af2SJacob Faibussowitsch #endif /* PETSC_CLANG_STATIC_ANALYZER */
1990