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