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