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