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