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