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