xref: /petsc/include/petscsys.h (revision 2ff79c18c26c94ed8cb599682f680f231dca6444)
1 /*
2    This is the main PETSc include file (for C and C++).  It is included by all
3    other PETSc include files, so it almost never has to be specifically included.
4    Portions of this code are under:
5    Copyright (c) 2022 Advanced Micro Devices, Inc. All rights reserved.
6 */
7 #pragma once
8 
9 /*MC
10    PeOP - indicates an argument to a PETSc function is optional and one can pass `NULL` instead. This is used by the Fortran API generator
11 
12    Level: developer
13 
14    Example:
15 .vb
16    PetscErrorCode XXXX(Vec v, PeOp PetscObject obj, PeOp PetscInt *idx, PeOp PetscInt *array[])
17 .ve
18 
19    Notes:
20    This is not part of the PETSc public API and should only be used in PETSc source code.
21 
22    Put this in the function declaration in front of each variable that is optional
23 
24    Developer Note:
25    Shortened form of PETSc optional
26 
27 .seealso: `PeNS`, `PeNSS`, `PeCtx`, `PetscInitialize()`
28 M*/
29 #define PeOp
30 
31 /*MC
32    PeNS - indicates a function that does not use the PETSc standard arguments which make it easy to generate automatic language stubs for other languages
33 
34    Level: developer
35 
36    Notes:
37    This is not part of the PETSc public API and should only be used in PETSc source code.
38 
39    Put this at the end of the function declaration closing parenthesis
40 
41    Developer Note:
42    Shortened form of PETSc non-standard
43 
44 .seealso: `PeOp`, `PeNSS`, `PeCtx`, `PetscInitialize()`
45 M*/
46 #define PeNS
47 
48 /*MC
49    PeNSS - indicates a function that needs a special treatment in the C-side stub when generating the binding for other languages
50 
51    Level: developer
52 
53    Notes:
54    This is not part of the PETSc public API and should only be used in PETSc source code.
55 
56    Put this at the end of the function declaration closing parenthesis
57 
58    It is similar to PeNS; in Fortran it will generate the Fortran interface definition automatically but not the C stub, which should be added manually under the appropriate `ftn-custom` directory
59 
60    Developer Note:
61    Shortened form of PETSc non-standard stub
62 
63 .seealso: `PeOp`, `PeNS`, `PeCtx`, `PetscInitialize()`
64 M*/
65 #define PeNSS
66 
67 /* ========================================================================== */
68 /*
69    petscconf.h is contained in ${PETSC_ARCH}/include/petscconf.h it is
70    found automatically by the compiler due to the -I${PETSC_DIR}/${PETSC_ARCH}/include that
71    PETSc's makefiles add to the compiler rules.
72    For --prefix installs the directory ${PETSC_ARCH} does not exist and petscconf.h is in the same
73    directory as the other PETSc include files.
74 */
75 #include <petscconf.h>
76 #include <petscpkg_version.h>
77 #include <petscconf_poison.h>
78 #include <petscfix.h>
79 #include <petscmacros.h>
80 
81 /* SUBMANSEC = Sys */
82 
83 #if defined(PETSC_DESIRE_FEATURE_TEST_MACROS)
84   /*
85    Feature test macros must be included before headers defined by IEEE Std 1003.1-2001
86    We only turn these in PETSc source files that require them by setting PETSC_DESIRE_FEATURE_TEST_MACROS
87 */
88   #if defined(PETSC__POSIX_C_SOURCE_200112L) && !defined(_POSIX_C_SOURCE)
89     #define _POSIX_C_SOURCE 200112L
90   #endif
91   #if defined(PETSC__BSD_SOURCE) && !defined(_BSD_SOURCE)
92     #define _BSD_SOURCE
93   #endif
94   #if defined(PETSC__DEFAULT_SOURCE) && !defined(_DEFAULT_SOURCE)
95     #define _DEFAULT_SOURCE
96   #endif
97   #if defined(PETSC__GNU_SOURCE) && !defined(_GNU_SOURCE)
98     #define _GNU_SOURCE
99   #endif
100 #endif
101 
102 #include <petscsystypes.h>
103 
104 /* ========================================================================== */
105 
106 /*
107     Defines the interface to MPI allowing the use of all MPI functions.
108 
109     PETSc does not use the C++ binding of MPI at ALL. The following flag
110     makes sure the C++ bindings are not included. The C++ bindings REQUIRE
111     putting mpi.h before ANY C++ include files, we cannot control this
112     with all PETSc users. Users who want to use the MPI C++ bindings can include
113     mpicxx.h directly in their code
114 */
115 #if !defined(MPICH_SKIP_MPICXX)
116   #define MPICH_SKIP_MPICXX 1
117 #endif
118 #if !defined(OMPI_SKIP_MPICXX)
119   #define OMPI_SKIP_MPICXX 1
120 #endif
121 #if defined(PETSC_HAVE_MPIUNI)
122   #include <petsc/mpiuni/mpi.h>
123 #else
124   #include <mpi.h>
125 #endif
126 
127 /*
128    Perform various sanity checks that the correct mpi.h is being included at compile time.
129    This usually happens because
130       * either an unexpected mpi.h is in the default compiler path (i.e. in /usr/include) or
131       * an extra include path -I/something (which contains the unexpected mpi.h) is being passed to the compiler
132    Note: with MPICH and OpenMPI, accept versions [x.y.z, x+1.0.0) as compatible
133 */
134 #if defined(PETSC_HAVE_MPIUNI)
135   #ifndef MPIUNI_H
136     #error "PETSc was configured with --with-mpi=0 but now appears to be compiling using a different mpi.h"
137   #endif
138 #elif defined(PETSC_HAVE_I_MPI)
139   #if !defined(I_MPI_NUMVERSION)
140     #error "PETSc was configured with I_MPI but now appears to be compiling using a non-I_MPI mpi.h"
141   #elif I_MPI_NUMVERSION != PETSC_PKG_I_MPI_NUMVERSION
142     #error "PETSc was configured with one I_MPI mpi.h version but now appears to be compiling using a different I_MPI mpi.h version"
143   #endif
144 #elif defined(PETSC_HAVE_MVAPICH2)
145   #if !defined(MVAPICH2_NUMVERSION)
146     #error "PETSc was configured with MVAPICH2 but now appears to be compiling using a non-MVAPICH2 mpi.h"
147   #elif MVAPICH2_NUMVERSION != PETSC_PKG_MVAPICH2_NUMVERSION
148     #error "PETSc was configured with one MVAPICH2 mpi.h version but now appears to be compiling using a different MVAPICH2 mpi.h version"
149   #endif
150 #elif defined(PETSC_HAVE_MPICH)
151   #if !defined(MPICH_NUMVERSION) || defined(MVAPICH2_NUMVERSION) || defined(I_MPI_NUMVERSION)
152     #error "PETSc was configured with MPICH but now appears to be compiling using a non-MPICH mpi.h"
153   #elif PETSC_PKG_MPICH_VERSION_GT(MPICH_NUMVERSION / 10000000, MPICH_NUMVERSION / 100000 % 100, MPICH_NUMVERSION / 1000 % 100)
154     #error "PETSc was configured with one MPICH mpi.h version but now appears to be compiling using an older MPICH mpi.h version"
155   #elif PETSC_PKG_MPICH_VERSION_LT(MPICH_NUMVERSION / 10000000, 0, 0)
156     #error "PETSc was configured with one MPICH mpi.h version but now appears to be compiling using a newer major MPICH mpi.h version"
157   #endif
158 #elif defined(PETSC_HAVE_OPENMPI)
159   #if !defined(OMPI_MAJOR_VERSION)
160     #error "PETSc was configured with Open MPI but now appears to be compiling using a non-Open MPI mpi.h"
161   #elif PETSC_PKG_OPENMPI_VERSION_GT(OMPI_MAJOR_VERSION, OMPI_MINOR_VERSION, OMPI_RELEASE_VERSION)
162     #error "PETSc was configured with one Open MPI mpi.h version but now appears to be compiling using an older Open MPI mpi.h version"
163   #elif PETSC_PKG_OPENMPI_VERSION_LT(OMPI_MAJOR_VERSION, 0, 0)
164     #error "PETSc was configured with one Open MPI mpi.h version but now appears to be compiling using a newer major Open MPI mpi.h version"
165   #endif
166 #elif defined(PETSC_HAVE_MSMPI_VERSION)
167   #if !defined(MSMPI_VER)
168     #error "PETSc was configured with MSMPI but now appears to be compiling using a non-MSMPI mpi.h"
169   #elif (MSMPI_VER != PETSC_HAVE_MSMPI_VERSION)
170     #error "PETSc was configured with one MSMPI mpi.h version but now appears to be compiling using a different MSMPI mpi.h version"
171   #endif
172 #elif defined(OMPI_MAJOR_VERSION) || defined(MPICH_NUMVERSION) || defined(MSMPI_VER)
173   #error "PETSc was configured with undetermined MPI - but now appears to be compiling using any of Open MPI, MS-MPI or a MPICH variant"
174 #endif
175 
176 /*
177     Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler
178     see the top of mpicxx.h in the MPICH2 distribution.
179 */
180 #include <stdio.h>
181 
182 /* MSMPI on 32-bit Microsoft Windows requires this yukky hack - that breaks MPI standard compliance */
183 #if !defined(MPIAPI)
184   #define MPIAPI
185 #endif
186 
187 PETSC_EXTERN MPI_Datatype MPIU_ENUM PETSC_ATTRIBUTE_MPI_TYPE_TAG(PetscEnum);
188 #define MPIU_BOOL MPI_C_BOOL PETSC_DEPRECATED_MACRO(3, 24, 0, "MPI_C_BOOL", )
189 
190 /*MC
191    MPIU_INT - Portable MPI datatype corresponding to `PetscInt` independent of the precision of `PetscInt`
192 
193    Level: beginner
194 
195    Note:
196    In MPI calls that require an MPI datatype that matches a `PetscInt` or array of `PetscInt` values, pass this value.
197 
198 .seealso: `PetscReal`, `PetscScalar`, `PetscComplex`, `PetscInt`, `MPIU_COUNT`, `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`
199 M*/
200 
201 PETSC_EXTERN MPI_Datatype MPIU_FORTRANADDR;
202 
203 #if defined(PETSC_USE_64BIT_INDICES)
204   #define MPIU_INT MPIU_INT64
205 #else
206   #define MPIU_INT MPI_INT
207 #endif
208 
209 /*MC
210    MPIU_COUNT - Portable MPI datatype corresponding to `PetscCount` independent of the precision of `PetscCount`
211 
212    Level: beginner
213 
214    Note:
215    In MPI calls that require an MPI datatype that matches a `PetscCount` or array of `PetscCount` values, pass this value.
216 
217   Developer Note:
218   It seems `MPI_AINT` is unsigned so this may be the wrong choice here since `PetscCount` is signed
219 
220 .seealso: `PetscReal`, `PetscScalar`, `PetscComplex`, `PetscInt`, `MPIU_INT`, `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`
221 M*/
222 #define MPIU_COUNT MPI_AINT
223 
224 /*
225     For the rare cases when one needs to send a size_t object with MPI
226 */
227 PETSC_EXTERN MPI_Datatype MPIU_SIZE_T PETSC_ATTRIBUTE_MPI_TYPE_TAG(size_t);
228 
229 /*
230       You can use PETSC_STDOUT as a replacement of stdout. You can also change
231     the value of PETSC_STDOUT to redirect all standard output elsewhere
232 */
233 PETSC_EXTERN FILE *PETSC_STDOUT;
234 
235 /*
236       You can use PETSC_STDERR as a replacement of stderr. You can also change
237     the value of PETSC_STDERR to redirect all standard error elsewhere
238 */
239 PETSC_EXTERN FILE *PETSC_STDERR;
240 
241 /*
242   Handle inclusion when using clang compiler with CUDA support
243   __float128 is not available for the device
244 */
245 #if defined(__clang__) && (defined(__CUDA_ARCH__) || defined(__HIPCC__))
246   #define PETSC_SKIP_REAL___FLOAT128
247 #endif
248 
249 /*
250     Declare extern C stuff after including external header files
251 */
252 
253 PETSC_EXTERN PetscBool PETSC_RUNNING_ON_VALGRIND;
254 /*
255     Defines elementary mathematics functions and constants.
256 */
257 #include <petscmath.h>
258 
259 /*MC
260    PETSC_IGNORE - same as `NULL`, means PETSc will ignore this argument
261 
262    Level: beginner
263 
264    Note:
265    Accepted by many PETSc functions to not set a parameter and instead use a default value
266 
267    Fortran Note:
268    Use `PETSC_NULL_INTEGER`, `PETSC_NULL_SCALAR` etc
269 
270 .seealso: `PETSC_DECIDE`, `PETSC_DEFAULT`, `PETSC_DETERMINE`
271 M*/
272 #define PETSC_IGNORE PETSC_NULLPTR
273 #define PETSC_NULL   PETSC_DEPRECATED_MACRO(3, 19, 0, "PETSC_NULLPTR", ) PETSC_NULLPTR
274 
275 /*MC
276    PETSC_UNLIMITED - standard way of passing an integer or floating point parameter to indicate PETSc there is no bound on the value allowed
277 
278    Level: beginner
279 
280    Example Usage:
281 .vb
282    KSPSetTolerances(ksp, PETSC_CURRENT, PETSC_CURRENT, PETSC_UNLIMITED, PETSC_UNLIMITED);
283 .ve
284   indicates that the solver is allowed to take any number of iterations and will not stop early no matter how the residual gets.
285 
286    Fortran Note:
287    Use `PETSC_UNLIMITED_INTEGER` or `PETSC_UNLIMITED_REAL`.
288 
289 .seealso: `PETSC_DEFAULT`, `PETSC_IGNORE`, `PETSC_DETERMINE`, `PETSC_DECIDE`
290 M*/
291 
292 /*MC
293    PETSC_DECIDE - standard way of passing an integer or floating point parameter to indicate PETSc should determine an appropriate value
294 
295    Level: beginner
296 
297    Example Usage:
298 .vb
299    VecSetSizes(ksp, PETSC_DECIDE, 10);
300 .ve
301   indicates that the global size of the vector is 10 and the local size will be automatically determined so that the sum of the
302   local sizes is the global size, see `PetscSplitOwnership()`.
303 
304    Fortran Note:
305    Use `PETSC_DECIDE_INTEGER` or `PETSC_DECIDE_REAL`.
306 
307 .seealso: `PETSC_DEFAULT`, `PETSC_IGNORE`, `PETSC_DETERMINE`, `PETSC_UNLIMITED'
308 M*/
309 
310 /*MC
311    PETSC_DETERMINE - standard way of passing an integer or floating point parameter to indicate PETSc should determine an appropriate value
312 
313    Level: beginner
314 
315     Example Usage:
316 .vb
317    VecSetSizes(ksp, 10, PETSC_DETERMINE);
318 .ve
319   indicates that the local size of the vector is 10 and the global size will be automatically summing up all the local sizes.
320 
321    Note:
322    Same as `PETSC_DECIDE`
323 
324    Fortran Note:
325    Use `PETSC_DETERMINE_INTEGER` or `PETSC_DETERMINE_REAL`.
326 
327    Developer Note:
328    I would like to use const `PetscInt` `PETSC_DETERMINE` = `PETSC_DECIDE`; but for
329    some reason this is not allowed by the standard even though `PETSC_DECIDE` is a constant value.
330 
331 .seealso: `PETSC_DECIDE`, `PETSC_DEFAULT`, `PETSC_IGNORE`, `VecSetSizes()`, `PETSC_UNLIMITED'
332 M*/
333 
334 /*MC
335    PETSC_CURRENT - standard way of indicating to an object not to change the current value of the parameter in the object
336 
337    Level: beginner
338 
339    Note:
340    Use `PETSC_DECIDE` to use the value that was set by PETSc when the object's type was set
341 
342    Fortran Note:
343    Use `PETSC_CURRENT_INTEGER` or `PETSC_CURRENT_REAL`.
344 
345 .seealso: `PETSC_DECIDE`, `PETSC_IGNORE`, `PETSC_DETERMINE`, `PETSC_DEFAULT`, `PETSC_UNLIMITED'
346 M*/
347 
348 /*MC
349    PETSC_DEFAULT - deprecated, see `PETSC_CURRENT` and `PETSC_DETERMINE`
350 
351    Level: beginner
352 
353    Note:
354    The name is confusing since it tells the object to continue to use the value it is using, not the default value when the object's type was set.
355 
356    Developer Note:
357    Unfortunately this was used for two different purposes in the past, to actually trigger the use of a default value or to continue the
358    use of currently set value (in, for example, `KSPSetTolerances()`.
359 
360 .seealso: `PETSC_DECIDE`, `PETSC_IGNORE`, `PETSC_DETERMINE`, `PETSC_CURRENT`, `PETSC_UNLIMITED'
361 M*/
362 
363 /* These MUST be preprocessor defines! see https://gitlab.com/petsc/petsc/-/issues/1370 */
364 #define PETSC_DECIDE    (-1)
365 #define PETSC_DETERMINE PETSC_DECIDE
366 #define PETSC_CURRENT   (-2)
367 #define PETSC_UNLIMITED (-3)
368 /*  PETSC_DEFAULT is deprecated in favor of PETSC_CURRENT for use in KSPSetTolerances() and similar functions */
369 #define PETSC_DEFAULT PETSC_CURRENT
370 
371 /*MC
372    PETSC_COMM_WORLD - the equivalent of the `MPI_COMM_WORLD` communicator which represents all the processes that PETSc knows about.
373 
374    Level: beginner
375 
376    Notes:
377    By default `PETSC_COMM_WORLD` and `MPI_COMM_WORLD` are identical unless you wish to
378    run PETSc on ONLY a subset of `MPI_COMM_WORLD`. In that case create your new (smaller)
379    communicator, call it, say comm, and set `PETSC_COMM_WORLD` = comm BEFORE calling
380    `PetscInitialize()`, but after `MPI_Init()` has been called.
381 
382    The value of `PETSC_COMM_WORLD` should never be used or accessed before `PetscInitialize()`
383    is called because it may not have a valid value yet.
384 
385 .seealso: `PETSC_COMM_SELF`
386 M*/
387 PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;
388 
389 /*MC
390    PETSC_COMM_SELF - This is always `MPI_COMM_SELF`
391 
392    Level: beginner
393 
394    Note:
395    Do not USE/access or set this variable before `PetscInitialize()` has been called.
396 
397 .seealso: `PETSC_COMM_WORLD`
398 M*/
399 #define PETSC_COMM_SELF MPI_COMM_SELF
400 
401 /*MC
402    PETSC_MPI_THREAD_REQUIRED - the required threading support used if PETSc initializes MPI with `MPI_Init_thread()`.
403 
404    No Fortran Support
405 
406    Level: beginner
407 
408    Note:
409    By default `PETSC_MPI_THREAD_REQUIRED` equals `MPI_THREAD_FUNNELED` when the MPI implementation provides `MPI_Init_thread()`, otherwise it equals `MPI_THREAD_SINGLE`
410 
411 .seealso: `PetscInitialize()`
412 M*/
413 PETSC_EXTERN PetscMPIInt PETSC_MPI_THREAD_REQUIRED;
414 
415 /*MC
416    PetscBeganMPI - indicates if PETSc initialized MPI using `MPI_Init()` during `PetscInitialize()` or if MPI was already initialized with `MPI_Init()`
417 
418    Synopsis:
419    #include <petscsys.h>
420    PetscBool PetscBeganMPI;
421 
422    No Fortran Support
423 
424    Level: developer
425 
426    Note:
427    `MPI_Init()` can never be called after `PetscInitialize()`
428 
429 .seealso: `PetscInitialize()`, `PetscInitializeCalled`
430 M*/
431 PETSC_EXTERN PetscBool PetscBeganMPI;
432 
433 PETSC_EXTERN PetscBool PetscErrorHandlingInitialized;
434 PETSC_EXTERN PetscBool PetscInitializeCalled;
435 PETSC_EXTERN PetscBool PetscFinalizeCalled;
436 PETSC_EXTERN PetscBool PetscViennaCLSynchronize;
437 
438 PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm), PetscErrorCode (*)(MPI_Comm));
439 PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm, MPI_Comm *, int *);
440 PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm *);
441 PETSC_EXTERN PetscErrorCode PetscCommGetComm(MPI_Comm, MPI_Comm *);
442 PETSC_EXTERN PetscErrorCode PetscCommRestoreComm(MPI_Comm, MPI_Comm *);
443 
444 #if defined(PETSC_HAVE_KOKKOS)
445 PETSC_EXTERN PetscErrorCode PetscKokkosInitializeCheck(void); /* Initialize Kokkos if not yet. */
446 #endif
447 
448 #if defined(PETSC_HAVE_NVSHMEM)
449 PETSC_EXTERN PetscBool      PetscBeganNvshmem;
450 PETSC_EXTERN PetscBool      PetscNvshmemInitialized;
451 PETSC_EXTERN PetscErrorCode PetscNvshmemFinalize(void);
452 #endif
453 
454 #if defined(PETSC_HAVE_ELEMENTAL)
455 PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
456 PETSC_EXTERN PetscErrorCode PetscElementalInitialized(PetscBool *);
457 PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
458 #endif
459 
460 /*MC
461    PetscMalloc - Allocates memory for use with PETSc. One should use `PetscNew()`, `PetscMalloc1()` or `PetscCalloc1()` usually instead of `PetscMalloc()`
462 
463    Synopsis:
464     #include <petscsys.h>
465    PetscErrorCode PetscMalloc(size_t m,void **result)
466 
467    Not Collective
468 
469    Input Parameter:
470 .  m - number of bytes to allocate
471 
472    Output Parameter:
473 .  result - memory allocated
474 
475    Level: beginner
476 
477    Notes:
478    Memory is always allocated at least double aligned
479 
480    It is safe to allocate with an m of 0 and pass the resulting pointer to `PetscFree()`.
481    However, the pointer should never be dereferenced or the program will crash.
482 
483    Developer Note:
484    All the `PetscMallocN()` routines actually call `PetscMalloc()` behind the scenes.
485 
486    Except for data structures that store information about the PETSc options database all memory allocated by PETSc is
487    obtained with `PetscMalloc()` or `PetscCalloc()`
488 
489 .seealso: `PetscFree()`, `PetscNew()`, `PetscCalloc()`
490 M*/
491 #define PetscMalloc(a, b) ((*PetscTrMalloc)((a), PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (void **)(b)))
492 
493 /*MC
494    PetscRealloc - Reallocates memory
495 
496    Synopsis:
497     #include <petscsys.h>
498    PetscErrorCode PetscRealloc(size_t m,void **result)
499 
500    Not Collective
501 
502    Input Parameters:
503 +  m      - number of bytes to allocate
504 -  result - previous memory
505 
506    Output Parameter:
507 .  result - new memory allocated
508 
509    Level: developer
510 
511    Notes:
512    `results` must have already been obtained with `PetscMalloc()`
513 
514    Memory is always allocated at least double aligned
515 
516 .seealso: `PetscMalloc()`, `PetscFree()`, `PetscNew()`
517 M*/
518 #define PetscRealloc(a, b) ((*PetscTrRealloc)((a), __LINE__, PETSC_FUNCTION_NAME, __FILE__, (void **)(b)))
519 
520 /*MC
521    PetscAddrAlign - Rounds up an address to `PETSC_MEMALIGN` alignment
522 
523    Synopsis:
524     #include <petscsys.h>
525    void *PetscAddrAlign(void *addr)
526 
527    Not Collective
528 
529    Input Parameter:
530 .  addr - address to align (any pointer type)
531 
532    Level: developer
533 
534 .seealso: `PetscMallocAlign()`
535 M*/
536 #define PetscAddrAlign(a) ((void *)((((PETSC_UINTPTR_T)(a)) + (PETSC_MEMALIGN - 1)) & ~(PETSC_MEMALIGN - 1)))
537 
538 /*MC
539    PetscCalloc - Allocates a cleared (zeroed) memory region aligned to `PETSC_MEMALIGN`, similar to `PetscMalloc()`
540 
541    Synopsis:
542     #include <petscsys.h>
543    PetscErrorCode PetscCalloc(size_t m,void **result)
544 
545    Not Collective
546 
547    Input Parameter:
548 .  m - number of bytes to allocate
549 
550    Output Parameter:
551 .  result - memory allocated
552 
553    Level: beginner
554 
555    Notes:
556    Memory is always allocated at least double aligned. This macro is useful in allocating memory pointed by void pointers
557 
558    It is safe to allocate with an m of 0 and pass the resulting pointer to `PetscFree()`.
559 
560    However, the pointer should never be dereferenced or the program will crash.
561 
562    Developer Note:
563    All `PetscCallocN()` routines call `PetscCalloc()` behind the scenes.
564 
565 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`
566 M*/
567 #define PetscCalloc(m, result) PetscMallocA(1, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)m), (result))
568 
569 /*MC
570    PetscMalloc1 - Allocates an array of memory aligned to `PETSC_MEMALIGN`
571 
572    Synopsis:
573     #include <petscsys.h>
574    PetscErrorCode PetscMalloc1(size_t m1,type **r1)
575 
576    Not Collective
577 
578    Input Parameter:
579 .  m1 - number of elements to allocate  (may be zero)
580 
581    Output Parameter:
582 .  r1 - memory allocated
583 
584    Level: beginner
585 
586    Note:
587    This uses `sizeof()` of the memory type requested to determine the total memory to be allocated; therefore, you should not
588    multiply the number of elements requested by the `sizeof()` the type. For example, use
589 .vb
590   PetscInt *id;
591   PetscMalloc1(10,&id);
592 .ve
593        not
594 .vb
595   PetscInt *id;
596   PetscMalloc1(10*sizeof(PetscInt),&id);
597 .ve
598 
599   Does not zero the memory allocated, use `PetscCalloc1()` to obtain memory that has been zeroed.
600 
601   The `PetscMalloc[N]()` and `PetscCalloc[N]()` take an argument of type `size_t`! However, most codes use `value`, computed via `int` or `PetscInt` variables. This can overflow in
602   32bit `int` computation - while computation in 64bit `size_t` would not overflow!
603   It's best if any arithmetic that is done for size computations is done with `size_t` type - avoiding arithmetic overflow!
604 
605   `PetscMalloc[N]()` and `PetscCalloc[N]()` attempt to work-around this by casting the first variable to `size_t`.
606   This works for most expressions, but not all, such as
607 .vb
608   PetscInt *id, a, b;
609   PetscMalloc1(use_a_squared ? a * a * b : a * b, &id); // use_a_squared is cast to size_t, but a and b are still PetscInt
610   PetscMalloc1(a + b * b, &id); // a is cast to size_t, but b * b is performed at PetscInt precision first due to order-of-operations
611 .ve
612 
613   These expressions should either be avoided, or appropriately cast variables to `size_t`:
614 .vb
615   PetscInt *id, a, b;
616   PetscMalloc1(use_a_squared ? (size_t)a * a * b : (size_t)a * b, &id); // Cast a to size_t before multiplication
617   PetscMalloc1(b * b + a, &id); // b is automatically cast to size_t and order-of-operations ensures size_t precision is maintained
618 .ve
619 
620 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc2()`
621 M*/
622 #define PetscMalloc1(m1, r1) PetscMallocA(1, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1))
623 
624 /*MC
625    PetscCalloc1 - Allocates a cleared (zeroed) array of memory aligned to `PETSC_MEMALIGN`
626 
627    Synopsis:
628     #include <petscsys.h>
629    PetscErrorCode PetscCalloc1(size_t m1,type **r1)
630 
631    Not Collective
632 
633    Input Parameter:
634 .  m1 - number of elements to allocate in 1st chunk  (may be zero)
635 
636    Output Parameter:
637 .  r1 - memory allocated
638 
639    Level: beginner
640 
641    Note:
642    See `PetscMalloc1()` for more details on usage.
643 
644 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc2()`
645 M*/
646 #define PetscCalloc1(m1, r1) PetscMallocA(1, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1))
647 
648 /*MC
649    PetscMalloc2 - Allocates 2 arrays of memory both aligned to `PETSC_MEMALIGN`
650 
651    Synopsis:
652     #include <petscsys.h>
653    PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)
654 
655    Not Collective
656 
657    Input Parameters:
658 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
659 -  m2 - number of elements to allocate in 2nd chunk  (may be zero)
660 
661    Output Parameters:
662 +  r1 - memory allocated in first chunk
663 -  r2 - memory allocated in second chunk
664 
665    Level: developer
666 
667 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc2()`
668 M*/
669 #define PetscMalloc2(m1, r1, m2, r2) PetscMallocA(2, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1), ((size_t)((size_t)m2) * sizeof(**(r2))), (r2))
670 
671 /*MC
672    PetscCalloc2 - Allocates 2 cleared (zeroed) arrays of memory both aligned to `PETSC_MEMALIGN`
673 
674    Synopsis:
675     #include <petscsys.h>
676    PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)
677 
678    Not Collective
679 
680    Input Parameters:
681 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
682 -  m2 - number of elements to allocate in 2nd chunk  (may be zero)
683 
684    Output Parameters:
685 +  r1 - memory allocated in first chunk
686 -  r2 - memory allocated in second chunk
687 
688    Level: developer
689 
690 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc2()`
691 M*/
692 #define PetscCalloc2(m1, r1, m2, r2) PetscMallocA(2, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1), ((size_t)((size_t)m2) * sizeof(**(r2))), (r2))
693 
694 /*MC
695    PetscMalloc3 - Allocates 3 arrays of memory, all aligned to `PETSC_MEMALIGN`
696 
697    Synopsis:
698     #include <petscsys.h>
699    PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
700 
701    Not Collective
702 
703    Input Parameters:
704 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
705 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
706 -  m3 - number of elements to allocate in 3rd chunk  (may be zero)
707 
708    Output Parameters:
709 +  r1 - memory allocated in first chunk
710 .  r2 - memory allocated in second chunk
711 -  r3 - memory allocated in third chunk
712 
713    Level: developer
714 
715 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc3()`, `PetscFree3()`
716 M*/
717 #define PetscMalloc3(m1, r1, m2, r2, m3, r3) \
718   PetscMallocA(3, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1), ((size_t)((size_t)m2) * sizeof(**(r2))), (r2), ((size_t)((size_t)m3) * sizeof(**(r3))), (r3))
719 
720 /*MC
721    PetscCalloc3 - Allocates 3 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
722 
723    Synopsis:
724     #include <petscsys.h>
725    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
726 
727    Not Collective
728 
729    Input Parameters:
730 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
731 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
732 -  m3 - number of elements to allocate in 3rd chunk  (may be zero)
733 
734    Output Parameters:
735 +  r1 - memory allocated in first chunk
736 .  r2 - memory allocated in second chunk
737 -  r3 - memory allocated in third chunk
738 
739    Level: developer
740 
741 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc2()`, `PetscMalloc3()`, `PetscFree3()`
742 M*/
743 #define PetscCalloc3(m1, r1, m2, r2, m3, r3) \
744   PetscMallocA(3, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1), ((size_t)((size_t)m2) * sizeof(**(r2))), (r2), ((size_t)((size_t)m3) * sizeof(**(r3))), (r3))
745 
746 /*MC
747    PetscMalloc4 - Allocates 4 arrays of memory, all aligned to `PETSC_MEMALIGN`
748 
749    Synopsis:
750     #include <petscsys.h>
751    PetscErrorCode PetscMalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
752 
753    Not Collective
754 
755    Input Parameters:
756 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
757 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
758 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
759 -  m4 - number of elements to allocate in 4th chunk  (may be zero)
760 
761    Output Parameters:
762 +  r1 - memory allocated in first chunk
763 .  r2 - memory allocated in second chunk
764 .  r3 - memory allocated in third chunk
765 -  r4 - memory allocated in fourth chunk
766 
767    Level: developer
768 
769 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc4()`, `PetscFree4()`
770 M*/
771 #define PetscMalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \
772   PetscMallocA(4, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1), ((size_t)((size_t)m2) * sizeof(**(r2))), (r2), ((size_t)((size_t)m3) * sizeof(**(r3))), (r3), ((size_t)((size_t)m4) * sizeof(**(r4))), (r4))
773 
774 /*MC
775    PetscCalloc4 - Allocates 4 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
776 
777    Synopsis:
778     #include <petscsys.h>
779    PetscErrorCode PetscCalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
780 
781    Not Collective
782 
783    Input Parameters:
784 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
785 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
786 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
787 -  m4 - number of elements to allocate in 4th chunk  (may be zero)
788 
789    Output Parameters:
790 +  r1 - memory allocated in first chunk
791 .  r2 - memory allocated in second chunk
792 .  r3 - memory allocated in third chunk
793 -  r4 - memory allocated in fourth chunk
794 
795    Level: developer
796 
797 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc4()`, `PetscFree4()`
798 M*/
799 #define PetscCalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \
800   PetscMallocA(4, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1), ((size_t)((size_t)m2) * sizeof(**(r2))), (r2), ((size_t)((size_t)m3) * sizeof(**(r3))), (r3), ((size_t)((size_t)m4) * sizeof(**(r4))), (r4))
801 
802 /*MC
803    PetscMalloc5 - Allocates 5 arrays of memory, all aligned to `PETSC_MEMALIGN`
804 
805    Synopsis:
806     #include <petscsys.h>
807    PetscErrorCode PetscMalloc5(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5)
808 
809    Not Collective
810 
811    Input Parameters:
812 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
813 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
814 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
815 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
816 -  m5 - number of elements to allocate in 5th chunk  (may be zero)
817 
818    Output Parameters:
819 +  r1 - memory allocated in first chunk
820 .  r2 - memory allocated in second chunk
821 .  r3 - memory allocated in third chunk
822 .  r4 - memory allocated in fourth chunk
823 -  r5 - memory allocated in fifth chunk
824 
825    Level: developer
826 
827 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc5()`, `PetscFree5()`
828 M*/
829 #define PetscMalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \
830   PetscMallocA(5, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1), ((size_t)((size_t)m2) * sizeof(**(r2))), (r2), ((size_t)((size_t)m3) * sizeof(**(r3))), (r3), ((size_t)((size_t)m4) * sizeof(**(r4))), (r4), ((size_t)((size_t)m5) * sizeof(**(r5))), (r5))
831 
832 /*MC
833    PetscCalloc5 - Allocates 5 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
834 
835    Synopsis:
836     #include <petscsys.h>
837    PetscErrorCode PetscCalloc5(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5)
838 
839    Not Collective
840 
841    Input Parameters:
842 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
843 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
844 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
845 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
846 -  m5 - number of elements to allocate in 5th chunk  (may be zero)
847 
848    Output Parameters:
849 +  r1 - memory allocated in first chunk
850 .  r2 - memory allocated in second chunk
851 .  r3 - memory allocated in third chunk
852 .  r4 - memory allocated in fourth chunk
853 -  r5 - memory allocated in fifth chunk
854 
855    Level: developer
856 
857 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc5()`, `PetscFree5()`
858 M*/
859 #define PetscCalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \
860   PetscMallocA(5, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1), ((size_t)((size_t)m2) * sizeof(**(r2))), (r2), ((size_t)((size_t)m3) * sizeof(**(r3))), (r3), ((size_t)((size_t)m4) * sizeof(**(r4))), (r4), ((size_t)((size_t)m5) * sizeof(**(r5))), (r5))
861 
862 /*MC
863    PetscMalloc6 - Allocates 6 arrays of memory, all aligned to `PETSC_MEMALIGN`
864 
865    Synopsis:
866     #include <petscsys.h>
867    PetscErrorCode PetscMalloc6(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5,size_t m6,type **r6)
868 
869    Not Collective
870 
871    Input Parameters:
872 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
873 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
874 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
875 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
876 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
877 -  m6 - number of elements to allocate in 6th chunk  (may be zero)
878 
879    Output Parameteasr:
880 +  r1 - memory allocated in first chunk
881 .  r2 - memory allocated in second chunk
882 .  r3 - memory allocated in third chunk
883 .  r4 - memory allocated in fourth chunk
884 .  r5 - memory allocated in fifth chunk
885 -  r6 - memory allocated in sixth chunk
886 
887    Level: developer
888 
889 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc6()`, `PetscFree3()`, `PetscFree4()`, `PetscFree5()`, `PetscFree6()`
890 M*/
891 #define PetscMalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \
892   PetscMallocA(6, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1), ((size_t)((size_t)m2) * sizeof(**(r2))), (r2), ((size_t)((size_t)m3) * sizeof(**(r3))), (r3), ((size_t)((size_t)m4) * sizeof(**(r4))), (r4), ((size_t)((size_t)m5) * sizeof(**(r5))), (r5), ((size_t)((size_t)m6) * sizeof(**(r6))), (r6))
893 
894 /*MC
895    PetscCalloc6 - Allocates 6 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
896 
897    Synopsis:
898     #include <petscsys.h>
899    PetscErrorCode PetscCalloc6(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5,size_t m6,type **r6)
900 
901    Not Collective
902 
903    Input Parameters:
904 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
905 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
906 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
907 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
908 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
909 -  m6 - number of elements to allocate in 6th chunk  (may be zero)
910 
911    Output Parameters:
912 +  r1 - memory allocated in first chunk
913 .  r2 - memory allocated in second chunk
914 .  r3 - memory allocated in third chunk
915 .  r4 - memory allocated in fourth chunk
916 .  r5 - memory allocated in fifth chunk
917 -  r6 - memory allocated in sixth chunk
918 
919    Level: developer
920 
921 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscMalloc6()`, `PetscFree6()`
922 M*/
923 #define PetscCalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \
924   PetscMallocA(6, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1), ((size_t)((size_t)m2) * sizeof(**(r2))), (r2), ((size_t)((size_t)m3) * sizeof(**(r3))), (r3), ((size_t)((size_t)m4) * sizeof(**(r4))), (r4), ((size_t)((size_t)m5) * sizeof(**(r5))), (r5), ((size_t)((size_t)m6) * sizeof(**(r6))), (r6))
925 
926 /*MC
927    PetscMalloc7 - Allocates 7 arrays of memory, all aligned to `PETSC_MEMALIGN`
928 
929    Synopsis:
930     #include <petscsys.h>
931    PetscErrorCode PetscMalloc7(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5,size_t m6,type **r6,size_t m7,type **r7)
932 
933    Not Collective
934 
935    Input Parameters:
936 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
937 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
938 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
939 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
940 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
941 .  m6 - number of elements to allocate in 6th chunk  (may be zero)
942 -  m7 - number of elements to allocate in 7th chunk  (may be zero)
943 
944    Output Parameters:
945 +  r1 - memory allocated in first chunk
946 .  r2 - memory allocated in second chunk
947 .  r3 - memory allocated in third chunk
948 .  r4 - memory allocated in fourth chunk
949 .  r5 - memory allocated in fifth chunk
950 .  r6 - memory allocated in sixth chunk
951 -  r7 - memory allocated in seventh chunk
952 
953    Level: developer
954 
955 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc7()`, `PetscFree7()`
956 M*/
957 #define PetscMalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \
958   PetscMallocA(7, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1), ((size_t)((size_t)m2) * sizeof(**(r2))), (r2), ((size_t)((size_t)m3) * sizeof(**(r3))), (r3), ((size_t)((size_t)m4) * sizeof(**(r4))), (r4), ((size_t)((size_t)m5) * sizeof(**(r5))), (r5), ((size_t)((size_t)m6) * sizeof(**(r6))), (r6), ((size_t)((size_t)m7) * sizeof(**(r7))), (r7))
959 
960 /*MC
961    PetscCalloc7 - Allocates 7 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
962 
963    Synopsis:
964     #include <petscsys.h>
965    PetscErrorCode PetscCalloc7(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5,size_t m6,type **r6,size_t m7,type **r7)
966 
967    Not Collective
968 
969    Input Parameters:
970 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
971 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
972 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
973 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
974 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
975 .  m6 - number of elements to allocate in 6th chunk  (may be zero)
976 -  m7 - number of elements to allocate in 7th chunk  (may be zero)
977 
978    Output Parameters:
979 +  r1 - memory allocated in first chunk
980 .  r2 - memory allocated in second chunk
981 .  r3 - memory allocated in third chunk
982 .  r4 - memory allocated in fourth chunk
983 .  r5 - memory allocated in fifth chunk
984 .  r6 - memory allocated in sixth chunk
985 -  r7 - memory allocated in seventh chunk
986 
987    Level: developer
988 
989 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscMalloc7()`, `PetscFree7()`
990 M*/
991 #define PetscCalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \
992   PetscMallocA(7, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1), ((size_t)((size_t)m2) * sizeof(**(r2))), (r2), ((size_t)((size_t)m3) * sizeof(**(r3))), (r3), ((size_t)((size_t)m4) * sizeof(**(r4))), (r4), ((size_t)((size_t)m5) * sizeof(**(r5))), (r5), ((size_t)((size_t)m6) * sizeof(**(r6))), (r6), ((size_t)((size_t)m7) * sizeof(**(r7))), (r7))
993 
994 /*MC
995    PetscNew - Allocates memory of a particular type, zeros the memory! Aligned to `PETSC_MEMALIGN`
996 
997    Synopsis:
998     #include <petscsys.h>
999    PetscErrorCode PetscNew(type **result)
1000 
1001    Not Collective
1002 
1003    Output Parameter:
1004 .  result - memory allocated, sized to match pointer `type`
1005 
1006    Level: beginner
1007 
1008    Developer Note:
1009    Calls `PetscCalloc()` with the appropriate memory size obtained from `type`
1010 
1011 .seealso: `PetscFree()`, `PetscMalloc()`, `PetscCall()`, `PetscCalloc1()`, `PetscMalloc1()`
1012 M*/
1013 #define PetscNew(b) PetscCalloc1(1, (b))
1014 
1015 #define PetscNewLog(o, b) PETSC_DEPRECATED_MACRO(3, 18, 0, "PetscNew()", ) PetscNew(b)
1016 
1017 /*MC
1018    PetscFree - Frees memory
1019 
1020    Synopsis:
1021     #include <petscsys.h>
1022    PetscErrorCode PetscFree(void *memory)
1023 
1024    Not Collective
1025 
1026    Input Parameter:
1027 .   memory - memory to free (the pointer is ALWAYS set to `NULL` upon success)
1028 
1029    Level: beginner
1030 
1031    Notes:
1032    Do not free memory obtained with `PetscMalloc2()`, `PetscCalloc2()` etc, they must be freed with `PetscFree2()` etc.
1033 
1034    It is safe to call `PetscFree()` on a `NULL` pointer.
1035 
1036 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc1()`
1037 M*/
1038 #define PetscFree(a) ((PetscErrorCode)((*PetscTrFree)((void *)(a), __LINE__, PETSC_FUNCTION_NAME, __FILE__) || ((a) = PETSC_NULLPTR, PETSC_SUCCESS)))
1039 
1040 /*MC
1041    PetscFree2 - Frees 2 chunks of memory obtained with `PetscMalloc2()`
1042 
1043    Synopsis:
1044     #include <petscsys.h>
1045    PetscErrorCode PetscFree2(void *memory1,void *memory2)
1046 
1047    Not Collective
1048 
1049    Input Parameters:
1050 +   memory1 - memory to free
1051 -   memory2 - 2nd memory to free
1052 
1053    Level: developer
1054 
1055    Notes:
1056     Memory must have been obtained with `PetscMalloc2()`
1057 
1058     The arguments need to be in the same order as they were in the call to `PetscMalloc2()`
1059 
1060 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`
1061 M*/
1062 #define PetscFree2(m1, m2) PetscFreeA(2, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2))
1063 
1064 /*MC
1065    PetscFree3 - Frees 3 chunks of memory obtained with `PetscMalloc3()`
1066 
1067    Synopsis:
1068     #include <petscsys.h>
1069    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)
1070 
1071    Not Collective
1072 
1073    Input Parameters:
1074 +   memory1 - memory to free
1075 .   memory2 - 2nd memory to free
1076 -   memory3 - 3rd memory to free
1077 
1078    Level: developer
1079 
1080    Notes:
1081     Memory must have been obtained with `PetscMalloc3()`
1082 
1083     The arguments need to be in the same order as they were in the call to `PetscMalloc3()`
1084 
1085 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`
1086 M*/
1087 #define PetscFree3(m1, m2, m3) PetscFreeA(3, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3))
1088 
1089 /*MC
1090    PetscFree4 - Frees 4 chunks of memory obtained with `PetscMalloc4()`
1091 
1092    Synopsis:
1093     #include <petscsys.h>
1094    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)
1095 
1096    Not Collective
1097 
1098    Input Parameters:
1099 +   m1 - memory to free
1100 .   m2 - 2nd memory to free
1101 .   m3 - 3rd memory to free
1102 -   m4 - 4th memory to free
1103 
1104    Level: developer
1105 
1106    Notes:
1107     Memory must have been obtained with `PetscMalloc4()`
1108 
1109     The arguments need to be in the same order as they were in the call to `PetscMalloc4()`
1110 
1111 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`
1112 M*/
1113 #define PetscFree4(m1, m2, m3, m4) PetscFreeA(4, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4))
1114 
1115 /*MC
1116    PetscFree5 - Frees 5 chunks of memory obtained with `PetscMalloc5()`
1117 
1118    Synopsis:
1119     #include <petscsys.h>
1120    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)
1121 
1122    Not Collective
1123 
1124    Input Parameters:
1125 +   m1 - memory to free
1126 .   m2 - 2nd memory to free
1127 .   m3 - 3rd memory to free
1128 .   m4 - 4th memory to free
1129 -   m5 - 5th memory to free
1130 
1131    Level: developer
1132 
1133    Notes:
1134     Memory must have been obtained with `PetscMalloc5()`
1135 
1136     The arguments need to be in the same order as they were in the call to `PetscMalloc5()`
1137 
1138 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`
1139 M*/
1140 #define PetscFree5(m1, m2, m3, m4, m5) PetscFreeA(5, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5))
1141 
1142 /*MC
1143    PetscFree6 - Frees 6 chunks of memory obtained with `PetscMalloc6()`
1144 
1145    Synopsis:
1146     #include <petscsys.h>
1147    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)
1148 
1149    Not Collective
1150 
1151    Input Parameters:
1152 +   m1 - memory to free
1153 .   m2 - 2nd memory to free
1154 .   m3 - 3rd memory to free
1155 .   m4 - 4th memory to free
1156 .   m5 - 5th memory to free
1157 -   m6 - 6th memory to free
1158 
1159    Level: developer
1160 
1161    Notes:
1162     Memory must have been obtained with `PetscMalloc6()`
1163 
1164     The arguments need to be in the same order as they were in the call to `PetscMalloc6()`
1165 
1166 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()`
1167 M*/
1168 #define PetscFree6(m1, m2, m3, m4, m5, m6) PetscFreeA(6, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6))
1169 
1170 /*MC
1171    PetscFree7 - Frees 7 chunks of memory obtained with `PetscMalloc7()`
1172 
1173    Synopsis:
1174     #include <petscsys.h>
1175    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)
1176 
1177    Not Collective
1178 
1179    Input Parameters:
1180 +   m1 - memory to free
1181 .   m2 - 2nd memory to free
1182 .   m3 - 3rd memory to free
1183 .   m4 - 4th memory to free
1184 .   m5 - 5th memory to free
1185 .   m6 - 6th memory to free
1186 -   m7 - 7th memory to free
1187 
1188    Level: developer
1189 
1190    Notes:
1191     Memory must have been obtained with `PetscMalloc7()`
1192 
1193     The arguments need to be in the same order as they were in the call to `PetscMalloc7()`
1194 
1195 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()`,
1196           `PetscMalloc7()`
1197 M*/
1198 #define PetscFree7(m1, m2, m3, m4, m5, m6, m7) PetscFreeA(7, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6), &(m7))
1199 
1200 PETSC_EXTERN PetscErrorCode PetscMallocA(int, PetscBool, int, const char *, const char *, size_t, void *, ...);
1201 PETSC_EXTERN PetscErrorCode PetscFreeA(int, int, const char *, const char *, void *, ...);
1202 PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t, PetscBool, int, const char[], const char[], void **);
1203 PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void *, int, const char[], const char[]);
1204 PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t, int, const char[], const char[], void **);
1205 PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1206 PETSC_EXTERN PetscErrorCode PetscMallocSet(PetscErrorCode (*)(size_t, PetscBool, int, const char[], const char[], void **), PetscErrorCode (*)(void *, int, const char[], const char[]), PetscErrorCode (*)(size_t, int, const char[], const char[], void **));
1207 PETSC_EXTERN PetscErrorCode PetscMallocClear(void);
1208 
1209 /*
1210   Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1211 */
1212 PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1213 PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);
1214 #if defined(PETSC_HAVE_CUDA)
1215 PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void);
1216 PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void);
1217 #endif
1218 #if defined(PETSC_HAVE_HIP)
1219 PETSC_EXTERN PetscErrorCode PetscMallocSetHIPHost(void);
1220 PETSC_EXTERN PetscErrorCode PetscMallocResetHIPHost(void);
1221 #endif
1222 
1223 #define MPIU_PETSCLOGDOUBLE  MPI_DOUBLE
1224 #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION
1225 
1226 /*
1227    Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1228 */
1229 PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1230 PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1231 PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1232 PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1233 PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1234 PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int, PetscLogDouble *);
1235 PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool, PetscBool);
1236 PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool *, PetscBool *, PetscBool *);
1237 PETSC_EXTERN PetscErrorCode PetscMallocValidate(int, const char[], const char[]);
1238 PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1239 PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool *);
1240 PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool);
1241 PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool *);
1242 
1243 PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType, MPI_Datatype *);
1244 PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype, PetscDataType *);
1245 PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType, size_t *);
1246 PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char *, PetscDataType *, PetscBool *);
1247 
1248 /*
1249    These are MPI operations for MPI_Allreduce() etc
1250 */
1251 PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1252 #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1253 PETSC_EXTERN MPI_Op MPIU_SUM;
1254 PETSC_EXTERN MPI_Op MPIU_MAX;
1255 PETSC_EXTERN MPI_Op MPIU_MIN;
1256 #else
1257   #define MPIU_SUM MPI_SUM
1258   #define MPIU_MAX MPI_MAX
1259   #define MPIU_MIN MPI_MIN
1260 #endif
1261 PETSC_EXTERN MPI_Op         Petsc_Garbage_SetIntersectOp;
1262 PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm, const PetscInt[], PetscInt *, PetscInt *);
1263 
1264 #if (defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128)) || (defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16))
1265 /*MC
1266    MPIU_SUM___FP16___FLOAT128 - MPI_Op that acts as a replacement for `MPI_SUM` with
1267    custom `MPI_Datatype` `MPIU___FLOAT128`, `MPIU___COMPLEX128`, and `MPIU___FP16`.
1268 
1269    Level: advanced
1270 
1271    Developer Note:
1272    This should be unified with `MPIU_SUM`
1273 
1274 .seealso: `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`
1275 M*/
1276 PETSC_EXTERN MPI_Op MPIU_SUM___FP16___FLOAT128;
1277 #endif
1278 
1279 /*
1280      These are so that in extern C code we can cast function pointers to non-extern C
1281    function pointers. Since the regular C++ code expects its function pointers to be C++
1282 */
1283 
1284 /*S
1285   PetscVoidFn - A prototype of a void (fn)(void) function
1286 
1287   Level: developer
1288 
1289   Notes:
1290   The deprecated `PetscVoidFunction` works as a replacement for `PetscVoidFn` *.
1291 
1292   The deprecated `PetscVoidStarFunction` works as a replacement for `PetscVoidFn` **.
1293 
1294 .seealso: `PetscObject`, `PetscObjectDestroy()`
1295 S*/
1296 PETSC_EXTERN_TYPEDEF typedef void PetscVoidFn(void);
1297 
1298 PETSC_EXTERN_TYPEDEF typedef PetscVoidFn  *PetscVoidFunction;
1299 PETSC_EXTERN_TYPEDEF typedef PetscVoidFn **PetscVoidStarFunction;
1300 
1301 /*S
1302   PetscErrorCodeFn - A prototype of a PetscErrorCode (fn)(void) function
1303 
1304   Level: developer
1305 
1306   Notes:
1307   The deprecated `PetscErrorCodeFunction` works as a replacement for `PetscErrorCodeFn` *.
1308 
1309 .seealso: `PetscObject`, `PetscObjectDestroy()`
1310 S*/
1311 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PetscErrorCodeFn(void);
1312 
1313 PETSC_EXTERN_TYPEDEF typedef PetscErrorCodeFn *PetscErrorCodeFunction;
1314 
1315 /*S
1316   PetscCtxDestroyFn - A prototype of a `PetscErrorCode (*)(void **)` function that is used to free user contexts
1317 
1318   Level: intermediate
1319 
1320   Note:
1321   Used in the prototype of functions such as `DMSetApplicationContextDestroy()`
1322 
1323 .seealso: `PetscObject`, `PetscCtxDestroyDefault()`, `PetscObjectDestroy()`, `DMSetApplicationContextDestroy()`
1324 S*/
1325 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PetscCtxDestroyFn(void **);
1326 
1327 PETSC_EXTERN PetscCtxDestroyFn PetscCtxDestroyDefault;
1328 
1329 /*
1330     Defines PETSc error handling.
1331 */
1332 #include <petscerror.h> // IWYU pragma: export
1333 
1334 PETSC_EXTERN PetscBool   PetscCIEnabled;                    /* code is running in the PETSc test harness CI */
1335 PETSC_EXTERN PetscBool   PetscCIEnabledPortableErrorOutput; /* error output is stripped to ensure portability of error messages across systems */
1336 PETSC_EXTERN const char *PetscCIFilename(const char *);
1337 PETSC_EXTERN int         PetscCILinenumber(int);
1338 
1339 #define PETSC_SMALLEST_CLASSID 1211211
1340 PETSC_EXTERN PetscClassId   PETSC_LARGEST_CLASSID;
1341 PETSC_EXTERN PetscClassId   PETSC_OBJECT_CLASSID;
1342 PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[], PetscClassId *);
1343 PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject, PetscObjectId *);
1344 PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject, PetscObjectId, PetscBool *);
1345 
1346 /*
1347    Routines that get memory usage information from the OS
1348 */
1349 PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1350 PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1351 PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1352 PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);
1353 
1354 PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);
1355 
1356 /*
1357    Initialization of PETSc
1358 */
1359 PETSC_EXTERN PetscErrorCode PetscInitialize(int *, char ***, const char[], const char[]);
1360 PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int, char *[], const char[], const char[]);
1361 PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1362 PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1363 PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1364 PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1365 PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1366 PETSC_EXTERN PetscErrorCode PetscGetArgs(int *, char ***);
1367 PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1368 PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);
1369 
1370 PETSC_EXTERN PetscErrorCode PetscEnd(void);
1371 PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);
1372 PETSC_EXTERN PetscErrorCode PetscSysFinalizePackage(void);
1373 
1374 PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[], const char[]);
1375 PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1376 PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1377 PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject, const char[]);
1378 
1379 PETSC_EXTERN PetscErrorCode PetscMonitorCompare(PetscErrorCode (*)(void), void *, PetscCtxDestroyFn *, PetscErrorCode (*)(void), void *, PetscCtxDestroyFn *, PetscBool *);
1380 
1381 /*
1382     Functions that can act on any PETSc object.
1383 */
1384 PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject *);
1385 PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject, MPI_Comm *);
1386 PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject, PetscClassId *);
1387 PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject, const char *[]);
1388 PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject, const char *[]);
1389 PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject, const char[]);
1390 PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject, const char *[]);
1391 PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject, PetscInt);
1392 PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject, PetscInt *);
1393 PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject, PetscObject, PetscInt);
1394 PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1395 PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject, PetscInt *);
1396 PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1397 PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject, PetscMPIInt *);
1398 PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject, const char[], PetscObject);
1399 PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject, const char[]);
1400 PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject, const char[], PetscObject *);
1401 PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject, const char[], void (*)(void));
1402 #define PetscObjectComposeFunction(a, b, ...) PetscObjectComposeFunction_Private((a), (b), (PetscVoidFn *)(__VA_ARGS__))
1403 PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1404 PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1405 PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1406 PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject, PetscObject);
1407 PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm, PetscMPIInt *);
1408 
1409 /*MC
1410    PetscObjectParameterSetDefault - sets a parameter default value in a `PetscObject` to a new default value.
1411    If the current value matches the old default value, then the current value is also set to the new value.
1412 
1413    No Fortran Support
1414 
1415    Synopsis:
1416    #include <petscsys.h>
1417    PetscBool PetscObjectParameterSetDefault(PetscObject obj, char* NAME, PetscReal value);
1418 
1419    Input Parameters:
1420 +  obj - the `PetscObject`
1421 .  NAME - the name of the parameter, unquoted
1422 -  value - the new value
1423 
1424    Level: developer
1425 
1426    Notes:
1427    The defaults for an object are the values set when the object's type is set.
1428 
1429    This should only be used in object constructors, such as, `SNESCreate_NGS()`.
1430 
1431    This only works for parameters that are declared in the struct with `PetscObjectParameterDeclare()`
1432 
1433 .seealso: `PetscObjectParameterDeclare()`, `PetscInitialize()`, `PetscFinalize()`, `PetscObject`, `SNESParametersInitialize()`
1434 M*/
1435 #define PetscObjectParameterSetDefault(obj, NAME, value) \
1436   do { \
1437     if (obj->NAME == obj->default_##NAME) obj->NAME = value; \
1438     obj->default_##NAME = value; \
1439   } while (0)
1440 
1441 /*MC
1442    PetscObjectParameterDeclare - declares a parameter in a `PetscObject` and a location to store its default
1443 
1444    No Fortran Support
1445 
1446    Synopsis:
1447    #include <petscsys.h>
1448    PetscBool PetscObjectParameterDeclare(type, char* NAME)
1449 
1450    Input Parameters:
1451 +  type - the type of the parameter, for example `PetscInt`
1452 -  NAME - the name of the parameter, unquoted
1453 
1454    Level: developer.
1455 
1456 .seealso: `PetscObjectParameterSetDefault()`, `PetscInitialize()`, `PetscFinalize()`, `PetscObject`, `SNESParametersInitialize()`
1457 M*/
1458 #define PetscObjectParameterDeclare(type, NAME)    type NAME, default_##NAME
1459 #define PetscObjectParameterDeclarePtr(type, NAME) type *NAME, *default_##NAME
1460 
1461 #include <petscviewertypes.h>
1462 #include <petscoptions.h>
1463 
1464 PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer, PetscBool, PetscLogDouble);
1465 PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool *);
1466 
1467 PETSC_EXTERN PetscErrorCode PetscObjectsListGetGlobalNumbering(MPI_Comm, PetscInt, PetscObject[], PetscInt *, PetscInt *);
1468 
1469 PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer, const char[]);
1470 PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject, PetscViewer);
1471 PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject, PetscViewer);
1472 #define PetscObjectQueryFunction(obj, name, fptr) PetscObjectQueryFunction_Private((obj), (name), (PetscVoidFn **)(fptr))
1473 PETSC_EXTERN PetscErrorCode PetscObjectHasFunction(PetscObject, const char[], PetscBool *);
1474 PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject, const char[], void (**)(void));
1475 PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject, const char[]);
1476 PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject, const char[]);
1477 PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject, const char[]);
1478 PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject, const char *[]);
1479 PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject, const char[]);
1480 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1481 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1482 PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject, PetscObject, const char[]);
1483 PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1484 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject, const char[], PetscBool *);
1485 PETSC_EXTERN PetscErrorCode PetscObjectObjectTypeCompare(PetscObject, PetscObject, PetscBool *);
1486 PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject, const char[], PetscBool *);
1487 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject, PetscBool *, const char[], ...);
1488 PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject, PetscBool *, const char[], ...);
1489 PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1490 PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);
1491 
1492 #if defined(PETSC_HAVE_SAWS)
1493 PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1494 PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1495 PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject, PetscBool);
1496 PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1497 PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1498 PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1499 PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1500 PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1501 PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1502 PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);
1503 
1504 #else
1505   #define PetscSAWsBlock()                  PETSC_SUCCESS
1506   #define PetscObjectSAWsViewOff(obj)       PETSC_SUCCESS
1507   #define PetscObjectSAWsSetBlock(obj, flg) PETSC_SUCCESS
1508   #define PetscObjectSAWsBlock(obj)         PETSC_SUCCESS
1509   #define PetscObjectSAWsGrantAccess(obj)   PETSC_SUCCESS
1510   #define PetscObjectSAWsTakeAccess(obj)    PETSC_SUCCESS
1511   #define PetscStackViewSAWs()              PETSC_SUCCESS
1512   #define PetscStackSAWsViewOff()           PETSC_SUCCESS
1513   #define PetscStackSAWsTakeAccess()
1514   #define PetscStackSAWsGrantAccess()
1515 
1516 #endif
1517 
1518 PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[], PetscDLMode, PetscDLHandle *);
1519 PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1520 PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle, const char[], void **);
1521 PETSC_EXTERN PetscErrorCode PetscDLAddr(void (*)(void), char *[]);
1522 #ifdef PETSC_HAVE_CXX
1523 PETSC_EXTERN PetscErrorCode PetscDemangleSymbol(const char *, char *[]);
1524 #endif
1525 
1526 PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void *, PetscStack **);
1527 
1528 PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE *, PetscBool);
1529 PETSC_EXTERN PetscErrorCode PetscObjectsView(PetscViewer);
1530 PETSC_EXTERN PetscErrorCode PetscObjectsGetObject(const char *, PetscObject *, const char *[]);
1531 PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList *);
1532 PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList, const char[], PetscObject *);
1533 PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList, PetscObject, const char *[], PetscBool *);
1534 PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *, const char[], PetscObject);
1535 PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *, const char[]);
1536 PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList, PetscObjectList *);
1537 
1538 /*
1539     Dynamic library lists. Lists of names of routines in objects or in dynamic
1540   link libraries that will be loaded as needed.
1541 */
1542 
1543 #define PetscFunctionListAdd(list, name, fptr) PetscFunctionListAdd_Private((list), (name), (PetscVoidFn *)(fptr))
1544 PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList *, const char[], PetscVoidFn *);
1545 PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList *);
1546 PETSC_EXTERN PetscErrorCode PetscFunctionListClear(PetscFunctionList);
1547 #define PetscFunctionListFind(list, name, fptr) PetscFunctionListFind_Private((list), (name), (PetscVoidFn **)(fptr))
1548 PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList, const char[], PetscVoidFn **);
1549 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm, FILE *, const char[], const char[], const char[], const char[], PetscFunctionList, const char[], const char[]);
1550 PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList, PetscFunctionList *);
1551 PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList, PetscViewer);
1552 PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList, const char ***, int *);
1553 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintNonEmpty(PetscFunctionList);
1554 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintAll(void);
1555 
1556 PETSC_EXTERN PetscDLLibrary PetscDLLibrariesLoaded;
1557 PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm, PetscDLLibrary *, const char[]);
1558 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm, PetscDLLibrary *, const char[]);
1559 PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm, PetscDLLibrary *, const char[], const char[], void **);
1560 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1561 PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm, const char[], char *, size_t, PetscBool *);
1562 PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm, const char[], PetscDLLibrary *);
1563 PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);
1564 
1565 /*
1566      Useful utility routines
1567 */
1568 PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm, PetscInt *, PetscInt *);
1569 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm, PetscInt, PetscInt *, PetscInt *);
1570 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm, PetscInt *, PetscInt *);
1571 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm, PetscMPIInt);
1572 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm, PetscMPIInt);
1573 PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1574 PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE *);
1575 PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm, const PetscInt[2], PetscInt[2]);
1576 PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm, const PetscReal[2], PetscReal[2]);
1577 
1578 /*MC
1579     PetscNot - negates a logical type value and returns result as a `PetscBool`
1580 
1581     Level: beginner
1582 
1583     Note:
1584     This is useful in cases like
1585 .vb
1586      int        *a;
1587      PetscBool  flag = PetscNot(a)
1588 .ve
1589      where !a would not return a `PetscBool` because we cannot provide a cast from int to `PetscBool` in C.
1590 
1591 .seealso: `PetscBool`, `PETSC_TRUE`, `PETSC_FALSE`
1592 M*/
1593 #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)
1594 
1595 /*MC
1596    PetscHelpPrintf - Prints help messages.
1597 
1598    Synopsis:
1599     #include <petscsys.h>
1600      PetscErrorCode (*PetscHelpPrintf)(MPI_Comm comm, const char format[],args);
1601 
1602    Not Collective, only applies on MPI rank 0; No Fortran Support
1603 
1604    Input Parameters:
1605 +  comm - the MPI communicator over which the help message is printed
1606 .  format - the usual printf() format string
1607 -  args - arguments to be printed
1608 
1609    Level: developer
1610 
1611    Notes:
1612    You can change how help messages are printed by replacing the function pointer with a function that does not simply write to stdout.
1613 
1614    To use, write your own function, for example,
1615 .vb
1616    PetscErrorCode mypetschelpprintf(MPI_Comm comm,const char format[],....)
1617    {
1618      PetscFunctionReturn(PETSC_SUCCESS);
1619    }
1620 .ve
1621 then do the assignment
1622 .vb
1623   PetscHelpPrintf = mypetschelpprintf;
1624 .ve
1625 
1626   You can do the assignment before `PetscInitialize()`.
1627 
1628   The default routine used is called `PetscHelpPrintfDefault()`.
1629 
1630 .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscErrorPrintf()`, `PetscHelpPrintfDefault()`
1631 M*/
1632 PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1633 
1634 /*
1635      Defines PETSc profiling.
1636 */
1637 #include <petsclog.h>
1638 
1639 /*
1640       Simple PETSc parallel IO for ASCII printing
1641 */
1642 PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[], char[]);
1643 PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm, const char[], const char[], FILE **);
1644 PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm, FILE *);
1645 PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1646 PETSC_EXTERN PetscErrorCode PetscFFlush(FILE *);
1647 PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1648 PETSC_EXTERN PetscErrorCode PetscSNPrintf(char *, size_t, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1649 PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char *, size_t, const char[], size_t *, ...) PETSC_ATTRIBUTE_FORMAT(3, 5);
1650 PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[], size_t, const char *, PetscInt, const PetscReal[]);
1651 
1652 PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1653 PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1654 PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1655 
1656 PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char *, size_t *);
1657 PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char *, char *);
1658 
1659 PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm, const char[], const char[], const char[], FILE **);
1660 PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm, FILE *);
1661 PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1662 
1663 PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1664 PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1665 PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm, FILE *);
1666 PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm, FILE *, size_t, char[]);
1667 PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm, const char[], const char[], FILE **);
1668 PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char *[]);
1669 
1670 /*MC
1671    PeCtx - indicates an argument that returns a pointer to a user defined C struct (or Fortran derived type)
1672 
1673    Level: developer
1674 
1675    Notes:
1676    This is not part of the PETSc public API and should only be used in PETSc source code.
1677 
1678    This should not be used for functions that return PETSc objects, or pointers to arrays of unknown type. Thus it is used for, for example,
1679    `KSPGetApplicationContext()` but not used for `DMNetworkGetComponent()`
1680 
1681    For pointers to arrays of unknown type and for functions that return PETSc internal objects that are opaque to users, such
1682    as `KSPMonitorDynamicToleranceCreate()` a `void **` should be used.
1683 
1684    Fortran Note:
1685    Should only be used with user defined Fortran datatypes
1686 .vb
1687    type(tUserType), pointer :: ctx
1688 .ve
1689 
1690    Developer Note:
1691    Put this in function declaration for the argument type instead of `void *`, or `void **`.
1692 
1693    C compilers generate a warning or error if one passes a pointer to a pointer to a specific type (instead of `void`), for example,
1694 .vb
1695    extern calledfunction(void **);
1696    SomeCtx *ctx;
1697    calledfunction(&ctx);   << warning that it is passing a pointer to a pointer to a SomeCtx instead of a void **
1698 .ve
1699    By using the common practice of prototyping the function as
1700 .vb
1701    extern calledfunction(void *);
1702 .ve
1703    the warning message is averted. `PeCtx` is used in PETSc source code so that the getAPI() code processor knows the argument is
1704    actually handled internally as `void **` so it can generate correct bindings for other languages.
1705 
1706 .seealso: `PeOp`, `PeNS`, `PetscInitialize()`
1707 M*/
1708 typedef void *PeCtx;
1709 
1710 PETSC_EXTERN PetscClassId   PETSC_CONTAINER_CLASSID;
1711 PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer, void *);
1712 PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer, void *);
1713 PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer *);
1714 PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm, PetscContainer *);
1715 PETSC_EXTERN PetscErrorCode PetscContainerSetCtxDestroy(PetscContainer, PetscCtxDestroyFn *);
1716 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 23, 0, "PetscContainerSetCtxDestroy()", ) PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void *));
1717 PETSC_EXTERN PetscErrorCode PetscObjectContainerCompose(PetscObject, const char *name, void *, PetscCtxDestroyFn *);
1718 PETSC_EXTERN PetscErrorCode PetscObjectContainerQuery(PetscObject, const char *, PeCtx);
1719 
1720 PETSC_DEPRECATED_FUNCTION(3, 23, 0, "PetscCtxDestroyDefault()", ) static inline PetscErrorCode PetscContainerCtxDestroyDefault(void **a)
1721 {
1722   return PetscCtxDestroyDefault(a);
1723 }
1724 
1725 /*
1726    For use in debuggers
1727 */
1728 PETSC_EXTERN PetscMPIInt    PetscGlobalRank;
1729 PETSC_EXTERN PetscMPIInt    PetscGlobalSize;
1730 PETSC_EXTERN PetscErrorCode PetscIntViewNumColumns(PetscInt, PetscInt, const PetscInt[], PetscViewer);
1731 PETSC_EXTERN PetscErrorCode PetscRealViewNumColumns(PetscInt, PetscInt, const PetscReal[], PetscViewer);
1732 PETSC_EXTERN PetscErrorCode PetscScalarViewNumColumns(PetscInt, PetscInt, const PetscScalar[], PetscViewer);
1733 PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt, const PetscInt[], PetscViewer);
1734 PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt, const PetscReal[], PetscViewer);
1735 PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt, const PetscScalar[], PetscViewer);
1736 
1737 /*
1738     Basic memory and string operations. These are usually simple wrappers
1739    around the basic Unix system calls, but a few of them have additional
1740    functionality and/or error checking.
1741 */
1742 #include <petscstring.h>
1743 
1744 #include <stddef.h>
1745 #include <stdlib.h>
1746 
1747 #if defined(PETSC_CLANG_STATIC_ANALYZER)
1748   #define PetscPrefetchBlock(a, b, c, d)
1749 #else
1750   /*MC
1751    PetscPrefetchBlock - Prefetches a block of memory
1752 
1753    Synopsis:
1754     #include <petscsys.h>
1755     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)
1756 
1757    Not Collective
1758 
1759    Input Parameters:
1760 +  a  - pointer to first element to fetch (any type but usually `PetscInt` or `PetscScalar`)
1761 .  n  - number of elements to fetch
1762 .  rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors)
1763 -  t  - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note
1764 
1765    Level: developer
1766 
1767    Notes:
1768    The last two arguments (`rw` and `t`) must be compile-time constants.
1769 
1770    Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality.  Not all architectures offer
1771    equivalent locality hints, but the following macros are always defined to their closest analogue.
1772 +  `PETSC_PREFETCH_HINT_NTA` - Non-temporal.  Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched).
1773 .  `PETSC_PREFETCH_HINT_T0`  - Fetch to all levels of cache and evict to the closest level.  Use this when the memory will be reused regularly despite necessary eviction from L1.
1774 .  `PETSC_PREFETCH_HINT_T1`  - Fetch to level 2 and higher (not L1).
1775 -  `PETSC_PREFETCH_HINT_T2`  - Fetch to high-level cache only.  (On many systems, T0 and T1 are equivalent.)
1776 
1777    This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid
1778    address).
1779 
1780 M*/
1781   #define PetscPrefetchBlock(a, n, rw, t) \
1782     do { \
1783       const char *_p = (const char *)(a), *_end = (const char *)((a) + (n)); \
1784       for (; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p, (rw), (t)); \
1785     } while (0)
1786 #endif
1787 /*
1788       Determine if some of the kernel computation routines use
1789    Fortran (rather than C) for the numerical calculations. On some machines
1790    and compilers (like complex numbers) the Fortran version of the routines
1791    is faster than the C/C++ versions. The flag --with-fortran-kernels
1792    should be used with ./configure to turn these on.
1793 */
1794 #if defined(PETSC_USE_FORTRAN_KERNELS)
1795 
1796   #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1797     #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1798   #endif
1799 
1800   #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1801     #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
1802   #endif
1803 
1804   #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1805     #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
1806   #endif
1807 
1808   #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
1809     #define PETSC_USE_FORTRAN_KERNEL_MAXPY
1810   #endif
1811 
1812   #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1813     #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
1814   #endif
1815 
1816   #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
1817     #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
1818   #endif
1819 
1820   #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1821     #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
1822   #endif
1823 
1824   #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
1825     #define PETSC_USE_FORTRAN_KERNEL_MDOT
1826   #endif
1827 
1828   #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
1829     #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
1830   #endif
1831 
1832   #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
1833     #define PETSC_USE_FORTRAN_KERNEL_AYPX
1834   #endif
1835 
1836   #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
1837     #define PETSC_USE_FORTRAN_KERNEL_WAXPY
1838   #endif
1839 
1840 #endif
1841 
1842 /*
1843     Macros for indicating code that should be compiled with a C interface,
1844    rather than a C++ interface. Any routines that are dynamically loaded
1845    (such as the PCCreate_XXX() routines) must be wrapped so that the name
1846    mangler does not change the functions symbol name. This just hides the
1847    ugly extern "C" {} wrappers.
1848 */
1849 #if defined(__cplusplus)
1850   #define EXTERN_C_BEGIN extern "C" {
1851   #define EXTERN_C_END   }
1852 #else
1853   #define EXTERN_C_BEGIN
1854   #define EXTERN_C_END
1855 #endif
1856 
1857 /*MC
1858    MPI_Comm - the basic object used by MPI to determine which processes are involved in a
1859    communication
1860 
1861    Level: beginner
1862 
1863    Note:
1864    This manual page is a place-holder because MPICH does not have a manual page for `MPI_Comm`
1865 
1866 .seealso: `PETSC_COMM_WORLD`, `PETSC_COMM_SELF`
1867 M*/
1868 
1869 #if defined(PETSC_HAVE_MPIIO)
1870 PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4);
1871 PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4);
1872 PETSC_EXTERN PetscErrorCode MPIU_File_write_at(MPI_File, MPI_Offset, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(3, 5);
1873 PETSC_EXTERN PetscErrorCode MPIU_File_read_at(MPI_File, MPI_Offset, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(3, 5);
1874 PETSC_EXTERN PetscErrorCode MPIU_File_write_at_all(MPI_File, MPI_Offset, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(3, 5);
1875 PETSC_EXTERN PetscErrorCode MPIU_File_read_at_all(MPI_File, MPI_Offset, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(3, 5);
1876 #endif
1877 
1878 #if defined(PETSC_HAVE_MPI_COUNT)
1879 typedef MPI_Count MPIU_Count;
1880 #else
1881 typedef PetscInt64 MPIU_Count;
1882 #endif
1883 
1884 /*@C
1885    PetscIntCast - casts a `MPI_Count`, `PetscInt64`, `PetscCount`, or `size_t` to a `PetscInt` (which may be 32-bits in size), generates an
1886    error if the `PetscInt` is not large enough to hold the number.
1887 
1888    Not Collective; No Fortran Support
1889 
1890    Input Parameter:
1891 .  a - the `PetscInt64` value
1892 
1893    Output Parameter:
1894 .  b - the resulting `PetscInt` value, or `NULL` if the result is not needed
1895 
1896    Level: advanced
1897 
1898    Note:
1899    If integers needed for the applications are too large to fit in 32-bit ints you can ./configure using `--with-64-bit-indices` to make `PetscInt` use 64-bit integers
1900 
1901 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscCIntCast()`, `PetscIntMultError()`, `PetscIntSumError()`
1902 @*/
1903 static inline PetscErrorCode PetscIntCast(MPIU_Count a, PetscInt *b)
1904 {
1905   PetscFunctionBegin;
1906   if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */
1907   PetscCheck(sizeof(MPIU_Count) <= sizeof(PetscInt) || (a <= (MPIU_Count)PETSC_INT_MAX && a >= (MPIU_Count)PETSC_INT_MIN), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt64_FMT " is too big for PetscInt, you may need to ./configure using --with-64-bit-indices", (PetscInt64)a);
1908   if (b) *b = (PetscInt)a;
1909   PetscFunctionReturn(PETSC_SUCCESS);
1910 }
1911 
1912 /*@C
1913    PetscBLASIntCast - casts a `MPI_Count`, `PetscInt`, `PetscCount` or `PetscInt64` to a `PetscBLASInt` (which may be 32-bits in size), generates an
1914    error if the `PetscBLASInt` is not large enough to hold the number.
1915 
1916    Not Collective; No Fortran Support
1917 
1918    Input Parameter:
1919 .  a - the `PetscInt` value
1920 
1921    Output Parameter:
1922 .  b - the resulting `PetscBLASInt` value, or `NULL` if the result is not needed
1923 
1924    Level: advanced
1925 
1926    Note:
1927    Errors if the integer is negative since PETSc calls to BLAS/LAPACK never need to cast negative integer inputs
1928 
1929 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscCIntCast()`, `PetscIntCast()`
1930 @*/
1931 static inline PetscErrorCode PetscBLASIntCast(MPIU_Count a, PetscBLASInt *b)
1932 {
1933   PetscFunctionBegin;
1934   if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */
1935   PetscCheck(sizeof(MPIU_Count) <= sizeof(PetscBLASInt) || a <= (MPIU_Count)PETSC_BLAS_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt64_FMT " is too big for BLAS/LAPACK, which is restricted to 32-bit integers. Either you have an invalidly large integer error in your code or you must ./configure PETSc with --with-64-bit-blas-indices for the case you are running", (PetscInt64)a);
1936   PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer to BLAS/LAPACK routine");
1937   if (b) *b = (PetscBLASInt)a;
1938   PetscFunctionReturn(PETSC_SUCCESS);
1939 }
1940 
1941 /*@C
1942    PetscCuBLASIntCast - like `PetscBLASIntCast()`, but for `PetscCuBLASInt`.
1943 
1944    Not Collective; No Fortran Support
1945 
1946    Input Parameter:
1947 .  a - the `PetscInt` value
1948 
1949    Output Parameter:
1950 .  b - the resulting `PetscCuBLASInt` value, or `NULL` if the result is not needed
1951 
1952    Level: advanced
1953 
1954    Note:
1955    Errors if the integer is negative since PETSc calls to cuBLAS and friends never need to cast negative integer inputs
1956 
1957 .seealso: `PetscCuBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscCIntCast()`, `PetscIntCast()`
1958 @*/
1959 static inline PetscErrorCode PetscCuBLASIntCast(MPIU_Count a, PetscCuBLASInt *b)
1960 {
1961   PetscFunctionBegin;
1962   if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */
1963   PetscCheck(sizeof(MPIU_Count) <= sizeof(PetscCuBLASInt) || a <= (MPIU_Count)PETSC_CUBLAS_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt64_FMT " is too big for cuBLAS, which is restricted to 32-bit integers.", (PetscInt64)a);
1964   PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer %" PetscInt64_FMT "to cuBLAS routine", (PetscInt64)a);
1965   if (b) *b = (PetscCuBLASInt)a;
1966   PetscFunctionReturn(PETSC_SUCCESS);
1967 }
1968 
1969 /*@C
1970    PetscHipBLASIntCast - like `PetscBLASIntCast()`, but for `PetscHipBLASInt`.
1971 
1972    Not Collective; No Fortran Support
1973 
1974    Input Parameter:
1975 .  a - the `PetscInt` value
1976 
1977    Output Parameter:
1978 .  b - the resulting `PetscHipBLASInt` value, or `NULL` if the result is not needed
1979 
1980    Level: advanced
1981 
1982    Note:
1983    Errors if the integer is negative since PETSc calls to hipBLAS and friends never need to cast negative integer inputs
1984 
1985 .seealso: `PetscHipBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscCIntCast()`, `PetscIntCast()`
1986 @*/
1987 static inline PetscErrorCode PetscHipBLASIntCast(MPIU_Count a, PetscHipBLASInt *b)
1988 {
1989   PetscFunctionBegin;
1990   if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */
1991   PetscCheck(sizeof(MPIU_Count) <= sizeof(PetscHipBLASInt) || a <= (MPIU_Count)PETSC_HIPBLAS_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt64_FMT " is too big for hipBLAS, which is restricted to 32-bit integers.", (PetscInt64)a);
1992   PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer %" PetscInt64_FMT "to hipBLAS routine", (PetscInt64)a);
1993   if (b) *b = (PetscHipBLASInt)a;
1994   PetscFunctionReturn(PETSC_SUCCESS);
1995 }
1996 
1997 /*@C
1998    PetscMPIIntCast - casts a `MPI_Count`, `PetscInt`, `PetscCount`, or `PetscInt64` to a `PetscMPIInt` (which is always 32-bits in size), generates an
1999    error if the `PetscMPIInt` is not large enough to hold the number.
2000 
2001    Not Collective; No Fortran Support
2002 
2003    Input Parameter:
2004 .  a - the `PetscInt` value
2005 
2006    Output Parameter:
2007 .  b - the resulting `PetscMPIInt` value, or `NULL` if the result is not needed
2008 
2009    Level: advanced
2010 
2011 .seealso: [](stylePetscCount), `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscIntCast()`
2012 @*/
2013 static inline PetscErrorCode PetscMPIIntCast(MPIU_Count a, PetscMPIInt *b)
2014 {
2015   PetscFunctionBegin;
2016   if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */
2017   PetscCheck(a <= (MPIU_Count)PETSC_MPI_INT_MAX && a >= (MPIU_Count)PETSC_MPI_INT_MIN, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt64_FMT " is too big for MPI buffer length. Maximum supported value is %d", (PetscInt64)a, PETSC_MPI_INT_MAX);
2018   if (b) *b = (PetscMPIInt)a;
2019   PetscFunctionReturn(PETSC_SUCCESS);
2020 }
2021 
2022 /*@C
2023    PetscCIntCast - casts a `MPI_Count`, `PetscInt`, `PetscCount`, or `PetscInt64` to a `int`, generates an error if the `int` is not large enough to hold the number.
2024 
2025    Not Collective; No Fortran Support
2026 
2027    Input Parameter:
2028 .  a - the `PetscInt` value
2029 
2030    Output Parameter:
2031 .  b - the resulting `int` value, or `NULL` if the result is not needed
2032 
2033    Level: advanced
2034 
2035 .seealso: [](stylePetscCount), `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscIntCast()`
2036 @*/
2037 static inline PetscErrorCode PetscCIntCast(MPIU_Count a, int *b)
2038 {
2039   PetscFunctionBegin;
2040   if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */
2041   PetscCheck(a <= INT_MAX && a >= INT_MIN, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt64_FMT " is too big to be casted to an int. Maximum supported value is %d", (PetscInt64)a, INT_MAX);
2042   if (b) *b = (int)a;
2043   PetscFunctionReturn(PETSC_SUCCESS);
2044 }
2045 
2046 /*MC
2047    PetscInt64Mult - Computes the product of two variables after casting them to `PetscInt64`.
2048 
2049    Not Collective; No Fortran Support
2050 
2051    Input Parameters:
2052 +  a - the first variable
2053 -  b - the second variable
2054 
2055    Level: advanced
2056 
2057 .seealso: [](stylePetscCount), `PetscIntMultError()`, `PetscIntMultTruncate()`
2058 M*/
2059 #if defined(PETSC_USE_64BIT_INDICES)
2060   #define PetscInt64Mult(a, b) ((a) * (b))
2061 #else
2062   #define PetscInt64Mult(a, b) (((PetscInt64)(a)) * ((PetscInt64)(b)))
2063 #endif
2064 
2065 /*@C
2066   PetscRealIntMultTruncate - Computes the product of a positive `PetscReal` and a positive
2067   `PetscInt` and truncates the value to slightly less than the maximal possible value.
2068 
2069   Not Collective; No Fortran Support
2070 
2071   Input Parameters:
2072 + a - The `PetscReal` value
2073 - b - The `PetscInt` value
2074 
2075   Level: advanced
2076 
2077   Notes:
2078   Returns the result as a `PetscInt` value.
2079 
2080   Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`.
2081 
2082   Use `PetscIntMultTruncate()` to compute the product of two positive `PetscInt` and truncate
2083   to fit a `PetscInt`.
2084 
2085   Use `PetscIntMultError()` to compute the product of two `PetscInt` if you wish to generate an
2086   error if the result will not fit in a `PetscInt`.
2087 
2088   Developer Notes:
2089   We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but
2090   requires many more checks.
2091 
2092   This is used where we compute approximate sizes for workspace and need to insure the
2093   workspace is index-able.
2094 
2095 .seealso: `PetscReal`, `PetscInt`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()`
2096 @*/
2097 static inline PetscInt PetscRealIntMultTruncate(PetscReal a, PetscInt b)
2098 {
2099   PetscInt64 r = (PetscInt64)(a * (PetscReal)b);
2100   if (r > PETSC_INT_MAX - 100) r = PETSC_INT_MAX - 100;
2101 #if defined(PETSC_USE_64BIT_INDICES)
2102   return r;
2103 #else
2104   return (PetscInt)r;
2105 #endif
2106 }
2107 
2108 /*@C
2109    PetscIntMultTruncate - Computes the product of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value
2110 
2111    Not Collective; No Fortran Support
2112 
2113    Input Parameters:
2114 +  a - the `PetscInt` value
2115 -  b - the second value
2116 
2117    Returns:
2118    The result as a `PetscInt` value
2119 
2120    Level: advanced
2121 
2122    Notes:
2123    Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`
2124 
2125    Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt`
2126 
2127    Use `PetscIntMultError()` to compute the product of two `PetscInt` if you wish to generate an error if the result will not fit in a `PetscInt`
2128 
2129    Developer Notes:
2130    We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but requires many more checks.
2131 
2132    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2133 
2134 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()`,
2135           `PetscIntSumTruncate()`
2136 @*/
2137 static inline PetscInt PetscIntMultTruncate(PetscInt a, PetscInt b)
2138 {
2139   PetscInt64 r = PetscInt64Mult(a, b);
2140   if (r > PETSC_INT_MAX - 100) r = PETSC_INT_MAX - 100;
2141 #if defined(PETSC_USE_64BIT_INDICES)
2142   return r;
2143 #else
2144   return (PetscInt)r;
2145 #endif
2146 }
2147 
2148 /*@C
2149    PetscIntSumTruncate - Computes the sum of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value
2150 
2151    Not Collective; No Fortran Support
2152 
2153    Input Parameters:
2154 +  a - the `PetscInt` value
2155 -  b - the second value
2156 
2157    Returns:
2158    The result as a `PetscInt` value
2159 
2160    Level: advanced
2161 
2162    Notes:
2163    Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`
2164 
2165    Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt`
2166 
2167    Use `PetscIntMultError()` to compute the product of two `PetscInt` if you wish to generate an error if the result will not fit in a `PetscInt`
2168 
2169    Developer Note:
2170    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2171 
2172 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`
2173 @*/
2174 static inline PetscInt PetscIntSumTruncate(PetscInt a, PetscInt b)
2175 {
2176   PetscInt64 r = a;
2177 
2178   r += b;
2179   if (r > PETSC_INT_MAX - 100) r = PETSC_INT_MAX - 100;
2180 #if defined(PETSC_USE_64BIT_INDICES)
2181   return r;
2182 #else
2183   return (PetscInt)r;
2184 #endif
2185 }
2186 
2187 /*@C
2188    PetscIntMultError - Computes the product of two positive `PetscInt` and generates an error with overflow.
2189 
2190    Not Collective; No Fortran Support
2191 
2192    Input Parameters:
2193 +  a - the `PetscInt` value
2194 -  b - the second value
2195 
2196    Output Parameter:
2197 .  result - the result as a `PetscInt` value, or `NULL` if you do not want the result, you just want to check if it overflows
2198 
2199    Level: advanced
2200 
2201    Notes:
2202    Use `PetscInt64Mult()` to compute the product of two `PetscInt` and store in a `PetscInt64`
2203 
2204    Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt`
2205 
2206    Developer Note:
2207    In most places in the source code we currently assume that `PetscInt` addition does not overflow, this is obviously wrong but requires many more checks.
2208    `PetscIntSumError()` can be used to check for this situation.
2209 
2210 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntSumError()`
2211 @*/
2212 static inline PetscErrorCode PetscIntMultError(PetscInt a, PetscInt b, PetscInt *result)
2213 {
2214   PetscInt64 r = PetscInt64Mult(a, b);
2215 
2216   PetscFunctionBegin;
2217 #if defined(PETSC_USE_64BIT_INDICES)
2218   if (result) *result = r;
2219 #else
2220   if (result) *result = (PetscInt)r;
2221 #endif
2222   if (!PetscDefined(USE_64BIT_INDICES)) {
2223     PetscCheck(r <= PETSC_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_SUP, "Product of two integers %" PetscInt_FMT " %" PetscInt_FMT " overflow, either you have an invalidly large integer error in your code or you must ./configure PETSc with --with-64-bit-indices for the case you are running", a, b);
2224   }
2225   PetscFunctionReturn(PETSC_SUCCESS);
2226 }
2227 
2228 /*@C
2229 
2230    PetscIntSumError - Computes the sum of two positive `PetscInt` and generates an error with overflow.
2231 
2232    Not Collective; No Fortran Support
2233 
2234    Input Parameters:
2235 +  a - the `PetscInt` value
2236 -  b - the second value
2237 
2238    Output Parameter:
2239 .  c - the result as a `PetscInt` value,  or `NULL` if you do not want the result, you just want to check if it overflows
2240 
2241    Level: advanced
2242 
2243    Notes:
2244    Use `PetscInt64Mult()` to compute the product of two 32-bit `PetscInt` and store in a `PetscInt64`
2245 
2246    Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt`
2247 
2248 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`
2249 @*/
2250 static inline PetscErrorCode PetscIntSumError(PetscInt a, PetscInt b, PetscInt *result)
2251 {
2252   PetscInt64 r = a;
2253 
2254   PetscFunctionBegin;
2255   r += b;
2256 #if defined(PETSC_USE_64BIT_INDICES)
2257   if (result) *result = r;
2258 #else
2259   if (result) *result = (PetscInt)r;
2260 #endif
2261   if (!PetscDefined(USE_64BIT_INDICES)) {
2262     PetscCheck(r <= PETSC_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_SUP, "Sum of two integers %" PetscInt_FMT " %" PetscInt_FMT " overflow, either you have an invalidly large integer error in your code or you must ./configure PETSc with --with-64-bit-indices for the case you are running", a, b);
2263   }
2264   PetscFunctionReturn(PETSC_SUCCESS);
2265 }
2266 
2267 /*
2268      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2269 */
2270 #if defined(hz)
2271   #undef hz
2272 #endif
2273 
2274 #if defined(PETSC_HAVE_SYS_TYPES_H)
2275   #include <sys/types.h>
2276 #endif
2277 
2278 /*MC
2279 
2280     PETSC_VERSION - This manual page provides information about how PETSc documents and uses its version information. This information is available to both C/C++
2281                     and Fortran compilers when `petscsys.h` is included.
2282 
2283     The current PETSc version and the API for accessing it are defined in <A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/include/petscversion.h.html">include/petscversion.html</A>
2284 
2285     The complete version number is given as the triple  PETSC_VERSION_MAJOR.PETSC_VERSION_MINOR.PETSC_VERSION_SUBMINOR (in short hand x.y.z)
2286 
2287     A change in the minor version number (y) indicates possible/likely changes in the PETSc API. Note this is different than with the semantic versioning convention
2288     where only a change in the major version number (x) indicates a change in the API.
2289 
2290     A subminor greater than zero indicates a patch release. Version x.y.z maintains source and binary compatibility with version x.y.w for all z and w
2291 
2292     Use the macros PETSC_VERSION_EQ(x,y,z), PETSC_VERSION_LT(x,y,z), PETSC_VERSION_LE(x,y,z), PETSC_VERSION_GT(x,y,z),
2293     PETSC_VERSION_GE(x,y,z) to determine if the current version is equal to, less than, less than or equal to, greater than or greater than or equal to a given
2294     version number (x.y.z).
2295 
2296     `PETSC_RELEASE_DATE` is the date the x.y version was released (i.e. the version before any patch releases)
2297 
2298     `PETSC_VERSION_DATE` is the date the x.y.z version was released
2299 
2300     `PETSC_VERSION_GIT` is the last git commit to the repository given in the form vx.y.z-wwwww
2301 
2302     `PETSC_VERSION_DATE_GIT` is the date of the last git commit to the repository
2303 
2304     `PETSC_VERSION_()` is deprecated and will eventually be removed.
2305 
2306     Level: intermediate
2307 M*/
2308 
2309 PETSC_EXTERN PetscErrorCode PetscGetArchType(char[], size_t);
2310 PETSC_EXTERN PetscErrorCode PetscGetHostName(char[], size_t);
2311 PETSC_EXTERN PetscErrorCode PetscGetUserName(char[], size_t);
2312 PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[], size_t);
2313 PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2314 PETSC_EXTERN PetscErrorCode PetscGetDate(char[], size_t);
2315 PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2316 PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt *, PetscInt *, PetscInt *, PetscInt *);
2317 
2318 PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscCount, const PetscInt[], PetscBool *);
2319 PETSC_EXTERN PetscErrorCode PetscSortedInt64(PetscCount, const PetscInt64[], PetscBool *);
2320 PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscCount, const PetscMPIInt[], PetscBool *);
2321 PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscCount, const PetscReal[], PetscBool *);
2322 PETSC_EXTERN PetscErrorCode PetscSortInt(PetscCount, PetscInt[]);
2323 PETSC_EXTERN PetscErrorCode PetscSortInt64(PetscCount, PetscInt64[]);
2324 PETSC_EXTERN PetscErrorCode PetscSortCount(PetscCount, PetscCount[]);
2325 PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscCount, PetscInt[]);
2326 PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *, PetscInt[]);
2327 PETSC_EXTERN PetscErrorCode PetscSortedCheckDupsInt(PetscCount, const PetscInt[], PetscBool *);
2328 PETSC_EXTERN PetscErrorCode PetscSortedCheckDupsCount(PetscCount, const PetscCount[], PetscBool *);
2329 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt *, PetscInt[]);
2330 PETSC_EXTERN PetscErrorCode PetscCheckDupsInt(PetscInt, const PetscInt[], PetscBool *);
2331 PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscCount, const PetscInt[], PetscInt *);
2332 PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscCount, const PetscMPIInt[], PetscInt *);
2333 PETSC_EXTERN PetscErrorCode PetscFindCount(PetscCount, PetscCount, const PetscCount[], PetscCount *);
2334 PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt, const PetscInt[], PetscInt[]);
2335 PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt, const char *[], PetscInt[]);
2336 PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscCount, PetscInt[], PetscInt[]);
2337 PETSC_EXTERN PetscErrorCode PetscSortIntWithCountArray(PetscCount, PetscInt[], PetscCount[]);
2338 PETSC_EXTERN PetscErrorCode PetscSortIntWithMPIIntArray(PetscCount, PetscInt[], PetscMPIInt[]);
2339 PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscCount, PetscInt[], PetscInt[], PetscInt[]);
2340 PETSC_EXTERN PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount, PetscInt[], PetscInt[], PetscCount[]);
2341 PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscCount, PetscMPIInt[]);
2342 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *, PetscMPIInt[]);
2343 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscCount, PetscMPIInt[], PetscMPIInt[]);
2344 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscCount, PetscMPIInt[], PetscInt[]);
2345 PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscCount, PetscInt[], PetscScalar[]);
2346 PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscCount, PetscInt[], void *, size_t, void *);
2347 PETSC_EXTERN PetscErrorCode PetscSortReal(PetscCount, PetscReal[]);
2348 PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscCount, PetscReal[], PetscInt[]);
2349 PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt, const PetscReal[], PetscInt[]);
2350 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt *, PetscReal[]);
2351 PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal, PetscCount, const PetscReal[], PetscReal, PetscInt *);
2352 PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt, PetscInt, PetscScalar[], PetscInt[]);
2353 PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt, PetscInt, PetscReal[], PetscInt[]);
2354 PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt, const PetscBool[], const PetscInt[], PetscInt *, PetscInt *[], PetscInt *[], PetscInt *[], PetscInt *[]);
2355 PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt, const PetscInt[], const PetscInt[], PetscInt, const PetscInt[], const PetscInt[], PetscInt *, PetscInt *[], PetscInt *[]);
2356 PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscInt *, PetscInt *[]);
2357 PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt, const PetscMPIInt[], PetscInt, const PetscMPIInt[], PetscInt *, PetscMPIInt *[]);
2358 PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *);
2359 
2360 PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt, void *, size_t, int (*)(const void *, const void *, void *), void *);
2361 PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt, PetscInt[]);
2362 PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt, PetscMPIInt[]);
2363 PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt, PetscReal[]);
2364 PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt, void *, size_t, void *, size_t, int (*)(const void *, const void *, void *), void *);
2365 PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt, PetscInt[], PetscInt[]);
2366 PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt, PetscMPIInt[], PetscMPIInt[]);
2367 PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt, PetscReal[], PetscInt[]);
2368 
2369 PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2370 PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[], size_t);
2371 
2372 /*J
2373     PetscRandomType - String with the name of a PETSc randomizer
2374 
2375    Level: beginner
2376 
2377    Note:
2378    To use `PETSCSPRNG` or `PETSCRANDOM123` you must have ./configure PETSc
2379    with the option `--download-sprng` or `--download-random123`. We recommend the default provided with PETSc.
2380 
2381 .seealso: `PetscRandomSetType()`, `PetscRandom`, `PetscRandomCreate()`
2382 J*/
2383 typedef const char *PetscRandomType;
2384 #define PETSCRAND      "rand"
2385 #define PETSCRAND48    "rand48"
2386 #define PETSCSPRNG     "sprng"
2387 #define PETSCRANDER48  "rander48"
2388 #define PETSCRANDOM123 "random123"
2389 #define PETSCCURAND    "curand"
2390 
2391 /* Logging support */
2392 PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;
2393 
2394 PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);
2395 PETSC_EXTERN PetscErrorCode PetscRandomFinalizePackage(void);
2396 
2397 /* Dynamic creation and loading functions */
2398 PETSC_EXTERN PetscFunctionList PetscRandomList;
2399 
2400 PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[], PetscErrorCode (*)(PetscRandom));
2401 PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2402 PETSC_EXTERN PetscErrorCode PetscRandomSetOptionsPrefix(PetscRandom, const char[]);
2403 PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2404 PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType *);
2405 PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom, PetscObject, const char[]);
2406 PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom, PetscViewer);
2407 
2408 PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm, PetscRandom *);
2409 PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom, PetscScalar *);
2410 PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom, PetscReal *);
2411 PETSC_EXTERN PetscErrorCode PetscRandomGetValues(PetscRandom, PetscInt, PetscScalar *);
2412 PETSC_EXTERN PetscErrorCode PetscRandomGetValuesReal(PetscRandom, PetscInt, PetscReal *);
2413 PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom, PetscScalar *, PetscScalar *);
2414 PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom, PetscScalar, PetscScalar);
2415 PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom, PetscInt64);
2416 PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom, PetscInt64 *);
2417 PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2418 PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom *);
2419 
2420 PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[], char[], size_t);
2421 PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[], char[], size_t);
2422 PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[], size_t);
2423 PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[], char[]);
2424 PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[], size_t);
2425 PETSC_EXTERN PetscErrorCode PetscTestFile(const char[], char, PetscBool *);
2426 PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[], char, PetscBool *);
2427 PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2428 PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2429 PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);
2430 
2431 /*MC
2432    PetscBinaryBigEndian - indicates if values in memory are stored with big endian format
2433 
2434    Synopsis:
2435    #include <petscsys.h>
2436    PetscBool PetscBinaryBigEndian(void);
2437 
2438    No Fortran Support
2439 
2440    Level: developer
2441 
2442 .seealso: `PetscInitialize()`, `PetscFinalize()`, `PetscInitializeCalled`
2443 M*/
2444 static inline PetscBool PetscBinaryBigEndian(void)
2445 {
2446   long _petsc_v = 1;
2447   return ((char *)&_petsc_v)[0] ? PETSC_FALSE : PETSC_TRUE;
2448 }
2449 
2450 PETSC_EXTERN PetscErrorCode PetscBinaryRead(int, void *, PetscCount, PetscInt *, PetscDataType);
2451 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm, int, void *, PetscInt, PetscInt *, PetscDataType);
2452 PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int, const void *, PetscCount, PetscDataType);
2453 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm, int, const void *, PetscInt, PetscDataType);
2454 PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[], PetscFileMode, int *);
2455 PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2456 PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm, PetscBool *);
2457 PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm, PetscBool *);
2458 PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm, char[], size_t);
2459 PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm, const char[], char[], size_t, PetscBool *);
2460 PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm, const char[], char[], size_t, PetscBool *);
2461 #if defined(PETSC_USE_SOCKET_VIEWER)
2462 PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[], int, int *);
2463 #endif
2464 
2465 PETSC_EXTERN PetscErrorCode PetscBinarySeek(int, off_t, PetscBinarySeekType, off_t *);
2466 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm, int, off_t, PetscBinarySeekType, off_t *);
2467 PETSC_EXTERN PetscErrorCode PetscByteSwap(void *, PetscDataType, PetscCount);
2468 
2469 PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2470 PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[], PetscBool);
2471 PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2472 PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char *);
2473 PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2474 PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2475 PETSC_EXTERN PetscErrorCode PetscWaitOnError(void);
2476 
2477 PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt *);
2478 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], PetscMPIInt *[], PetscMPIInt *[]);
2479 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt *[], PetscMPIInt *[], PetscMPIInt *[]);
2480 PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscInt ***, MPI_Request **);
2481 PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscScalar ***, MPI_Request **);
2482 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm, PetscMPIInt, MPI_Datatype, PetscMPIInt, const PetscMPIInt[], const void *, PetscMPIInt *, PetscMPIInt *[], void *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(6, 3);
2483 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm, PetscMPIInt, MPI_Datatype, PetscMPIInt, const PetscMPIInt[], const void *, PetscMPIInt *, PetscMPIInt **, void *, PetscMPIInt, PetscErrorCode (*send)(MPI_Comm, const PetscMPIInt[], PetscMPIInt, PetscMPIInt, void *, MPI_Request[], void *), PetscErrorCode (*recv)(MPI_Comm, const PetscMPIInt[], PetscMPIInt, void *, MPI_Request[], void *), void *ctx) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(6, 3);
2484 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm, PetscMPIInt, MPI_Datatype, PetscMPIInt, const PetscMPIInt[], const void *, PetscMPIInt *, PetscMPIInt **, void *, PetscMPIInt, MPI_Request **, MPI_Request **, PetscErrorCode (*send)(MPI_Comm, const PetscMPIInt[], PetscMPIInt, PetscMPIInt, void *, MPI_Request[], void *), PetscErrorCode (*recv)(MPI_Comm, const PetscMPIInt[], PetscMPIInt, void *, MPI_Request[], void *), void *ctx) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(6, 3);
2485 
2486 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm, PetscBuildTwoSidedType);
2487 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm, PetscBuildTwoSidedType *);
2488 
2489 PETSC_DEPRECATED_FUNCTION(3, 24, 0, "PetscSSEIsEnabled()", ) static inline PetscErrorCode PetscSSEIsEnabled(PETSC_UNUSED MPI_Comm comm, PetscBool *lflag, PetscBool *gflag)
2490 {
2491   if (lflag) *lflag = PETSC_FALSE;
2492   if (gflag) *gflag = PETSC_FALSE;
2493   return PETSC_SUCCESS;
2494 }
2495 
2496 PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);
2497 
2498 struct _n_PetscSubcomm {
2499   MPI_Comm         parent;    /* parent communicator */
2500   MPI_Comm         dupparent; /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2501   MPI_Comm         child;     /* the sub-communicator */
2502   PetscMPIInt      n;         /* num of subcommunicators under the parent communicator */
2503   PetscMPIInt      color;     /* color of processors belong to this communicator */
2504   PetscMPIInt     *subsize;   /* size of subcommunicator[color] */
2505   PetscSubcommType type;
2506   char            *subcommprefix;
2507 };
2508 
2509 static inline MPI_Comm PetscSubcommParent(PetscSubcomm scomm)
2510 {
2511   return scomm->parent;
2512 }
2513 static inline MPI_Comm PetscSubcommChild(PetscSubcomm scomm)
2514 {
2515   return scomm->child;
2516 }
2517 static inline MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm)
2518 {
2519   return scomm->dupparent;
2520 }
2521 PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm, PetscSubcomm *);
2522 PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm *);
2523 PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm, PetscInt);
2524 PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm, PetscSubcommType);
2525 PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm, PetscMPIInt, PetscMPIInt);
2526 PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm, PetscViewer);
2527 PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2528 PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm, const char[]);
2529 PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm, MPI_Comm *);
2530 PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm, MPI_Comm *);
2531 PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm, MPI_Comm *);
2532 
2533 PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt, PetscHeap *);
2534 PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap, PetscInt, PetscInt);
2535 PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap, PetscInt *, PetscInt *);
2536 PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap, PetscInt *, PetscInt *);
2537 PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap, PetscInt, PetscInt);
2538 PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2539 PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap *);
2540 PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap, PetscViewer);
2541 
2542 PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2543 PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm, PetscShmComm *);
2544 PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm, PetscMPIInt, PetscMPIInt *);
2545 PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm, PetscMPIInt, PetscMPIInt *);
2546 PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm, MPI_Comm *);
2547 
2548 /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2549 PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm, PetscInt, PetscOmpCtrl *);
2550 PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl, MPI_Comm *, MPI_Comm *, PetscBool *);
2551 PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl *);
2552 PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2553 PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2554 PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);
2555 
2556 PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t, PetscCount, PetscSegBuffer *);
2557 PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer *);
2558 PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer, PetscCount, void *);
2559 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer, void *);
2560 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer, void *);
2561 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer, void *);
2562 PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer, PetscCount *);
2563 PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer, PetscCount);
2564 
2565 /*MC
2566   PetscSegBufferGetInts - access an array of `PetscInt` from a `PetscSegBuffer`
2567 
2568   Synopsis:
2569   #include <petscsys.h>
2570   PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg, size_t count, PetscInt *PETSC_RESTRICT *slot);
2571 
2572   No Fortran Support
2573 
2574   Input Parameters:
2575 + seg   - `PetscSegBuffer` buffer
2576 - count - number of entries needed
2577 
2578   Output Parameter:
2579 . buf - address of new buffer for contiguous data
2580 
2581   Level: intermediate
2582 
2583   Developer Note:
2584   Type-safe wrapper to encourage use of PETSC_RESTRICT. Does not use PetscFunctionBegin because the error handling
2585   prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as
2586   possible.
2587 
2588 .seealso: `PetscSegBuffer`, `PetscSegBufferGet()`, `PetscInitialize()`, `PetscFinalize()`, `PetscInitializeCalled`
2589 M*/
2590 static inline PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg, PetscCount count, PetscInt *PETSC_RESTRICT *slot)
2591 {
2592   return PetscSegBufferGet(seg, count, (void **)slot);
2593 }
2594 
2595 extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2596 PETSC_EXTERN PetscErrorCode    PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted *);
2597 PETSC_EXTERN PetscErrorCode    PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted *);
2598 PETSC_EXTERN PetscErrorCode    PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted, const char *, const char *, PetscBool *);
2599 
2600 #include <stdarg.h>
2601 PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char *, size_t, const char[], size_t *, va_list);
2602 PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE *, const char[], va_list);
2603 
2604 PETSC_EXTERN PetscSegBuffer PetscCitationsList;
2605 
2606 /*@
2607      PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code.
2608 
2609      Not Collective; No Fortran Support
2610 
2611      Input Parameters:
2612 +    cite - the bibtex item, formatted to displayed on multiple lines nicely
2613 -    set - a boolean variable initially set to `PETSC_FALSE`; this is used to insure only a single registration of the citation
2614 
2615      Options Database Key:
2616 .     -citations [filename]   - print out the bibtex entries for the given computation
2617 
2618      Level: intermediate
2619 @*/
2620 static inline PetscErrorCode PetscCitationsRegister(const char cit[], PetscBool *set)
2621 {
2622   size_t len;
2623   char  *vstring;
2624 
2625   PetscFunctionBegin;
2626   if (set && *set) PetscFunctionReturn(PETSC_SUCCESS);
2627   PetscCall(PetscStrlen(cit, &len));
2628   PetscCall(PetscSegBufferGet(PetscCitationsList, (PetscCount)len, &vstring));
2629   PetscCall(PetscArraycpy(vstring, cit, len));
2630   if (set) *set = PETSC_TRUE;
2631   PetscFunctionReturn(PETSC_SUCCESS);
2632 }
2633 
2634 PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm, char[], char[], size_t);
2635 PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm, const char[], char[], size_t);
2636 PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm, const char[], const char[]);
2637 
2638 PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm, char[], char[], size_t);
2639 PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm, const char[], char[], char[], size_t);
2640 PETSC_EXTERN PetscErrorCode PetscBoxUpload(MPI_Comm, const char[], const char[]);
2641 
2642 PETSC_EXTERN PetscErrorCode PetscGlobusGetTransfers(MPI_Comm, const char[], char[], size_t);
2643 PETSC_EXTERN PetscErrorCode PetscGlobusAuthorize(MPI_Comm, char[], size_t);
2644 PETSC_EXTERN PetscErrorCode PetscGlobusUpload(MPI_Comm, const char[], const char[]);
2645 
2646 PETSC_EXTERN PetscErrorCode PetscPullJSONValue(const char[], const char[], char[], size_t, PetscBool *);
2647 PETSC_EXTERN PetscErrorCode PetscPushJSONValue(char[], const char[], const char[], size_t);
2648 
2649 #if !defined(PETSC_HAVE_MPI_LARGE_COUNT)
2650   /*
2651    Cast PetscCount <a> to PetscMPIInt <b>, where <a> is likely used for the 'count' argument in MPI routines.
2652    It is similar to PetscMPIIntCast() except that here it returns an MPI error code.
2653 */
2654   #define PetscMPIIntCast_Internal(a, b) \
2655     do { \
2656       *b = 0; \
2657       if (PetscUnlikely(a > (MPIU_Count)PETSC_MPI_INT_MAX)) return MPI_ERR_COUNT; \
2658       *b = (PetscMPIInt)a; \
2659     } while (0)
2660 
2661 static inline PetscMPIInt MPIU_Get_count(MPI_Status *status, MPI_Datatype dtype, PetscCount *count)
2662 {
2663   PetscMPIInt count2, err;
2664 
2665   *count = 0; /* to prevent incorrect warnings of uninitialized variables */
2666   err    = MPI_Get_count(status, dtype, &count2);
2667   *count = count2;
2668   return err;
2669 }
2670 
2671 static inline PetscMPIInt MPIU_Send(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt dest, PetscMPIInt tag, MPI_Comm comm)
2672 {
2673   PetscMPIInt count2, err;
2674 
2675   PetscMPIIntCast_Internal(count, &count2);
2676   err = MPI_Send((void *)buf, count2, dtype, dest, tag, comm);
2677   return err;
2678 }
2679 
2680 static inline PetscMPIInt MPIU_Send_init(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt dest, PetscMPIInt tag, MPI_Comm comm, MPI_Request *request)
2681 {
2682   PetscMPIInt count2, err;
2683 
2684   PetscMPIIntCast_Internal(count, &count2);
2685   err = MPI_Send_init((void *)buf, count2, dtype, dest, tag, comm, request);
2686   return err;
2687 }
2688 
2689 static inline PetscMPIInt MPIU_Isend(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt dest, PetscMPIInt tag, MPI_Comm comm, MPI_Request *request)
2690 {
2691   PetscMPIInt count2, err;
2692 
2693   PetscMPIIntCast_Internal(count, &count2);
2694   err = MPI_Isend((void *)buf, count2, dtype, dest, tag, comm, request);
2695   return err;
2696 }
2697 
2698 static inline PetscMPIInt MPIU_Recv(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt source, PetscMPIInt tag, MPI_Comm comm, MPI_Status *status)
2699 {
2700   PetscMPIInt count2, err;
2701 
2702   PetscMPIIntCast_Internal(count, &count2);
2703   err = MPI_Recv((void *)buf, count2, dtype, source, tag, comm, status);
2704   return err;
2705 }
2706 
2707 static inline PetscMPIInt MPIU_Recv_init(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt source, PetscMPIInt tag, MPI_Comm comm, MPI_Request *request)
2708 {
2709   PetscMPIInt count2, err;
2710 
2711   PetscMPIIntCast_Internal(count, &count2);
2712   err = MPI_Recv_init((void *)buf, count2, dtype, source, tag, comm, request);
2713   return err;
2714 }
2715 
2716 static inline PetscMPIInt MPIU_Irecv(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt source, PetscMPIInt tag, MPI_Comm comm, MPI_Request *request)
2717 {
2718   PetscMPIInt count2, err;
2719 
2720   PetscMPIIntCast_Internal(count, &count2);
2721   err = MPI_Irecv((void *)buf, count2, dtype, source, tag, comm, request);
2722   return err;
2723 }
2724 
2725 static inline PetscMPIInt MPIU_Reduce(const void *inbuf, void *outbuf, MPIU_Count count, MPI_Datatype dtype, MPI_Op op, PetscMPIInt root, MPI_Comm comm)
2726 {
2727   PetscMPIInt count2, err;
2728 
2729   PetscMPIIntCast_Internal(count, &count2);
2730   err = MPI_Reduce((void *)inbuf, outbuf, count2, dtype, op, root, comm);
2731   return err;
2732 }
2733 
2734   #if defined(PETSC_HAVE_MPI_REDUCE_LOCAL)
2735 static inline PetscMPIInt MPIU_Reduce_local(const void *inbuf, void *inoutbuf, MPIU_Count count, MPI_Datatype dtype, MPI_Op op)
2736 {
2737   PetscMPIInt count2, err;
2738 
2739   PetscMPIIntCast_Internal(count, &count2);
2740   err = MPI_Reduce_local((void *)inbuf, inoutbuf, count2, dtype, op);
2741   return err;
2742 }
2743   #endif
2744 
2745   #if !defined(PETSC_USE_64BIT_INDICES)
2746     #define MPIU_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm) MPI_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm)
2747     #define MPIU_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm)  MPI_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm)
2748   #else
2749     #define MPIU_Scatterv(sendbuf, sendcount, displs, sendtype, recvbuf, recvcount, recvtype, root, comm) \
2750       ((void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_SUP, PETSC_ERROR_INITIAL, "Must have MPI 4 support for MPI_Scatterv_c() for this functionality, upgrade your MPI"), MPI_ERR_COUNT)
2751     #define MPIU_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm) \
2752       ((void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_SUP, PETSC_ERROR_INITIAL, "Must have MPI 4 support for MPI_Scatterv_c() for this functionality, upgrade your MPI"), MPI_ERR_COUNT)
2753   #endif
2754 
2755 #else
2756 
2757   /* on 32 bit systems MPI_Count maybe 64-bit while PetscCount is 32-bit */
2758   #define PetscCountCast_Internal(a, b) \
2759     do { \
2760       *b = 0; \
2761       if (PetscUnlikely(a > (MPI_Count)PETSC_COUNT_MAX)) return MPI_ERR_COUNT; \
2762       *b = (PetscMPIInt)a; \
2763     } while (0)
2764 
2765 static inline PetscMPIInt MPIU_Get_count(MPI_Status *status, MPI_Datatype dtype, PetscCount *count)
2766 {
2767   MPI_Count   count2;
2768   PetscMPIInt err;
2769 
2770   *count = 0; /* to prevent incorrect warnings of uninitialized variables */
2771   err    = MPI_Get_count_c(status, dtype, &count2);
2772   if (err) return err;
2773   PetscCountCast_Internal(count2, count);
2774   return MPI_SUCCESS;
2775 }
2776 
2777   #define MPIU_Reduce(inbuf, outbuf, count, dtype, op, root, comm)      MPI_Reduce_c(inbuf, outbuf, (MPI_Count)(count), dtype, op, root, comm)
2778   #define MPIU_Send(buf, count, dtype, dest, tag, comm)                 MPI_Send_c(buf, (MPI_Count)(count), dtype, dest, tag, comm)
2779   #define MPIU_Send_init(buf, count, dtype, dest, tag, comm, request)   MPI_Send_init_c(buf, (MPI_Count)(count), dtype, dest, tag, comm, request)
2780   #define MPIU_Isend(buf, count, dtype, dest, tag, comm, request)       MPI_Isend_c(buf, (MPI_Count)(count), dtype, dest, tag, comm, request)
2781   #define MPIU_Recv(buf, count, dtype, source, tag, comm, status)       MPI_Recv_c(buf, (MPI_Count)(count), dtype, source, tag, comm, status)
2782   #define MPIU_Recv_init(buf, count, dtype, source, tag, comm, request) MPI_Recv_init_c(buf, (MPI_Count)(count), dtype, source, tag, comm, request)
2783   #define MPIU_Irecv(buf, count, dtype, source, tag, comm, request)     MPI_Irecv_c(buf, (MPI_Count)(count), dtype, source, tag, comm, request)
2784   #if defined(PETSC_HAVE_MPI_REDUCE_LOCAL)
2785     #define MPIU_Reduce_local(inbuf, inoutbuf, count, dtype, op) MPI_Reduce_local_c(inbuf, inoutbuf, (MPI_Count)(count), dtype, op)
2786   #endif
2787 
2788 /*MC
2789   MPIU_Scatterv - A replacement for `MPI_Scatterv()` that can be called with `PetscInt` types when PETSc is built for either 32-bit indices or 64-bit indices.
2790 
2791   Synopsis:
2792   #include <petscsys.h>
2793   PetscMPIInt MPIU_Scatterv(const void *sendbuf, const PetscInt sendcounts[], const PetscInt displs[], MPI_Datatype sendtype, void *recvbuf, PetscInt recvcount, MPI_Datatype recvtype, PetscMPIInt root, MPI_Comm comm)
2794 
2795   Collective
2796 
2797   Input Parameters:
2798 + sendbuf    - address of send buffer
2799 . sendcounts - non-negative `PetscInt` array (of length `comm` group size) specifying the number of elements to send to each MPI process
2800 . displs     - `PetscInt` array (of length `comm` group size). Entry i specifies the displacement (relative to `sendbuf`) from which to take the outgoing data to process i
2801 . sendtype   - data type of `sendbuf` elements
2802 . recvcount  - number of elements in `recvbuf` (non-negative integer)
2803 . recvtype   - data type of `recvbuf` elements
2804 . root       - Rank of the MPI root process, which will dispatch the data to scatter
2805 - comm       - `MPI_Comm` communicator
2806 
2807   Output Parameter:
2808 . recvbuf - the resulting scattered values on this MPI process
2809 
2810   Level: developer
2811 
2812   Notes:
2813   Should be wrapped with `PetscCallMPI()` for error checking
2814 
2815   This is different than most of the `MPIU_` wrappers in that all the count arguments are in `PetscInt`
2816 
2817 .seealso: [](stylePetscCount), `MPI_Allreduce()`, `MPIU_Gatherv()`
2818 M*/
2819 
2820   #if !defined(PETSC_USE_64BIT_INDICES)
2821     #define MPIU_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm) MPI_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm)
2822     #define MPIU_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm)  MPI_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm)
2823   #else
2824     #define MPIU_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm) MPI_Scatterv_c(sendbuf, (const MPI_Count *)(sendcounts), (const MPI_Aint *)(displs), sendtype, recvbuf, recvcount, recvtype, root, comm)
2825     #define MPIU_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm)  MPI_Gatherv_c(sendbuf, sendcount, sendtype, recvbuf, (const MPI_Count *)(recvcounts), (const MPI_Aint *)(displs), recvtype, root, comm)
2826   #endif
2827 
2828 #endif
2829 
2830 PETSC_EXTERN PetscMPIInt    MPIU_Allreduce_Private(const void *, void *, MPIU_Count, MPI_Datatype, MPI_Op, MPI_Comm);
2831 PETSC_EXTERN PetscErrorCode PetscCheckAllreduceSameLineAndCount_Private(MPI_Comm, const char *, PetscMPIInt, PetscMPIInt);
2832 
2833 #if defined(PETSC_USE_DEBUG)
2834 static inline unsigned int PetscStrHash(const char *str)
2835 {
2836   unsigned int c, hash = 5381;
2837 
2838   while ((c = (unsigned int)*str++)) hash = ((hash << 5) + hash) + c; /* hash * 33 + c */
2839   return hash;
2840 }
2841 #endif
2842 
2843 /*MC
2844   MPIU_Allreduce - A replacement for `MPI_Allreduce()` that (1) performs single-count `MPIU_INT` operations in `PetscInt64` to detect
2845                    integer overflows and (2) tries to determine if the call from all the MPI ranks occur in the
2846                    same place in the PETSc code. This helps to detect bugs where different MPI ranks follow different code paths
2847                    resulting in inconsistent and incorrect calls to `MPI_Allreduce()`.
2848 
2849   Synopsis:
2850   #include <petscsys.h>
2851   PetscMPIInt MPIU_Allreduce(void *indata,void *outdata,PetscCount count,MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
2852 
2853   Collective
2854 
2855   Input Parameters:
2856 + a     - pointer to the input data to be reduced
2857 . count - the number of MPI data items in `a` and `b`
2858 . dtype - the MPI datatype, for example `MPI_INT`
2859 . op    - the MPI operation, for example `MPI_SUM`
2860 - comm   - the MPI communicator on which the operation occurs
2861 
2862   Output Parameter:
2863 . b - the reduced values
2864 
2865   Level: developer
2866 
2867   Note:
2868   Should be wrapped with `PetscCallMPI()` for error checking
2869 
2870 .seealso: [](stylePetscCount), `MPI_Allreduce()`
2871 M*/
2872 #if defined(PETSC_USE_DEBUG)
2873   #define MPIU_Allreduce(a, b, count, dtype, op, comm) \
2874     PetscMacroReturnStandard( \
2875     PetscCall(PetscCheckAllreduceSameLineAndCount_Private((comm), __FILE__, (PetscMPIInt)__LINE__, (PetscMPIInt)(count))); \
2876     PetscCallMPI(MPIU_Allreduce_Private((a), (b), (count), (dtype), (op), (comm)));)
2877 #else
2878   #define MPIU_Allreduce(a, b, count, dtype, op, comm) MPIU_Allreduce_Private((a), (b), (count), (dtype), (op), (comm))
2879 #endif
2880 
2881 /* this is a vile hack */
2882 #if defined(PETSC_HAVE_NECMPI)
2883   #if !defined(PETSC_NECMPI_VERSION_MAJOR) || !defined(PETSC_NECMPI_VERSION_MINOR) || PETSC_NECMPI_VERSION_MAJOR < 2 || (PETSC_NECMPI_VERSION_MAJOR == 2 && PETSC_NECMPI_VERSION_MINOR < 18)
2884     #define MPI_Type_free(a) (*(a) = MPI_DATATYPE_NULL, 0);
2885   #endif
2886 #endif
2887 
2888 /*
2889     List of external packages and queries on it
2890 */
2891 PETSC_EXTERN PetscErrorCode PetscHasExternalPackage(const char[], PetscBool *);
2892 
2893 /* this cannot go here because it may be in a different shared library */
2894 PETSC_EXTERN PetscErrorCode PCMPIServerBegin(void);
2895 PETSC_EXTERN PetscErrorCode PCMPIServerEnd(void);
2896 PETSC_EXTERN PetscBool      PCMPIServerActive;
2897 PETSC_EXTERN PetscBool      PCMPIServerInSolve;
2898 PETSC_EXTERN PetscBool      PCMPIServerUseShmget;
2899 PETSC_EXTERN PetscErrorCode PetscShmgetAllocateArray(size_t, size_t, void **);
2900 PETSC_EXTERN PetscErrorCode PetscShmgetDeallocateArray(void **);
2901 PETSC_EXTERN PetscErrorCode PetscShmgetMapAddresses(MPI_Comm, PetscInt, const void **, void **);
2902 PETSC_EXTERN PetscErrorCode PetscShmgetUnmapAddresses(PetscInt, void **);
2903 PETSC_EXTERN PetscErrorCode PetscShmgetAddressesFinalize(void);
2904 
2905 typedef struct {
2906   PetscInt n;
2907   void    *addr[3];
2908 } PCMPIServerAddresses;
2909 PETSC_EXTERN PetscCtxDestroyFn PCMPIServerAddressesDestroy;
2910 
2911 #define PETSC_HAVE_FORTRAN PETSC_DEPRECATED_MACRO(3, 20, 0, "PETSC_USE_FORTRAN_BINDINGS", ) PETSC_USE_FORTRAN_BINDINGS
2912 
2913 PETSC_EXTERN PetscErrorCode PetscBLASSetNumThreads(PetscInt);
2914 PETSC_EXTERN PetscErrorCode PetscBLASGetNumThreads(PetscInt *);
2915 
2916 /*MC
2917    PetscSafePointerPlusOffset - Checks that a pointer is not `NULL` before applying an offset
2918 
2919    Level: beginner
2920 
2921    Note:
2922    This is needed to avoid errors with undefined-behavior sanitizers such as
2923    UBSan, assuming PETSc has been configured with `-fsanitize=undefined` as part of the compiler flags
2924 M*/
2925 #define PetscSafePointerPlusOffset(ptr, offset) ((ptr) ? (ptr) + (offset) : NULL)
2926 
2927 /* this is required to force PetscDevice to be visible at the system level for the Fortran interface */
2928 #include <petscdevicetypes.h>
2929 
2930 #if PetscDefined(USE_DEBUG) && !PetscDefined(HAVE_THREADSAFETY)
2931 PETSC_EXTERN PetscErrorCode PetscStackView(FILE *);
2932 #else
2933   #define PetscStackView(file) PETSC_SUCCESS
2934 #endif
2935