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