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