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