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