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