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