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