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