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