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