xref: /petsc/include/petscsys.h (revision fef353a4b7d6ea90a49c0912b64f64e31ef80aa8)
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 NULL
214 
215 /* This is deprecated */
216 #define PETSC_NULL NULL
217 
218 /*MC
219     PETSC_DECIDE - standard way of passing in integer or floating point parameter
220        where you wish PETSc to use the default.
221 
222    Level: beginner
223 
224 .seealso: `PETSC_DEFAULT`, `PETSC_IGNORE`, `PETSC_DETERMINE`
225 
226 M*/
227 #define PETSC_DECIDE -1
228 
229 /*MC
230     PETSC_DETERMINE - standard way of passing in integer or floating point parameter
231        where you wish PETSc to compute the required value.
232 
233    Level: beginner
234 
235    Developer Note:
236    I would like to use const `PetscInt` `PETSC_DETERMINE` = `PETSC_DECIDE`; but for
237      some reason this is not allowed by the standard even though `PETSC_DECIDE` is a constant value.
238 
239 .seealso: `PETSC_DECIDE`, `PETSC_DEFAULT`, `PETSC_IGNORE`, `VecSetSizes()`
240 
241 M*/
242 #define PETSC_DETERMINE PETSC_DECIDE
243 
244 /*MC
245     PETSC_DEFAULT - standard way of passing in integer or floating point parameter
246        where you wish PETSc to use the default.
247 
248    Level: beginner
249 
250    Fortran Note:
251    You need to use `PETSC_DEFAULT_INTEGER` or `PETSC_DEFAULT_REAL`.
252 
253 .seealso: `PETSC_DECIDE`, `PETSC_IGNORE`, `PETSC_DETERMINE`
254 
255 M*/
256 #define PETSC_DEFAULT -2
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)((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) 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))
563 
564 /*MC
565    PetscCalloc3 - Allocates 3 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
566 
567    Synopsis:
568     #include <petscsys.h>
569    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
570 
571    Not Collective
572 
573    Input Parameters:
574 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
575 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
576 -  m3 - number of elements to allocate in 3rd chunk  (may be zero)
577 
578    Output Parameters:
579 +  r1 - memory allocated in first chunk
580 .  r2 - memory allocated in second chunk
581 -  r3 - memory allocated in third chunk
582 
583    Level: developer
584 
585 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc2()`, `PetscMalloc3()`, `PetscFree3()`
586 
587 M*/
588 #define PetscCalloc3(m1, r1, m2, r2, m3, r3) PetscMallocA(3, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (size_t)((size_t)m1) * sizeof(**(r1)), (r1), (size_t)((size_t)m2) * sizeof(**(r2)), (r2), (size_t)((size_t)m3) * sizeof(**(r3)), (r3))
589 
590 /*MC
591    PetscMalloc4 - Allocates 4 arrays of memory, all aligned to `PETSC_MEMALIGN`
592 
593    Synopsis:
594     #include <petscsys.h>
595    PetscErrorCode PetscMalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
596 
597    Not Collective
598 
599    Input Parameters:
600 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
601 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
602 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
603 -  m4 - number of elements to allocate in 4th chunk  (may be zero)
604 
605    Output Parameters:
606 +  r1 - memory allocated in first chunk
607 .  r2 - memory allocated in second chunk
608 .  r3 - memory allocated in third chunk
609 -  r4 - memory allocated in fourth chunk
610 
611    Level: developer
612 
613 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc4()`, `PetscFree4()`
614 
615 M*/
616 #define PetscMalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \
617   PetscMallocA(4, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (size_t)((size_t)m1) * sizeof(**(r1)), (r1), (size_t)((size_t)m2) * sizeof(**(r2)), (r2), (size_t)((size_t)m3) * sizeof(**(r3)), (r3), (size_t)((size_t)m4) * sizeof(**(r4)), (r4))
618 
619 /*MC
620    PetscCalloc4 - Allocates 4 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
621 
622    Synopsis:
623     #include <petscsys.h>
624    PetscErrorCode PetscCalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
625 
626    Not Collective
627 
628    Input Parameters:
629 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
630 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
631 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
632 -  m4 - number of elements to allocate in 4th chunk  (may be zero)
633 
634    Output Parameters:
635 +  r1 - memory allocated in first chunk
636 .  r2 - memory allocated in second chunk
637 .  r3 - memory allocated in third chunk
638 -  r4 - memory allocated in fourth chunk
639 
640    Level: developer
641 
642 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc4()`, `PetscFree4()`
643 
644 M*/
645 #define PetscCalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \
646   PetscMallocA(4, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (size_t)((size_t)m1) * sizeof(**(r1)), (r1), (size_t)((size_t)m2) * sizeof(**(r2)), (r2), (size_t)((size_t)m3) * sizeof(**(r3)), (r3), (size_t)((size_t)m4) * sizeof(**(r4)), (r4))
647 
648 /*MC
649    PetscMalloc5 - Allocates 5 arrays of memory, all aligned to `PETSC_MEMALIGN`
650 
651    Synopsis:
652     #include <petscsys.h>
653    PetscErrorCode PetscMalloc5(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5)
654 
655    Not Collective
656 
657    Input Parameters:
658 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
659 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
660 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
661 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
662 -  m5 - number of elements to allocate in 5th chunk  (may be zero)
663 
664    Output Parameters:
665 +  r1 - memory allocated in first chunk
666 .  r2 - memory allocated in second chunk
667 .  r3 - memory allocated in third chunk
668 .  r4 - memory allocated in fourth chunk
669 -  r5 - memory allocated in fifth chunk
670 
671    Level: developer
672 
673 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc5()`, `PetscFree5()`
674 
675 M*/
676 #define PetscMalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \
677   PetscMallocA(5, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (size_t)((size_t)m1) * sizeof(**(r1)), (r1), (size_t)((size_t)m2) * sizeof(**(r2)), (r2), (size_t)((size_t)m3) * sizeof(**(r3)), (r3), (size_t)((size_t)m4) * sizeof(**(r4)), (r4), (size_t)((size_t)m5) * sizeof(**(r5)), (r5))
678 
679 /*MC
680    PetscCalloc5 - Allocates 5 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
681 
682    Synopsis:
683     #include <petscsys.h>
684    PetscErrorCode PetscCalloc5(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5)
685 
686    Not Collective
687 
688    Input Parameters:
689 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
690 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
691 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
692 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
693 -  m5 - number of elements to allocate in 5th chunk  (may be zero)
694 
695    Output Parameters:
696 +  r1 - memory allocated in first chunk
697 .  r2 - memory allocated in second chunk
698 .  r3 - memory allocated in third chunk
699 .  r4 - memory allocated in fourth chunk
700 -  r5 - memory allocated in fifth chunk
701 
702    Level: developer
703 
704 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc5()`, `PetscFree5()`
705 
706 M*/
707 #define PetscCalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \
708   PetscMallocA(5, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (size_t)((size_t)m1) * sizeof(**(r1)), (r1), (size_t)((size_t)m2) * sizeof(**(r2)), (r2), (size_t)((size_t)m3) * sizeof(**(r3)), (r3), (size_t)((size_t)m4) * sizeof(**(r4)), (r4), (size_t)((size_t)m5) * sizeof(**(r5)), (r5))
709 
710 /*MC
711    PetscMalloc6 - Allocates 6 arrays of memory, all aligned to `PETSC_MEMALIGN`
712 
713    Synopsis:
714     #include <petscsys.h>
715    PetscErrorCode PetscMalloc6(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5,size_t m6,type **r6)
716 
717    Not Collective
718 
719    Input Parameters:
720 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
721 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
722 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
723 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
724 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
725 -  m6 - number of elements to allocate in 6th chunk  (may be zero)
726 
727    Output Parameteasr:
728 +  r1 - memory allocated in first chunk
729 .  r2 - memory allocated in second chunk
730 .  r3 - memory allocated in third chunk
731 .  r4 - memory allocated in fourth chunk
732 .  r5 - memory allocated in fifth chunk
733 -  r6 - memory allocated in sixth chunk
734 
735    Level: developer
736 
737 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc6()`, `PetscFree3()`, `PetscFree4()`, `PetscFree5()`, `PetscFree6()`
738 
739 M*/
740 #define PetscMalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \
741   PetscMallocA(6, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (size_t)((size_t)m1) * sizeof(**(r1)), (r1), (size_t)((size_t)m2) * sizeof(**(r2)), (r2), (size_t)((size_t)m3) * sizeof(**(r3)), (r3), (size_t)((size_t)m4) * sizeof(**(r4)), (r4), (size_t)((size_t)m5) * sizeof(**(r5)), (r5), (size_t)((size_t)m6) * sizeof(**(r6)), (r6))
742 
743 /*MC
744    PetscCalloc6 - Allocates 6 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
745 
746    Synopsis:
747     #include <petscsys.h>
748    PetscErrorCode PetscCalloc6(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5,size_t m6,type **r6)
749 
750    Not Collective
751 
752    Input Parameters:
753 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
754 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
755 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
756 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
757 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
758 -  m6 - number of elements to allocate in 6th chunk  (may be zero)
759 
760    Output Parameters:
761 +  r1 - memory allocated in first chunk
762 .  r2 - memory allocated in second chunk
763 .  r3 - memory allocated in third chunk
764 .  r4 - memory allocated in fourth chunk
765 .  r5 - memory allocated in fifth chunk
766 -  r6 - memory allocated in sixth chunk
767 
768    Level: developer
769 
770 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscMalloc6()`, `PetscFree6()`
771 
772 M*/
773 #define PetscCalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \
774   PetscMallocA(6, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (size_t)((size_t)m1) * sizeof(**(r1)), (r1), (size_t)((size_t)m2) * sizeof(**(r2)), (r2), (size_t)((size_t)m3) * sizeof(**(r3)), (r3), (size_t)((size_t)m4) * sizeof(**(r4)), (r4), (size_t)((size_t)m5) * sizeof(**(r5)), (r5), (size_t)((size_t)m6) * sizeof(**(r6)), (r6))
775 
776 /*MC
777    PetscMalloc7 - Allocates 7 arrays of memory, all aligned to `PETSC_MEMALIGN`
778 
779    Synopsis:
780     #include <petscsys.h>
781    PetscErrorCode PetscMalloc7(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5,size_t m6,type **r6,size_t m7,type **r7)
782 
783    Not Collective
784 
785    Input Parameters:
786 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
787 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
788 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
789 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
790 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
791 .  m6 - number of elements to allocate in 6th chunk  (may be zero)
792 -  m7 - number of elements to allocate in 7th chunk  (may be zero)
793 
794    Output Parameters:
795 +  r1 - memory allocated in first chunk
796 .  r2 - memory allocated in second chunk
797 .  r3 - memory allocated in third chunk
798 .  r4 - memory allocated in fourth chunk
799 .  r5 - memory allocated in fifth chunk
800 .  r6 - memory allocated in sixth chunk
801 -  r7 - memory allocated in seventh chunk
802 
803    Level: developer
804 
805 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc7()`, `PetscFree7()`
806 
807 M*/
808 #define PetscMalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \
809   PetscMallocA(7, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (size_t)((size_t)m1) * sizeof(**(r1)), (r1), (size_t)((size_t)m2) * sizeof(**(r2)), (r2), (size_t)((size_t)m3) * sizeof(**(r3)), (r3), (size_t)((size_t)m4) * sizeof(**(r4)), (r4), (size_t)((size_t)m5) * sizeof(**(r5)), (r5), (size_t)((size_t)m6) * sizeof(**(r6)), (r6), (size_t)((size_t)m7) * sizeof(**(r7)), (r7))
810 
811 /*MC
812    PetscCalloc7 - Allocates 7 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
813 
814    Synopsis:
815     #include <petscsys.h>
816    PetscErrorCode PetscCalloc7(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5,size_t m6,type **r6,size_t m7,type **r7)
817 
818    Not Collective
819 
820    Input Parameters:
821 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
822 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
823 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
824 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
825 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
826 .  m6 - number of elements to allocate in 6th chunk  (may be zero)
827 -  m7 - number of elements to allocate in 7th chunk  (may be zero)
828 
829    Output Parameters:
830 +  r1 - memory allocated in first chunk
831 .  r2 - memory allocated in second chunk
832 .  r3 - memory allocated in third chunk
833 .  r4 - memory allocated in fourth chunk
834 .  r5 - memory allocated in fifth chunk
835 .  r6 - memory allocated in sixth chunk
836 -  r7 - memory allocated in seventh chunk
837 
838    Level: developer
839 
840 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscMalloc7()`, `PetscFree7()`
841 
842 M*/
843 #define PetscCalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \
844   PetscMallocA(7, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (size_t)((size_t)m1) * sizeof(**(r1)), (r1), (size_t)((size_t)m2) * sizeof(**(r2)), (r2), (size_t)((size_t)m3) * sizeof(**(r3)), (r3), (size_t)((size_t)m4) * sizeof(**(r4)), (r4), (size_t)((size_t)m5) * sizeof(**(r5)), (r5), (size_t)((size_t)m6) * sizeof(**(r6)), (r6), (size_t)((size_t)m7) * sizeof(**(r7)), (r7))
845 
846 /*MC
847    PetscNew - Allocates memory of a particular type, zeros the memory! Aligned to `PETSC_MEMALIGN`
848 
849    Synopsis:
850     #include <petscsys.h>
851    PetscErrorCode PetscNew(type **result)
852 
853    Not Collective
854 
855    Output Parameter:
856 .  result - memory allocated, sized to match pointer type
857 
858    Level: beginner
859 
860 .seealso: `PetscFree()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc1()`
861 
862 M*/
863 #define PetscNew(b) PetscCalloc1(1, (b))
864 
865 #define PetscNewLog(o, b) PETSC_DEPRECATED_MACRO("GCC warning \"PetscNewLog is deprecated, use PetscNew() instead (since version 3.18)\"") PetscNew((b))
866 
867 /*MC
868    PetscFree - Frees memory
869 
870    Synopsis:
871     #include <petscsys.h>
872    PetscErrorCode PetscFree(void *memory)
873 
874    Not Collective
875 
876    Input Parameter:
877 .   memory - memory to free (the pointer is ALWAYS set to NULL upon success)
878 
879    Level: beginner
880 
881    Notes:
882    Do not free memory obtained with `PetscMalloc2()`, `PetscCalloc2()` etc, they must be freed with `PetscFree2()` etc.
883 
884    It is safe to call `PetscFree()` on a NULL pointer.
885 
886 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc1()`
887 
888 M*/
889 #define PetscFree(a) ((*PetscTrFree)((void *)(a), __LINE__, PETSC_FUNCTION_NAME, __FILE__) || ((a) = NULL, 0))
890 
891 /*MC
892    PetscFree2 - Frees 2 chunks of memory obtained with `PetscMalloc2()`
893 
894    Synopsis:
895     #include <petscsys.h>
896    PetscErrorCode PetscFree2(void *memory1,void *memory2)
897 
898    Not Collective
899 
900    Input Parameters:
901 +   memory1 - memory to free
902 -   memory2 - 2nd memory to free
903 
904    Level: developer
905 
906    Notes:
907     Memory must have been obtained with `PetscMalloc2()`
908 
909 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`
910 
911 M*/
912 #define PetscFree2(m1, m2) PetscFreeA(2, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2))
913 
914 /*MC
915    PetscFree3 - Frees 3 chunks of memory obtained with `PetscMalloc3()`
916 
917    Synopsis:
918     #include <petscsys.h>
919    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)
920 
921    Not Collective
922 
923    Input Parameters:
924 +   memory1 - memory to free
925 .   memory2 - 2nd memory to free
926 -   memory3 - 3rd memory to free
927 
928    Level: developer
929 
930    Notes:
931     Memory must have been obtained with `PetscMalloc3()`
932 
933 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`
934 
935 M*/
936 #define PetscFree3(m1, m2, m3) PetscFreeA(3, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3))
937 
938 /*MC
939    PetscFree4 - Frees 4 chunks of memory obtained with `PetscMalloc4()`
940 
941    Synopsis:
942     #include <petscsys.h>
943    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)
944 
945    Not Collective
946 
947    Input Parameters:
948 +   m1 - memory to free
949 .   m2 - 2nd memory to free
950 .   m3 - 3rd memory to free
951 -   m4 - 4th memory to free
952 
953    Level: developer
954 
955    Notes:
956     Memory must have been obtained with `PetscMalloc4()`
957 
958 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`
959 
960 M*/
961 #define PetscFree4(m1, m2, m3, m4) PetscFreeA(4, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4))
962 
963 /*MC
964    PetscFree5 - Frees 5 chunks of memory obtained with `PetscMalloc5()`
965 
966    Synopsis:
967     #include <petscsys.h>
968    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)
969 
970    Not Collective
971 
972    Input Parameters:
973 +   m1 - memory to free
974 .   m2 - 2nd memory to free
975 .   m3 - 3rd memory to free
976 .   m4 - 4th memory to free
977 -   m5 - 5th memory to free
978 
979    Level: developer
980 
981    Notes:
982     Memory must have been obtained with `PetscMalloc5()`
983 
984 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`
985 
986 M*/
987 #define PetscFree5(m1, m2, m3, m4, m5) PetscFreeA(5, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5))
988 
989 /*MC
990    PetscFree6 - Frees 6 chunks of memory obtained with `PetscMalloc6()`
991 
992    Synopsis:
993     #include <petscsys.h>
994    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)
995 
996    Not Collective
997 
998    Input Parameters:
999 +   m1 - memory to free
1000 .   m2 - 2nd memory to free
1001 .   m3 - 3rd memory to free
1002 .   m4 - 4th memory to free
1003 .   m5 - 5th memory to free
1004 -   m6 - 6th memory to free
1005 
1006    Level: developer
1007 
1008    Notes:
1009     Memory must have been obtained with `PetscMalloc6()`
1010 
1011 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()`
1012 
1013 M*/
1014 #define PetscFree6(m1, m2, m3, m4, m5, m6) PetscFreeA(6, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6))
1015 
1016 /*MC
1017    PetscFree7 - Frees 7 chunks of memory obtained with `PetscMalloc7()`
1018 
1019    Synopsis:
1020     #include <petscsys.h>
1021    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)
1022 
1023    Not Collective
1024 
1025    Input Parameters:
1026 +   m1 - memory to free
1027 .   m2 - 2nd memory to free
1028 .   m3 - 3rd memory to free
1029 .   m4 - 4th memory to free
1030 .   m5 - 5th memory to free
1031 .   m6 - 6th memory to free
1032 -   m7 - 7th memory to free
1033 
1034    Level: developer
1035 
1036    Notes:
1037     Memory must have been obtained with `PetscMalloc7()`
1038 
1039 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()`,
1040           `PetscMalloc7()`
1041 
1042 M*/
1043 #define PetscFree7(m1, m2, m3, m4, m5, m6, m7) PetscFreeA(7, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6), &(m7))
1044 
1045 PETSC_EXTERN PetscErrorCode PetscMallocA(int, PetscBool, int, const char *, const char *, size_t, void *, ...);
1046 PETSC_EXTERN PetscErrorCode PetscFreeA(int, int, const char *, const char *, void *, ...);
1047 PETSC_EXTERN                PetscErrorCode (*PetscTrMalloc)(size_t, PetscBool, int, const char[], const char[], void **);
1048 PETSC_EXTERN                PetscErrorCode (*PetscTrFree)(void *, int, const char[], const char[]);
1049 PETSC_EXTERN                PetscErrorCode (*PetscTrRealloc)(size_t, int, const char[], const char[], void **);
1050 PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1051 PETSC_EXTERN PetscErrorCode PetscMallocSet(PetscErrorCode (*)(size_t, PetscBool, int, const char[], const char[], void **), PetscErrorCode (*)(void *, int, const char[], const char[]), PetscErrorCode (*)(size_t, int, const char[], const char[], void **));
1052 PETSC_EXTERN PetscErrorCode PetscMallocClear(void);
1053 
1054 /*
1055   Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1056 */
1057 PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1058 PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);
1059 #if defined(PETSC_HAVE_CUDA)
1060 PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void);
1061 PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void);
1062 #endif
1063 #if defined(PETSC_HAVE_HIP)
1064 PETSC_EXTERN PetscErrorCode PetscMallocSetHIPHost(void);
1065 PETSC_EXTERN PetscErrorCode PetscMallocResetHIPHost(void);
1066 #endif
1067 
1068 #define MPIU_PETSCLOGDOUBLE  MPI_DOUBLE
1069 #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION
1070 
1071 /*
1072    Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1073 */
1074 PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1075 PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1076 PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1077 PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1078 PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1079 PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int, PetscLogDouble *);
1080 PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool, PetscBool);
1081 PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool *, PetscBool *, PetscBool *);
1082 PETSC_EXTERN PetscErrorCode PetscMallocValidate(int, const char[], const char[]);
1083 PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1084 PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool *);
1085 PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool);
1086 PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool *);
1087 
1088 PETSC_EXTERN const char *const PetscDataTypes[];
1089 PETSC_EXTERN PetscErrorCode    PetscDataTypeToMPIDataType(PetscDataType, MPI_Datatype *);
1090 PETSC_EXTERN PetscErrorCode    PetscMPIDataTypeToPetscDataType(MPI_Datatype, PetscDataType *);
1091 PETSC_EXTERN PetscErrorCode    PetscDataTypeGetSize(PetscDataType, size_t *);
1092 PETSC_EXTERN PetscErrorCode    PetscDataTypeFromString(const char *, PetscDataType *, PetscBool *);
1093 
1094 /*
1095     Basic memory and string operations. These are usually simple wrappers
1096    around the basic Unix system calls, but a few of them have additional
1097    functionality and/or error checking.
1098 */
1099 PETSC_EXTERN PetscErrorCode PetscMemcmp(const void *, const void *, size_t, PetscBool *);
1100 PETSC_EXTERN PetscErrorCode PetscStrlen(const char[], size_t *);
1101 PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[], char, int *, char ***);
1102 PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int, char **);
1103 PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[], const char[], PetscBool *);
1104 PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[], const char[], PetscBool *);
1105 PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[], const char[], PetscBool *);
1106 PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[], const char[], size_t, PetscBool *);
1107 PETSC_EXTERN PetscErrorCode PetscStrcpy(char[], const char[]);
1108 PETSC_EXTERN PetscErrorCode PetscStrcat(char[], const char[]);
1109 PETSC_EXTERN PetscErrorCode PetscStrlcat(char[], const char[], size_t);
1110 PETSC_EXTERN PetscErrorCode PetscStrncpy(char[], const char[], size_t);
1111 PETSC_EXTERN PetscErrorCode PetscStrchr(const char[], char, char *[]);
1112 PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1113 PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1114 PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[], char, char *[]);
1115 PETSC_EXTERN PetscErrorCode PetscStrstr(const char[], const char[], char *[]);
1116 PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[], const char[], char *[]);
1117 PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[], const char[], PetscBool *);
1118 PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[], const char[], PetscBool *);
1119 PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[], const char *const *, PetscInt *);
1120 PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[], char *[]);
1121 PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const *, char ***);
1122 PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char ***);
1123 PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt, const char *const *, char ***);
1124 PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt, char ***);
1125 PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm, const char[], char[], size_t);
1126 
1127 PETSC_EXTERN void PetscStrcmpNoError(const char[], const char[], PetscBool *);
1128 
1129 PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[], const char, PetscToken *);
1130 PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken, char *[]);
1131 PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken *);
1132 
1133 PETSC_EXTERN PetscErrorCode PetscStrInList(const char[], const char[], char, PetscBool *);
1134 PETSC_EXTERN const char    *PetscBasename(const char[]);
1135 PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt, const char *const *, const char *, PetscInt *, PetscBool *);
1136 PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const *, const char *, PetscEnum *, PetscBool *);
1137 
1138 /*
1139    These are MPI operations for MPI_Allreduce() etc
1140 */
1141 PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1142 #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1143 PETSC_EXTERN MPI_Op MPIU_SUM;
1144 PETSC_EXTERN MPI_Op MPIU_MAX;
1145 PETSC_EXTERN MPI_Op MPIU_MIN;
1146 #else
1147   #define MPIU_SUM MPI_SUM
1148   #define MPIU_MAX MPI_MAX
1149   #define MPIU_MIN MPI_MIN
1150 #endif
1151 PETSC_EXTERN MPI_Op         Petsc_Garbage_SetIntersectOp;
1152 PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm, const PetscInt[], PetscInt *, PetscInt *);
1153 
1154 #if defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
1155 /*MC
1156     MPIU_SUM___FP16___FLOAT128 - MPI_Op that acts as a replacement for MPI_SUM with
1157     custom MPI_Datatype `MPIU___FLOAT128`, `MPIU___COMPLEX128`, and `MPIU___FP16`.
1158 
1159    Level: advanced
1160 
1161    Developer Note:
1162    This should be unified with `MPIU_SUM`
1163 
1164 .seealso: `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`
1165 M*/
1166 PETSC_EXTERN MPI_Op MPIU_SUM___FP16___FLOAT128;
1167 #endif
1168 PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm, const PetscInt[], PetscInt *, PetscInt *);
1169 
1170 PETSC_EXTERN PetscErrorCode MPIULong_Send(void *, PetscInt, MPI_Datatype, PetscMPIInt, PetscMPIInt, MPI_Comm) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(1, 3);
1171 PETSC_EXTERN PetscErrorCode MPIULong_Recv(void *, PetscInt, MPI_Datatype, PetscMPIInt, PetscMPIInt, MPI_Comm) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(1, 3);
1172 
1173 PETSC_EXTERN const char *const PetscFileModes[];
1174 
1175 /*
1176     Defines PETSc error handling.
1177 */
1178 #include <petscerror.h>
1179 
1180 PETSC_EXTERN PetscBool   PetscCIEnabled;                    /* code is running in the PETSc test harness CI */
1181 PETSC_EXTERN PetscBool   PetscCIEnabledPortableErrorOutput; /* error output is stripped to ensure portability of error messages across systems */
1182 PETSC_EXTERN const char *PetscCIFilename(const char *);
1183 PETSC_EXTERN int         PetscCILinenumber(int);
1184 
1185 #define PETSC_SMALLEST_CLASSID 1211211
1186 PETSC_EXTERN PetscClassId   PETSC_LARGEST_CLASSID;
1187 PETSC_EXTERN PetscClassId   PETSC_OBJECT_CLASSID;
1188 PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[], PetscClassId *);
1189 PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject, PetscObjectId *);
1190 PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject, PetscObjectId, PetscBool *);
1191 
1192 /*
1193    Routines that get memory usage information from the OS
1194 */
1195 PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1196 PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1197 PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1198 PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);
1199 
1200 PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);
1201 
1202 /*
1203    Initialization of PETSc
1204 */
1205 PETSC_EXTERN PetscErrorCode PetscInitialize(int *, char ***, const char[], const char[]);
1206 PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int, char **, const char[], const char[]);
1207 PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1208 PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1209 PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1210 PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1211 PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1212 PETSC_EXTERN PetscErrorCode PetscGetArgs(int *, char ***);
1213 PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1214 PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);
1215 
1216 PETSC_EXTERN PetscErrorCode PetscEnd(void);
1217 PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);
1218 
1219 PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[], const char[]);
1220 PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1221 PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1222 PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject, const char[]);
1223 
1224 PETSC_EXTERN PetscErrorCode PetscMonitorCompare(PetscErrorCode (*)(void), void *, PetscErrorCode (*)(void **), PetscErrorCode (*)(void), void *, PetscErrorCode (*)(void **), PetscBool *);
1225 
1226 /*
1227      These are so that in extern C code we can caste function pointers to non-extern C
1228    function pointers. Since the regular C++ code expects its function pointers to be C++
1229 */
1230 PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1231 PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1232 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);
1233 
1234 /*
1235     Functions that can act on any PETSc object.
1236 */
1237 PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject *);
1238 PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject, MPI_Comm *);
1239 PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject, PetscClassId *);
1240 PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject, const char *[]);
1241 PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject, const char *[]);
1242 PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject, const char[]);
1243 PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject, const char *[]);
1244 PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject, PetscInt);
1245 PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject, PetscInt *);
1246 PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject, PetscObject, PetscInt);
1247 PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1248 PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject, PetscInt *);
1249 PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1250 PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject, PetscMPIInt *);
1251 PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject, const char[], PetscObject);
1252 PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject, const char[]);
1253 PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject, const char[], PetscObject *);
1254 PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject, const char[], void (*)(void));
1255 #define PetscObjectComposeFunction(a, b, d) PetscObjectComposeFunction_Private(a, b, (PetscVoidFunction)(d))
1256 PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1257 PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1258 PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1259 PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject, PetscObject);
1260 PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm, PetscMPIInt *);
1261 
1262 #include <petscviewertypes.h>
1263 #include <petscoptions.h>
1264 
1265 PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer, PetscBool, PetscLogDouble);
1266 PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool *);
1267 
1268 PETSC_EXTERN PetscErrorCode PetscObjectsListGetGlobalNumbering(MPI_Comm, PetscInt, PetscObject *, PetscInt *, PetscInt *);
1269 
1270 PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer, const char[]);
1271 PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject, PetscViewer);
1272 PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject, PetscViewer);
1273 #define PetscObjectQueryFunction(obj, name, fptr) PetscObjectQueryFunction_Private((obj), (name), (PetscVoidFunction *)(fptr))
1274 PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject, const char[], void (**)(void));
1275 PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject, const char[]);
1276 PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject, const char[]);
1277 PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject, const char[]);
1278 PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject, const char *[]);
1279 PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject, const char[]);
1280 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1281 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1282 PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject, PetscObject, const char[]);
1283 PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1284 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject, const char[], PetscBool *);
1285 PETSC_EXTERN PetscErrorCode PetscObjectObjectTypeCompare(PetscObject, PetscObject, PetscBool *);
1286 PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject, const char[], PetscBool *);
1287 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject, PetscBool *, const char[], ...);
1288 PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject, PetscBool *, const char[], ...);
1289 PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1290 PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);
1291 
1292 #if defined(PETSC_HAVE_SAWS)
1293 PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1294 PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1295 PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject, PetscBool);
1296 PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1297 PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1298 PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1299 PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1300 PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1301 PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1302 PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);
1303 
1304 #else
1305   #define PetscSAWsBlock()                  0
1306   #define PetscObjectSAWsViewOff(obj)       0
1307   #define PetscObjectSAWsSetBlock(obj, flg) 0
1308   #define PetscObjectSAWsBlock(obj)         0
1309   #define PetscObjectSAWsGrantAccess(obj)   0
1310   #define PetscObjectSAWsTakeAccess(obj)    0
1311   #define PetscStackViewSAWs()              0
1312   #define PetscStackSAWsViewOff()           0
1313   #define PetscStackSAWsTakeAccess()
1314   #define PetscStackSAWsGrantAccess()
1315 
1316 #endif
1317 
1318 PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[], PetscDLMode, PetscDLHandle *);
1319 PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1320 PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle, const char[], void **);
1321 PETSC_EXTERN PetscErrorCode PetscDLAddr(void (*)(void), char **);
1322 #ifdef PETSC_HAVE_CXX
1323 PETSC_EXTERN PetscErrorCode PetscDemangleSymbol(const char *, char **);
1324 #endif
1325 
1326 PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void *, PetscStack **);
1327 
1328 PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE *, PetscBool);
1329 PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList *);
1330 PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList, const char[], PetscObject *);
1331 PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList, PetscObject, char **, PetscBool *);
1332 PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *, const char[], PetscObject);
1333 PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *, const char[]);
1334 PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList, PetscObjectList *);
1335 
1336 /*
1337     Dynamic library lists. Lists of names of routines in objects or in dynamic
1338   link libraries that will be loaded as needed.
1339 */
1340 
1341 #define PetscFunctionListAdd(list, name, fptr) PetscFunctionListAdd_Private((list), (name), (PetscVoidFunction)(fptr))
1342 PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList *, const char[], void (*)(void));
1343 PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList *);
1344 PETSC_EXTERN PetscErrorCode PetscFunctionListClear(PetscFunctionList);
1345 #define PetscFunctionListFind(list, name, fptr) PetscFunctionListFind_Private((list), (name), (PetscVoidFunction *)(fptr))
1346 PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList, const char[], void (**)(void));
1347 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm, FILE *, const char[], const char[], const char[], const char[], PetscFunctionList, const char[], const char[]);
1348 PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList, PetscFunctionList *);
1349 PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList, PetscViewer);
1350 PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList, const char ***, int *);
1351 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintNonEmpty(PetscFunctionList);
1352 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintAll(void);
1353 
1354 PETSC_EXTERN PetscDLLibrary PetscDLLibrariesLoaded;
1355 PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm, PetscDLLibrary *, const char[]);
1356 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm, PetscDLLibrary *, const char[]);
1357 PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm, PetscDLLibrary *, const char[], const char[], void **);
1358 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1359 PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm, const char[], char *, size_t, PetscBool *);
1360 PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm, const char[], PetscDLLibrary *);
1361 PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);
1362 
1363 /*
1364      Useful utility routines
1365 */
1366 PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm, PetscInt *, PetscInt *);
1367 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm, PetscInt, PetscInt *, PetscInt *);
1368 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm, PetscInt *, PetscInt *);
1369 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm, PetscMPIInt);
1370 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm, PetscMPIInt);
1371 PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1372 PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE *);
1373 PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm, const PetscInt[2], PetscInt[2]);
1374 PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm, const PetscReal[2], PetscReal[2]);
1375 
1376 /*MC
1377     PetscNot - negates a logical type value and returns result as a `PetscBool`
1378 
1379     Level: beginner
1380 
1381     Note:
1382     This is useful in cases like
1383 .vb
1384      int        *a;
1385      PetscBool  flag = PetscNot(a)
1386 .ve
1387      where !a would not return a `PetscBool` because we cannot provide a cast from int to `PetscBool` in C.
1388 
1389 .seealso: `PetscBool`, `PETSC_TRUE`, `PETSC_FALSE`
1390 M*/
1391 #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)
1392 
1393 /*MC
1394     PetscHelpPrintf - Prints help messages.
1395 
1396     Not Collective, only applies on rank 0; No Fortran Support
1397 
1398    Synopsis:
1399     #include <petscsys.h>
1400      PetscErrorCode (*PetscHelpPrintf)(MPI_Comm comm, const char format[],args);
1401 
1402     Input Parameters:
1403 +  comm - the MPI communicator over which the help message is printed
1404 .  format - the usual printf() format string
1405 -  args - arguments to be printed
1406 
1407    Level: developer
1408 
1409    Note:
1410      You can change how help messages are printed by replacing the function pointer with a function that does not simply write to stdout.
1411 
1412       To use, write your own function, for example,
1413 $PetscErrorCode mypetschelpprintf(MPI_Comm comm,const char format[],....)
1414 ${
1415 $ PetscFunctionReturn(0);
1416 $}
1417 then do the assignment
1418 $    PetscHelpPrintf = mypetschelpprintf;
1419    You can do the assignment before `PetscInitialize()`.
1420 
1421   The default routine used is called `PetscHelpPrintfDefault()`.
1422 
1423 .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscErrorPrintf()`
1424 M*/
1425 PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1426 
1427 /*
1428      Defines PETSc profiling.
1429 */
1430 #include <petsclog.h>
1431 
1432 /*
1433       Simple PETSc parallel IO for ASCII printing
1434 */
1435 PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[], char[]);
1436 PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm, const char[], const char[], FILE **);
1437 PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm, FILE *);
1438 PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1439 PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1440 PETSC_EXTERN PetscErrorCode PetscSNPrintf(char *, size_t, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1441 PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char *, size_t, const char[], size_t *, ...) PETSC_ATTRIBUTE_FORMAT(3, 5);
1442 PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[], size_t, const char *, PetscInt, const PetscReal[]);
1443 
1444 PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1445 PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1446 PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1447 
1448 PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char *, size_t *);
1449 PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char *, char *);
1450 
1451 #if defined(PETSC_HAVE_POPEN)
1452 PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm, const char[], const char[], const char[], FILE **);
1453 PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm, FILE *);
1454 PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1455 #endif
1456 
1457 PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1458 PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1459 PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm, FILE *);
1460 PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm, FILE *, size_t, char[]);
1461 PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm, const char[], const char[], FILE **);
1462 PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm, const char[], const char[], FILE **);
1463 PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char *[]);
1464 
1465 PETSC_EXTERN PetscClassId   PETSC_CONTAINER_CLASSID;
1466 PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer, void **);
1467 PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer, void *);
1468 PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer *);
1469 PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm, PetscContainer *);
1470 PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void *));
1471 PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void *);
1472 
1473 /*
1474    For use in debuggers
1475 */
1476 PETSC_EXTERN PetscMPIInt    PetscGlobalRank;
1477 PETSC_EXTERN PetscMPIInt    PetscGlobalSize;
1478 PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt, const PetscInt[], PetscViewer);
1479 PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt, const PetscReal[], PetscViewer);
1480 PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt, const PetscScalar[], PetscViewer);
1481 
1482 #include <stddef.h>
1483 #include <string.h> /* for memcpy, memset */
1484 #include <stdlib.h>
1485 
1486 #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1487   #include <xmmintrin.h>
1488 #endif
1489 
1490 /*@C
1491    PetscMemmove - Copies n bytes, beginning at location b, to the space
1492    beginning at location a. Copying  between regions that overlap will
1493    take place correctly. Use `PetscMemcpy()` if the locations do not overlap
1494 
1495    Not Collective
1496 
1497    Input Parameters:
1498 +  b - pointer to initial memory space
1499 .  a - pointer to copy space
1500 -  n - length (in bytes) of space to copy
1501 
1502    Level: intermediate
1503 
1504    Notes:
1505    `PetscArraymove()` is preferred
1506 
1507    This routine is analogous to `memmove()`.
1508 
1509    Developers Note:
1510    This is inlined for performance
1511 
1512 .seealso: `PetscMemcpy()`, `PetscMemcmp()`, `PetscArrayzero()`, `PetscMemzero()`, `PetscArraycmp()`, `PetscArraycpy()`, `PetscStrallocpy()`,
1513           `PetscArraymove()`
1514 @*/
1515 static inline PetscErrorCode PetscMemmove(void *a, const void *b, size_t n)
1516 {
1517   PetscFunctionBegin;
1518   if (n > 0) {
1519     PetscCheck(a, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Trying to copy to null pointer");
1520     PetscCheck(b, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Trying to copy from a null pointer");
1521   }
1522 #if !defined(PETSC_HAVE_MEMMOVE)
1523   if (a < b) {
1524     if (a <= b - n) memcpy(a, b, n);
1525     else {
1526       memcpy(a, b, (int)(b - a));
1527       PetscMemmove(b, b + (int)(b - a), n - (int)(b - a));
1528     }
1529   } else {
1530     if (b <= a - n) memcpy(a, b, n);
1531     else {
1532       memcpy(b + n, b + (n - (int)(a - b)), (int)(a - b));
1533       PetscMemmove(a, b, n - (int)(a - b));
1534     }
1535   }
1536 #else
1537   memmove((char *)(a), (char *)(b), n);
1538 #endif
1539   PetscFunctionReturn(0);
1540 }
1541 
1542 /*@C
1543    PetscMemcpy - Copies n bytes, beginning at location b, to the space
1544    beginning at location a. The two memory regions CANNOT overlap, use
1545    `PetscMemmove()` in that case.
1546 
1547    Not Collective
1548 
1549    Input Parameters:
1550 +  b - pointer to initial memory space
1551 -  n - length (in bytes) of space to copy
1552 
1553    Output Parameter:
1554 .  a - pointer to copy space
1555 
1556    Level: intermediate
1557 
1558    Compile Option:
1559     `PETSC_PREFER_DCOPY_FOR_MEMCPY` will cause the BLAS dcopy() routine to be used
1560                                   for memory copies on double precision values.
1561     `PETSC_PREFER_COPY_FOR_MEMCPY` will cause C code to be used
1562                                   for memory copies on double precision values.
1563     `PETSC_PREFER_FORTRAN_FORMEMCPY` will cause Fortran code to be used
1564                                   for memory copies on double precision values.
1565 
1566    Notes:
1567    Prefer `PetscArraycpy()`
1568 
1569    This routine is analogous to `memcpy()`.
1570 
1571    Not available from Fortran
1572 
1573    Developer Note:
1574    This is inlined for fastest performance
1575 
1576 .seealso: `PetscMemzero()`, `PetscMemcmp()`, `PetscArrayzero()`, `PetscArraycmp()`, `PetscArraycpy()`, `PetscMemmove()`, `PetscStrallocpy()`
1577 
1578 @*/
1579 static inline PetscErrorCode PetscMemcpy(void *a, const void *b, size_t n)
1580 {
1581 #if defined(PETSC_USE_DEBUG)
1582   size_t al = (size_t)a, bl = (size_t)b;
1583   size_t nl = (size_t)n;
1584   PetscFunctionBegin;
1585   if (n > 0) {
1586     PetscCheck(b, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Trying to copy from a null pointer");
1587     PetscCheck(a, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Trying to copy to a null pointer");
1588   }
1589 #else
1590   PetscFunctionBegin;
1591 #endif
1592   if (a != b && n > 0) {
1593 #if defined(PETSC_USE_DEBUG)
1594     PetscCheck(!((al > bl && (al - bl) < nl) || (bl - al) < nl), PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Memory regions overlap: either use PetscMemmov()\n\
1595               or make sure your copy regions and lengths are correct. \n\
1596               Length (bytes) %zu first address %zu second address %zu",
1597                nl, al, bl);
1598 #endif
1599 #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1600     if (!(a % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1601       size_t len = n / sizeof(PetscScalar);
1602   #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1603       PetscBLASInt one = 1, blen;
1604       PetscCall(PetscBLASIntCast(len, &blen));
1605       PetscCallBLAS("BLAScopy", BLAScopy_(&blen, (PetscScalar *)b, &one, (PetscScalar *)a, &one));
1606   #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1607       fortrancopy_(&len, (PetscScalar *)b, (PetscScalar *)a);
1608   #else
1609       PetscScalar *x = (PetscScalar *)b, *y = (PetscScalar *)a;
1610       for (size_t i = 0; i < len; i++) y[i] = x[i];
1611   #endif
1612     } else {
1613       memcpy((char *)(a), (char *)(b), n);
1614     }
1615 #else
1616     memcpy((char *)(a), (char *)(b), n);
1617 #endif
1618   }
1619   PetscFunctionReturn(0);
1620 }
1621 
1622 /*@C
1623    PetscMemzero - Zeros the specified memory.
1624 
1625    Not Collective
1626 
1627    Input Parameters:
1628 +  a - pointer to beginning memory location
1629 -  n - length (in bytes) of memory to initialize
1630 
1631    Level: intermediate
1632 
1633    Compile Option:
1634    `PETSC_PREFER_BZERO` - on certain machines (the IBM RS6000) the bzero() routine happens
1635   to be faster than the memset() routine. This flag causes the bzero() routine to be used.
1636 
1637    Notes:
1638    Not available from Fortran
1639 
1640    Prefer `PetscArrayzero()`
1641 
1642    Developer Note:
1643    This is inlined for fastest performance
1644 
1645 .seealso: `PetscMemcpy()`, `PetscMemcmp()`, `PetscArrayzero()`, `PetscArraycmp()`, `PetscArraycpy()`, `PetscMemmove()`, `PetscStrallocpy()`
1646 @*/
1647 static inline PetscErrorCode PetscMemzero(void *a, size_t n)
1648 {
1649   if (n > 0) {
1650 #if defined(PETSC_USE_DEBUG)
1651     PetscCheck(a, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Trying to zero at a null pointer with %zu bytes", n);
1652 #endif
1653 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1654     if (!(((long)a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1655       size_t       i, len = n / sizeof(PetscScalar);
1656       PetscScalar *x = (PetscScalar *)a;
1657       for (i = 0; i < len; i++) x[i] = 0.0;
1658     } else {
1659 #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1660     if (!(((long)a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1661       PetscInt len = n / sizeof(PetscScalar);
1662       fortranzero_(&len, (PetscScalar *)a);
1663     } else {
1664 #endif
1665 #if defined(PETSC_PREFER_BZERO)
1666       bzero((char *)a, n);
1667 #else
1668       memset((char *)a, 0, n);
1669 #endif
1670 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1671     }
1672 #endif
1673   }
1674   return 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 /* the following petsc_static_inline require petscerror.h */
1924 
1925 /* Limit MPI to 32-bits */
1926 #define PETSC_MPI_INT_MAX 2147483647
1927 #define PETSC_MPI_INT_MIN -2147483647
1928 /* Limit BLAS to 32-bits */
1929 #define PETSC_BLAS_INT_MAX    2147483647
1930 #define PETSC_BLAS_INT_MIN    -2147483647
1931 #define PETSC_CUBLAS_INT_MAX  2147483647
1932 #define PETSC_HIPBLAS_INT_MAX 2147483647
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 #if !defined(PETSC_USE_64BIT_INDICES)
1959   if (a > PETSC_MAX_INT) {
1960     *b = 0;
1961     SETERRQ(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);
1962   }
1963 #endif
1964   *b = (PetscInt)(a);
1965   PetscFunctionReturn(0);
1966 }
1967 
1968 /*@C
1969     PetscCountCast - casts a `PetscCount` to a `PetscInt` (which may be 32 bits in size), generates an
1970          error if the `PetscInt` is not large enough to hold the number.
1971 
1972    Not Collective
1973 
1974    Input Parameter:
1975 .     a - the `PetscCount` value
1976 
1977    Output Parameter:
1978 .     b - the resulting `PetscInt` value
1979 
1980    Level: advanced
1981 
1982    Notes:
1983      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
1984 
1985    Not available from Fortran
1986 
1987 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscIntMultError()`, `PetscIntSumError()`, `PetscIntCast()`
1988 @*/
1989 static inline PetscErrorCode PetscCountCast(PetscCount a, PetscInt *b)
1990 {
1991   PetscFunctionBegin;
1992   if (sizeof(PetscCount) > sizeof(PetscInt) && a > PETSC_MAX_INT) {
1993     *b = 0;
1994     SETERRQ(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);
1995   }
1996   *b = (PetscInt)(a);
1997   PetscFunctionReturn(0);
1998 }
1999 
2000 /*@C
2001     PetscBLASIntCast - casts a `PetscInt` (which may be 64 bits in size) to a `PetscBLASInt` (which may be 32 bits in size), generates an
2002          error if the `PetscBLASInt` is not large enough to hold the number.
2003 
2004    Not Collective
2005 
2006    Input Parameter:
2007 .     a - the `PetscInt` value
2008 
2009    Output Parameter:
2010 .     b - the resulting `PetscBLASInt` value
2011 
2012    Level: advanced
2013 
2014    Notes:
2015    Not available from Fortran
2016 
2017    Errors if the integer is negative since PETSc calls to BLAS/LAPACK never need to cast negative integer inputs
2018 
2019 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscIntCast()`, `PetscCountCast()`
2020 @*/
2021 static inline PetscErrorCode PetscBLASIntCast(PetscInt a, PetscBLASInt *b)
2022 {
2023   PetscFunctionBegin;
2024   *b = 0;
2025   if (PetscDefined(USE_64BIT_INDICES) && !PetscDefined(HAVE_64BIT_BLAS_INDICES)) {
2026     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);
2027   }
2028   PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer to BLAS/LAPACK routine");
2029   *b = (PetscBLASInt)(a);
2030   PetscFunctionReturn(0);
2031 }
2032 
2033 /*@C
2034     PetscCuBLASIntCast - like `PetscBLASIntCast()`, but for `PetscCuBLASInt`.
2035 
2036    Not Collective
2037 
2038    Input Parameter:
2039 .     a - the `PetscInt` value
2040 
2041    Output Parameter:
2042 .     b - the resulting `PetscCuBLASInt` value
2043 
2044    Level: advanced
2045 
2046    Notes:
2047       Errors if the integer is negative since PETSc calls to cuBLAS and friends never need to cast negative integer inputs
2048 
2049 .seealso: `PetscCuBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscIntCast()`
2050 @*/
2051 static inline PetscErrorCode PetscCuBLASIntCast(PetscInt a, PetscCuBLASInt *b)
2052 {
2053   PetscFunctionBegin;
2054   *b = 0;
2055   if (PetscDefined(USE_64BIT_INDICES)) 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);
2056   PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer to cuBLAS routine");
2057   *b = (PetscCuBLASInt)(a);
2058   PetscFunctionReturn(0);
2059 }
2060 
2061 /*@C
2062     PetscHipBLASIntCast - like `PetscBLASIntCast()`, but for `PetscHipBLASInt`.
2063 
2064    Not Collective
2065 
2066    Input Parameter:
2067 .     a - the `PetscInt` value
2068 
2069    Output Parameter:
2070 .     b - the resulting `PetscHipBLASInt` value
2071 
2072    Level: advanced
2073 
2074    Notes:
2075       Errors if the integer is negative since PETSc calls to hipBLAS and friends never need to cast negative integer inputs
2076 
2077 .seealso: `PetscHipBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscIntCast()`
2078 @*/
2079 static inline PetscErrorCode PetscHipBLASIntCast(PetscInt a, PetscHipBLASInt *b)
2080 {
2081   PetscFunctionBegin;
2082   *b = 0;
2083   if (PetscDefined(USE_64BIT_INDICES)) 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);
2084   PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer to hipBLAS routine");
2085   *b = (PetscHipBLASInt)(a);
2086   PetscFunctionReturn(0);
2087 }
2088 
2089 /*@C
2090     PetscMPIIntCast - casts a PetscInt (which may be 64 bits in size) to a PetscMPIInt (which may be 32 bits in size), generates an
2091          error if the PetscMPIInt is not large enough to hold the number.
2092 
2093    Not Collective
2094 
2095    Input Parameter:
2096 .     a - the `PetscInt` value
2097 
2098    Output Parameter:
2099 .     b - the resulting `PetscMPIInt` value
2100 
2101    Level: advanced
2102 
2103    Not available from Fortran
2104 
2105 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscIntCast()`
2106 @*/
2107 static inline PetscErrorCode PetscMPIIntCast(PetscInt a, PetscMPIInt *b)
2108 {
2109   PetscFunctionBegin;
2110   *b = 0;
2111   if (PetscDefined(USE_64BIT_INDICES)) PetscCheck(a <= PETSC_MPI_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " is too big for MPI buffer length. We currently only support 32 bit integers", a);
2112   *b = (PetscMPIInt)(a);
2113   PetscFunctionReturn(0);
2114 }
2115 
2116 #define PetscInt64Mult(a, b) ((PetscInt64)(a)) * ((PetscInt64)(b))
2117 
2118 /*@C
2119 
2120    PetscRealIntMultTruncate - Computes the product of a positive `PetscReal` and a positive `PetscInt` and truncates the value to slightly less than the maximal possible value
2121 
2122    Not Collective
2123 
2124    Input Parameters:
2125 +     a - the `PetscReal` value
2126 -     b - the second value
2127 
2128    Returns:
2129       the result as a `PetscInt` value
2130 
2131    Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`
2132    Use `PetscIntMultTruncate()` to compute the product of two positive `PetscInt` and truncate to fit a `PetscInt`
2133    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`
2134 
2135    Developers Note:
2136    We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but requires many more checks.
2137 
2138    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2139 
2140    Not available from Fortran
2141 
2142    Level: advanced
2143 
2144 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()`
2145 @*/
2146 static inline PetscInt PetscRealIntMultTruncate(PetscReal a, PetscInt b)
2147 {
2148   PetscInt64 r = (PetscInt64)(a * (PetscReal)b);
2149   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2150   return (PetscInt)r;
2151 }
2152 
2153 /*@C
2154 
2155    PetscIntMultTruncate - Computes the product of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value
2156 
2157    Not Collective
2158 
2159    Input Parameters:
2160 +     a - the PetscInt value
2161 -     b - the second value
2162 
2163    Returns:
2164       the result as a `PetscInt` value
2165 
2166    Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`
2167    Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt`
2168    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`
2169 
2170    Not available from Fortran
2171 
2172    Developers Note:
2173    We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but requires many more checks.
2174 
2175    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2176 
2177    Level: advanced
2178 
2179 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()`
2180 @*/
2181 static inline PetscInt PetscIntMultTruncate(PetscInt a, PetscInt b)
2182 {
2183   PetscInt64 r = PetscInt64Mult(a, b);
2184   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2185   return (PetscInt)r;
2186 }
2187 
2188 /*@C
2189 
2190    PetscIntSumTruncate - Computes the sum of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value
2191 
2192    Not Collective
2193 
2194    Input Parameters:
2195 +     a - the `PetscInt` value
2196 -     b - the second value
2197 
2198    Returns:
2199      the result as a `PetscInt` value
2200 
2201    Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`
2202    Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt`
2203    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`
2204 
2205    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2206 
2207    Not available from Fortran
2208 
2209    Level: advanced
2210 
2211 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`
2212 @*/
2213 static inline PetscInt PetscIntSumTruncate(PetscInt a, PetscInt b)
2214 {
2215   PetscInt64 r = ((PetscInt64)a) + ((PetscInt64)b);
2216   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2217   return (PetscInt)r;
2218 }
2219 
2220 /*@C
2221 
2222    PetscIntMultError - Computes the product of two positive `PetscInt` and generates an error with overflow.
2223 
2224    Not Collective
2225 
2226    Input Parameters:
2227 +     a - the `PetscInt` value
2228 -     b - the second value
2229 
2230    Output Parameter:
2231 .     result - the result as a `PetscInt` value, or NULL if you do not want the result, you just want to check if it overflows
2232 
2233    Use `PetscInt64Mult()` to compute the product of two `PetscInt` and store in a `PetscInt64`
2234    Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt`
2235 
2236    Not available from Fortran
2237 
2238    Developers Note:
2239    We currently assume that `PetscInt` addition does not overflow, this is obviously wrong but requires many more checks.
2240 
2241    Level: advanced
2242 
2243 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscIntMult64()`, `PetscIntSumError()`
2244 @*/
2245 static inline PetscErrorCode PetscIntMultError(PetscInt a, PetscInt b, PetscInt *result)
2246 {
2247   PetscInt64 r;
2248 
2249   PetscFunctionBegin;
2250   r = PetscInt64Mult(a, b);
2251 #if !defined(PETSC_USE_64BIT_INDICES)
2252   PetscCheck(r <= PETSC_MAX_INT, PETSC_COMM_SELF, PETSC_ERR_SUP, "Product of two integers %d %d 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);
2253 #endif
2254   if (result) *result = (PetscInt)r;
2255   PetscFunctionReturn(0);
2256 }
2257 
2258 /*@C
2259 
2260    PetscIntSumError - Computes the sum of two positive `PetscInt` and generates an error with overflow.
2261 
2262    Not Collective
2263 
2264    Input Parameters:
2265 +     a - the `PetscInt` value
2266 -     b - the second value
2267 
2268    Output Parameter:
2269 .     c - the result as a `PetscInt` value,  or NULL if you do not want the result, you just want to check if it overflows
2270 
2271    Use `PetscInt64Mult()` to compute the product of two 32 bit PetscInt and store in a `PetscInt64`
2272    Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt`
2273 
2274    Not available from Fortran
2275 
2276    Level: advanced
2277 
2278 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`
2279 @*/
2280 static inline PetscErrorCode PetscIntSumError(PetscInt a, PetscInt b, PetscInt *result)
2281 {
2282   PetscInt64 r;
2283 
2284   PetscFunctionBegin;
2285   r = ((PetscInt64)a) + ((PetscInt64)b);
2286 #if !defined(PETSC_USE_64BIT_INDICES)
2287   PetscCheck(r <= PETSC_MAX_INT, PETSC_COMM_SELF, PETSC_ERR_SUP, "Sum of two integers %d %d 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);
2288 #endif
2289   if (result) *result = (PetscInt)r;
2290   PetscFunctionReturn(0);
2291 }
2292 
2293 /*
2294      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2295 */
2296 #if defined(hz)
2297   #undef hz
2298 #endif
2299 
2300 #include <limits.h>
2301 
2302 /* The number of bits in a byte */
2303 
2304 #define PETSC_BITS_PER_BYTE CHAR_BIT
2305 
2306 #if defined(PETSC_HAVE_SYS_TYPES_H)
2307   #include <sys/types.h>
2308 #endif
2309 
2310 /*MC
2311 
2312     PETSC_VERSION - This manual page provides information about how PETSc documents and uses its version information. This information is available to both C/C++
2313                     and Fortran compilers when petscsys.h is included.
2314 
2315     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>
2316 
2317     The complete version number is given as the triple  PETSC_VERSION_MAJOR.PETSC_VERSION_MINOR.PETSC_VERSION_SUBMINOR (in short hand x.y.z)
2318 
2319     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
2320     where only a change in the major version number (x) indicates a change in the API.
2321 
2322     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
2323 
2324     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),
2325     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
2326     version number (x.y.z).
2327 
2328     `PETSC_RELEASE_DATE` is the date the x.y version was released (i.e. the version before any patch releases)
2329 
2330     `PETSC_VERSION_DATE` is the date the x.y.z version was released
2331 
2332     `PETSC_VERSION_GIT` is the last git commit to the repository given in the form vx.y.z-wwwww
2333 
2334     `PETSC_VERSION_DATE_GIT` is the date of the last git commit to the repository
2335 
2336     `PETSC_VERSION_()` is deprecated and will eventually be removed.
2337 
2338     Level: intermediate
2339 
2340 M*/
2341 
2342 PETSC_EXTERN PetscErrorCode PetscGetArchType(char[], size_t);
2343 PETSC_EXTERN PetscErrorCode PetscGetHostName(char[], size_t);
2344 PETSC_EXTERN PetscErrorCode PetscGetUserName(char[], size_t);
2345 PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[], size_t);
2346 PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2347 PETSC_EXTERN PetscErrorCode PetscGetDate(char[], size_t);
2348 PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2349 PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt *, PetscInt *, PetscInt *, PetscInt *);
2350 
2351 PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt, const PetscInt[], PetscBool *);
2352 PETSC_EXTERN PetscErrorCode PetscSortedInt64(PetscInt, const PetscInt64[], PetscBool *);
2353 PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt, const PetscMPIInt[], PetscBool *);
2354 PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt, const PetscReal[], PetscBool *);
2355 PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt, PetscInt[]);
2356 PETSC_EXTERN PetscErrorCode PetscSortInt64(PetscInt, PetscInt64[]);
2357 PETSC_EXTERN PetscErrorCode PetscSortCount(PetscInt, PetscCount[]);
2358 PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt, PetscInt[]);
2359 PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *, PetscInt[]);
2360 PETSC_EXTERN PetscErrorCode PetscSortedCheckDupsInt(PetscInt, const PetscInt[], PetscBool *);
2361 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt *, PetscInt[]);
2362 PETSC_EXTERN PetscErrorCode PetscCheckDupsInt(PetscInt, const PetscInt[], PetscBool *);
2363 PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt *);
2364 PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt *);
2365 PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt, const PetscInt[], PetscInt[]);
2366 PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt, const char *[], PetscInt[]);
2367 PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt, PetscInt[], PetscInt[]);
2368 PETSC_EXTERN PetscErrorCode PetscSortIntWithCountArray(PetscCount, PetscInt[], PetscCount[]);
2369 PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt, PetscInt[], PetscInt[], PetscInt[]);
2370 PETSC_EXTERN PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount, PetscInt[], PetscInt[], PetscCount[]);
2371 PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt, PetscMPIInt[]);
2372 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *, PetscMPIInt[]);
2373 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt, PetscMPIInt[], PetscMPIInt[]);
2374 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt, PetscMPIInt[], PetscInt[]);
2375 PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt, PetscInt[], PetscScalar[]);
2376 PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt, PetscInt[], void *, size_t, void *);
2377 PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt, PetscReal[]);
2378 PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt, PetscReal[], PetscInt[]);
2379 PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt, const PetscReal[], PetscInt[]);
2380 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt *, PetscReal[]);
2381 PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal, PetscInt, const PetscReal[], PetscReal, PetscInt *);
2382 PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt, PetscInt, PetscScalar[], PetscInt[]);
2383 PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt, PetscInt, PetscReal[], PetscInt[]);
2384 PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt, const PetscBool[], const PetscInt[], PetscInt *, PetscInt **, PetscInt **, PetscInt **, PetscInt **);
2385 PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt, const PetscInt[], const PetscInt[], PetscInt, const PetscInt[], const PetscInt[], PetscInt *, PetscInt **, PetscInt **);
2386 PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscInt *, PetscInt **);
2387 PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt, const PetscMPIInt[], PetscInt, const PetscMPIInt[], PetscInt *, PetscMPIInt **);
2388 PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *);
2389 
2390 PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt, void *, size_t, int (*)(const void *, const void *, void *), void *);
2391 PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt, PetscInt[]);
2392 PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt, PetscMPIInt[]);
2393 PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt, PetscReal[]);
2394 PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt, void *, size_t, void *, size_t, int (*)(const void *, const void *, void *), void *);
2395 PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt, PetscInt[], PetscInt[]);
2396 PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt, PetscMPIInt[], PetscMPIInt[]);
2397 PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt, PetscReal[], PetscInt[]);
2398 
2399 PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2400 PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[], size_t);
2401 
2402 /*J
2403     PetscRandomType - String with the name of a PETSc randomizer
2404 
2405    Level: beginner
2406 
2407    Notes:
2408    To use `PETSCSPRNG` or `PETSCRANDOM123` you must have ./configure PETSc
2409    with the option --download-sprng or --download-random123
2410 
2411 .seealso: `PetscRandomSetType()`, `PetscRandom`, `PetscRandomCreate()`
2412 J*/
2413 typedef const char *PetscRandomType;
2414 #define PETSCRAND      "rand"
2415 #define PETSCRAND48    "rand48"
2416 #define PETSCSPRNG     "sprng"
2417 #define PETSCRANDER48  "rander48"
2418 #define PETSCRANDOM123 "random123"
2419 #define PETSCCURAND    "curand"
2420 
2421 /* Logging support */
2422 PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;
2423 
2424 PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);
2425 
2426 /* Dynamic creation and loading functions */
2427 PETSC_EXTERN PetscFunctionList PetscRandomList;
2428 
2429 PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[], PetscErrorCode (*)(PetscRandom));
2430 PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2431 PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2432 PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType *);
2433 PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom, PetscObject, const char[]);
2434 PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom, PetscViewer);
2435 
2436 PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm, PetscRandom *);
2437 PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom, PetscScalar *);
2438 PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom, PetscReal *);
2439 PETSC_EXTERN PetscErrorCode PetscRandomGetValues(PetscRandom, PetscInt, PetscScalar *);
2440 PETSC_EXTERN PetscErrorCode PetscRandomGetValuesReal(PetscRandom, PetscInt, PetscReal *);
2441 PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom, PetscScalar *, PetscScalar *);
2442 PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom, PetscScalar, PetscScalar);
2443 PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom, unsigned long);
2444 PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom, unsigned long *);
2445 PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2446 PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom *);
2447 
2448 PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[], char[], size_t);
2449 PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[], char[], size_t);
2450 PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[], size_t);
2451 PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[], char[]);
2452 PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[], size_t);
2453 PETSC_EXTERN PetscErrorCode PetscTestFile(const char[], char, PetscBool *);
2454 PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[], char, PetscBool *);
2455 PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2456 PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2457 PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);
2458 
2459 static inline PetscBool PetscBinaryBigEndian(void)
2460 {
2461   long _petsc_v = 1;
2462   return ((char *)&_petsc_v)[0] ? PETSC_FALSE : PETSC_TRUE;
2463 }
2464 
2465 PETSC_EXTERN PetscErrorCode PetscBinaryRead(int, void *, PetscInt, PetscInt *, PetscDataType);
2466 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm, int, void *, PetscInt, PetscInt *, PetscDataType);
2467 PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int, const void *, PetscInt, PetscDataType);
2468 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm, int, const void *, PetscInt, PetscDataType);
2469 PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[], PetscFileMode, int *);
2470 PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2471 PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm, PetscBool *);
2472 PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm, PetscBool *);
2473 PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm, char[], size_t);
2474 PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm, const char[], char[], size_t, PetscBool *);
2475 PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm, const char[], char[], size_t, PetscBool *);
2476 #if defined(PETSC_USE_SOCKET_VIEWER)
2477 PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[], int, int *);
2478 #endif
2479 
2480 PETSC_EXTERN PetscErrorCode PetscBinarySeek(int, off_t, PetscBinarySeekType, off_t *);
2481 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm, int, off_t, PetscBinarySeekType, off_t *);
2482 PETSC_EXTERN PetscErrorCode PetscByteSwap(void *, PetscDataType, PetscInt);
2483 
2484 PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2485 PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[], PetscBool);
2486 PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2487 PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char *);
2488 PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2489 PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2490 PETSC_EXTERN PetscErrorCode PetscWaitOnError(void);
2491 
2492 PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt *);
2493 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], PetscMPIInt **, PetscMPIInt **);
2494 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt **, PetscMPIInt **, PetscMPIInt **);
2495 PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscInt ***, MPI_Request **);
2496 PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscScalar ***, MPI_Request **);
2497 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);
2498 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);
2499 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);
2500 
2501 PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2502 PETSC_EXTERN PetscErrorCode    PetscCommBuildTwoSidedSetType(MPI_Comm, PetscBuildTwoSidedType);
2503 PETSC_EXTERN PetscErrorCode    PetscCommBuildTwoSidedGetType(MPI_Comm, PetscBuildTwoSidedType *);
2504 
2505 PETSC_EXTERN PetscErrorCode PetscSSEIsEnabled(MPI_Comm, PetscBool *, PetscBool *);
2506 
2507 PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);
2508 
2509 PETSC_EXTERN const char *const PetscSubcommTypes[];
2510 
2511 struct _n_PetscSubcomm {
2512   MPI_Comm         parent;    /* parent communicator */
2513   MPI_Comm         dupparent; /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2514   MPI_Comm         child;     /* the sub-communicator */
2515   PetscMPIInt      n;         /* num of subcommunicators under the parent communicator */
2516   PetscMPIInt      color;     /* color of processors belong to this communicator */
2517   PetscMPIInt     *subsize;   /* size of subcommunicator[color] */
2518   PetscSubcommType type;
2519   char            *subcommprefix;
2520 };
2521 
2522 static inline MPI_Comm PetscSubcommParent(PetscSubcomm scomm)
2523 {
2524   return scomm->parent;
2525 }
2526 static inline MPI_Comm PetscSubcommChild(PetscSubcomm scomm)
2527 {
2528   return scomm->child;
2529 }
2530 static inline MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm)
2531 {
2532   return scomm->dupparent;
2533 }
2534 PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm, PetscSubcomm *);
2535 PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm *);
2536 PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm, PetscInt);
2537 PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm, PetscSubcommType);
2538 PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm, PetscMPIInt, PetscMPIInt);
2539 PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm, PetscViewer);
2540 PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2541 PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm, const char[]);
2542 PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm, MPI_Comm *);
2543 PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm, MPI_Comm *);
2544 PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm, MPI_Comm *);
2545 
2546 PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt, PetscHeap *);
2547 PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap, PetscInt, PetscInt);
2548 PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap, PetscInt *, PetscInt *);
2549 PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap, PetscInt *, PetscInt *);
2550 PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap, PetscInt, PetscInt);
2551 PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2552 PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap *);
2553 PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap, PetscViewer);
2554 
2555 PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2556 PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm, PetscShmComm *);
2557 PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm, PetscMPIInt, PetscMPIInt *);
2558 PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm, PetscMPIInt, PetscMPIInt *);
2559 PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm, MPI_Comm *);
2560 
2561 /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2562 PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm, PetscInt, PetscOmpCtrl *);
2563 PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl, MPI_Comm *, MPI_Comm *, PetscBool *);
2564 PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl *);
2565 PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2566 PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2567 PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);
2568 
2569 PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t, size_t, PetscSegBuffer *);
2570 PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer *);
2571 PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer, size_t, void *);
2572 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer, void *);
2573 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer, void *);
2574 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer, void *);
2575 PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer, size_t *);
2576 PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer, size_t);
2577 
2578 /* Type-safe wrapper to encourage use of PETSC_RESTRICT. Does not use PetscFunctionBegin because the error handling
2579  * prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as
2580  * possible. */
2581 static inline PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg, size_t count, PetscInt *PETSC_RESTRICT *slot)
2582 {
2583   return PetscSegBufferGet(seg, count, (void **)slot);
2584 }
2585 
2586 extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2587 PETSC_EXTERN PetscErrorCode    PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted *);
2588 PETSC_EXTERN PetscErrorCode    PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted *);
2589 PETSC_EXTERN PetscErrorCode    PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted, const char *, const char *, PetscBool *);
2590 
2591 #include <stdarg.h>
2592 PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char *, size_t, const char[], size_t *, va_list);
2593 PETSC_EXTERN                PetscErrorCode (*PetscVFPrintf)(FILE *, const char[], va_list);
2594 
2595 PETSC_EXTERN PetscSegBuffer PetscCitationsList;
2596 
2597 /*@C
2598       PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code.
2599 
2600      Not Collective - only what is registered on rank 0 of `PETSC_COMM_WORLD` will be printed
2601 
2602      Input Parameters:
2603 +      cite - the bibtex item, formatted to displayed on multiple lines nicely
2604 -      set - a boolean variable initially set to `PETSC_FALSE`; this is used to insure only a single registration of the citation
2605 
2606      Options Database: Key
2607 .     -citations [filename]   - print out the bibtex entries for the given computation
2608 
2609      Level: intermediate
2610 
2611      Fortran Note:
2612      Not available from Fortran
2613 
2614 @*/
2615 static inline PetscErrorCode PetscCitationsRegister(const char cit[], PetscBool *set)
2616 {
2617   size_t len;
2618   char  *vstring;
2619 
2620   PetscFunctionBegin;
2621   if (set && *set) PetscFunctionReturn(0);
2622   PetscCall(PetscStrlen(cit, &len));
2623   PetscCall(PetscSegBufferGet(PetscCitationsList, len, &vstring));
2624   PetscCall(PetscArraycpy(vstring, cit, len));
2625   if (set) *set = PETSC_TRUE;
2626   PetscFunctionReturn(0);
2627 }
2628 
2629 PETSC_EXTERN                PETSC_DEPRECATED_FUNCTION("Google has discontinued its URL shortener service") PetscErrorCode PetscURLShorten(const char[], char[], size_t c);
2630 PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm, char[], char[], size_t);
2631 PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm, const char[], char[], size_t);
2632 PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm, const char[], const char[]);
2633 
2634 PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm, char[], char[], size_t);
2635 PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm, const char[], char[], char[], size_t);
2636 
2637 PETSC_EXTERN PetscErrorCode PetscGlobusGetTransfers(MPI_Comm, const char[], char[], size_t);
2638 
2639 PETSC_EXTERN PetscErrorCode PetscTextBelt(MPI_Comm, const char[], const char[], PetscBool *);
2640 PETSC_EXTERN PetscErrorCode PetscTellMyCell(MPI_Comm, const char[], const char[], PetscBool *);
2641 
2642 PETSC_EXTERN PetscErrorCode PetscPullJSONValue(const char[], const char[], char[], size_t, PetscBool *);
2643 PETSC_EXTERN PetscErrorCode PetscPushJSONValue(char[], const char[], const char[], size_t);
2644 
2645 #if defined(PETSC_USE_DEBUG)
2646 static inline unsigned int PetscStrHash(const char *str)
2647 {
2648   unsigned int c, hash = 5381;
2649 
2650   while ((c = (unsigned int)*str++)) hash = ((hash << 5) + hash) + c; /* hash * 33 + c */
2651   return hash;
2652 }
2653 
2654   /*MC
2655    MPIU_Allreduce - a PETSc replacement for `MPI_Allreduce()` that tries to determine if the call from all the MPI ranks occur from the
2656                     same place in the PETSc code. This helps to detect bugs where different MPI ranks follow different code paths
2657                     resulting in inconsistent and incorrect calls to `MPI_Allreduce()`.
2658 
2659    Collective
2660 
2661    Synopsis:
2662      #include <petscsys.h>
2663      PetscErrorCode MPIU_Allreduce(void *indata,void *outdata,PetscMPIInt count,MPI_Datatype datatype, MPI_Op op, MPI_Comm comm);
2664 
2665    Input Parameters:
2666 +  a - pointer to the input data to be reduced
2667 .  c - the number of MPI data items in a and b
2668 .  d - the MPI datatype, for example `MPI_INT`
2669 .  e - the MPI operation, for example `MPI_SUM`
2670 -  fcomm - the MPI communicator on which the operation occurs
2671 
2672    Output Parameter:
2673 .  b - the reduced values
2674 
2675    Notes:
2676      In optimized mode this directly calls `MPI_Allreduce()`
2677 
2678      This is defined as a macro that can return error codes internally so it cannot be used in a subroutine that returns void.
2679 
2680      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
2681 
2682    Level: developer
2683 
2684 .seealso: `MPI_Allreduce()`
2685 M*/
2686   #define MPIU_Allreduce(a, b, c, d, e, fcomm) \
2687     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)));)
2688 #else
2689   #define MPIU_Allreduce(a, b, c, d, e, fcomm) PetscMacroReturnStandard(PetscCallMPI(MPI_Allreduce((a), (b), (c), d, e, (fcomm))))
2690 #endif
2691 
2692 #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
2693 PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint, PetscMPIInt, MPI_Info, MPI_Comm, void *, MPI_Win *);
2694 PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win, PetscMPIInt, MPI_Aint *, PetscMPIInt *, void *);
2695 #endif
2696 
2697 /* this is a vile hack */
2698 #if defined(PETSC_HAVE_NECMPI)
2699   #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)
2700     #define MPI_Type_free(a) (*(a) = MPI_DATATYPE_NULL, 0);
2701   #endif
2702 #endif
2703 
2704 /*
2705     List of external packages and queries on it
2706 */
2707 PETSC_EXTERN PetscErrorCode PetscHasExternalPackage(const char[], PetscBool *);
2708 
2709 /*
2710  OpenMP support
2711 */
2712 #if defined(_OPENMP)
2713   #define PetscPragmaOMP(...) _Pragma(PetscStringize(omp __VA_ARGS__))
2714 #else // no OpenMP so no threads
2715   #define PetscPragmaOMP(...)
2716 #endif
2717 
2718 /* this cannot go here because it may be in a different shared library */
2719 PETSC_EXTERN PetscErrorCode PCMPIServerBegin(void);
2720 PETSC_EXTERN PetscErrorCode PCMPIServerEnd(void);
2721 PETSC_EXTERN PetscErrorCode PCMPICommsDestroy(void);
2722 #endif
2723