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