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