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