xref: /petsc/include/petscsys.h (revision 1795a4d16c893ec2fc06bbbc6c5ce592a2de75d4)
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     Defines some elementary mathematics functions and constants.
335 */
336 #include <petscmath.h>
337 
338 /*
339     Declare extern C stuff after including external header files
340 */
341 
342 
343 /*
344        Basic PETSc constants
345 */
346 
347 /*E
348     PetscBool  - Logical variable. Actually an int in C and a logical in Fortran.
349 
350    Level: beginner
351 
352    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
353       boolean values. It is not easy to have a simple macro that that will work properly in all circumstances with all three mechanisms.
354 
355 E*/
356 typedef enum { PETSC_FALSE,PETSC_TRUE } PetscBool;
357 PETSC_EXTERN const char *const PetscBools[];
358 PETSC_EXTERN MPI_Datatype MPIU_BOOL PetscAttrMPITypeTag(PetscBool);
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_DOUBLE_PRECISION.
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: Memory is always allocated at least double aligned
520 
521           If you request memory of zero size it will allocate no space and assign the pointer to 0; PetscFree() will
522           properly handle not freeing the null pointer.
523 
524 .seealso: PetscFree(), PetscNew()
525 
526   Concepts: memory allocation
527 
528 M*/
529 #define PetscMalloc(a,b)  ((a != 0) ? (*PetscTrMalloc)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__,(void**)(b)) : (*(b) = 0,0) )
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   (PetscMalloc6(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(*(r6),(m6)*sizeof(**(r6))))
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(struct type,((type *))result)
1021 
1022    Not Collective
1023 
1024    Input Parameter:
1025 .  type - structure name of space to be allocated. Memory of size sizeof(type) is allocated
1026 
1027    Output Parameter:
1028 .  result - memory allocated
1029 
1030    Level: beginner
1031 
1032 .seealso: PetscFree(), PetscMalloc(), PetscNewLog()
1033 
1034   Concepts: memory allocation
1035 
1036 M*/
1037 #define PetscNew(A,b)      (PetscMalloc(sizeof(A),(b)) || PetscMemzero(*(b),sizeof(A)))
1038 
1039 /*MC
1040    PetscNewLog - Allocates memory of a particular type, zeros the memory! Aligned to PETSC_MEMALIGN. Associates the memory allocated
1041          with the given object using PetscLogObjectMemory().
1042 
1043    Synopsis:
1044     #include "petscsys.h"
1045    PetscErrorCode PetscNewLog(PetscObject obj,struct type,((type *))result)
1046 
1047    Not Collective
1048 
1049    Input Parameter:
1050 +  obj - object memory is logged to
1051 -  type - structure name of space to be allocated. Memory of size sizeof(type) is allocated
1052 
1053    Output Parameter:
1054 .  result - memory allocated
1055 
1056    Level: developer
1057 
1058 .seealso: PetscFree(), PetscMalloc(), PetscNew(), PetscLogObjectMemory()
1059 
1060   Concepts: memory allocation
1061 
1062 M*/
1063 #define PetscNewLog(o,A,b) (PetscNew(A,b) || ((o) ? PetscLogObjectMemory((PetscObject)o,sizeof(A)) : 0))
1064 
1065 /*MC
1066    PetscFree - Frees memory
1067 
1068    Synopsis:
1069     #include "petscsys.h"
1070    PetscErrorCode PetscFree(void *memory)
1071 
1072    Not Collective
1073 
1074    Input Parameter:
1075 .   memory - memory to free (the pointer is ALWAYS set to 0 upon sucess)
1076 
1077    Level: beginner
1078 
1079    Notes: Memory must have been obtained with PetscNew() or PetscMalloc()
1080 
1081 .seealso: PetscNew(), PetscMalloc(), PetscFreeVoid()
1082 
1083   Concepts: memory allocation
1084 
1085 M*/
1086 #define PetscFree(a)   ((a) && ((*PetscTrFree)((void*)(a),__LINE__,PETSC_FUNCTION_NAME,__FILE__) || ((a) = 0,0)))
1087 
1088 /*MC
1089    PetscFreeVoid - Frees memory
1090 
1091    Synopsis:
1092     #include "petscsys.h"
1093    void PetscFreeVoid(void *memory)
1094 
1095    Not Collective
1096 
1097    Input Parameter:
1098 .   memory - memory to free
1099 
1100    Level: beginner
1101 
1102    Notes: This is different from PetscFree() in that no error code is returned
1103 
1104 .seealso: PetscFree(), PetscNew(), PetscMalloc()
1105 
1106   Concepts: memory allocation
1107 
1108 M*/
1109 #define PetscFreeVoid(a) ((*PetscTrFree)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__),(a) = 0)
1110 
1111 
1112 /*MC
1113    PetscFree2 - Frees 2 chunks of memory obtained with PetscMalloc2()
1114 
1115    Synopsis:
1116     #include "petscsys.h"
1117    PetscErrorCode PetscFree2(void *memory1,void *memory2)
1118 
1119    Not Collective
1120 
1121    Input Parameter:
1122 +   memory1 - memory to free
1123 -   memory2 - 2nd memory to free
1124 
1125    Level: developer
1126 
1127    Notes: Memory must have been obtained with PetscMalloc2()
1128 
1129 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree()
1130 
1131   Concepts: memory allocation
1132 
1133 M*/
1134 #if defined(PETSC_USE_DEBUG)
1135 #define PetscFree2(m1,m2)   (PetscFree(m2) || PetscFree(m1))
1136 #else
1137 #define PetscFree2(m1,m2)   ((m2)=0, PetscFree(m1))
1138 #endif
1139 
1140 /*MC
1141    PetscFree3 - Frees 3 chunks of memory obtained with PetscMalloc3()
1142 
1143    Synopsis:
1144     #include "petscsys.h"
1145    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)
1146 
1147    Not Collective
1148 
1149    Input Parameter:
1150 +   memory1 - memory to free
1151 .   memory2 - 2nd memory to free
1152 -   memory3 - 3rd memory to free
1153 
1154    Level: developer
1155 
1156    Notes: Memory must have been obtained with PetscMalloc3()
1157 
1158 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3()
1159 
1160   Concepts: memory allocation
1161 
1162 M*/
1163 #if defined(PETSC_USE_DEBUG)
1164 #define PetscFree3(m1,m2,m3)   (PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1165 #else
1166 #define PetscFree3(m1,m2,m3)   ((m3)=0,(m2)=0,PetscFree(m1))
1167 #endif
1168 
1169 /*MC
1170    PetscFree4 - Frees 4 chunks of memory obtained with PetscMalloc4()
1171 
1172    Synopsis:
1173     #include "petscsys.h"
1174    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)
1175 
1176    Not Collective
1177 
1178    Input Parameter:
1179 +   m1 - memory to free
1180 .   m2 - 2nd memory to free
1181 .   m3 - 3rd memory to free
1182 -   m4 - 4th memory to free
1183 
1184    Level: developer
1185 
1186    Notes: Memory must have been obtained with PetscMalloc4()
1187 
1188 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4()
1189 
1190   Concepts: memory allocation
1191 
1192 M*/
1193 #if defined(PETSC_USE_DEBUG)
1194 #define PetscFree4(m1,m2,m3,m4)   (PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1195 #else
1196 #define PetscFree4(m1,m2,m3,m4)   ((m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
1197 #endif
1198 
1199 /*MC
1200    PetscFree5 - Frees 5 chunks of memory obtained with PetscMalloc5()
1201 
1202    Synopsis:
1203     #include "petscsys.h"
1204    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)
1205 
1206    Not Collective
1207 
1208    Input Parameter:
1209 +   m1 - memory to free
1210 .   m2 - 2nd memory to free
1211 .   m3 - 3rd memory to free
1212 .   m4 - 4th memory to free
1213 -   m5 - 5th memory to free
1214 
1215    Level: developer
1216 
1217    Notes: Memory must have been obtained with PetscMalloc5()
1218 
1219 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5()
1220 
1221   Concepts: memory allocation
1222 
1223 M*/
1224 #if defined(PETSC_USE_DEBUG)
1225 #define PetscFree5(m1,m2,m3,m4,m5)   (PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1226 #else
1227 #define PetscFree5(m1,m2,m3,m4,m5)   ((m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
1228 #endif
1229 
1230 
1231 /*MC
1232    PetscFree6 - Frees 6 chunks of memory obtained with PetscMalloc6()
1233 
1234    Synopsis:
1235     #include "petscsys.h"
1236    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)
1237 
1238    Not Collective
1239 
1240    Input Parameter:
1241 +   m1 - memory to free
1242 .   m2 - 2nd memory to free
1243 .   m3 - 3rd memory to free
1244 .   m4 - 4th memory to free
1245 .   m5 - 5th memory to free
1246 -   m6 - 6th memory to free
1247 
1248 
1249    Level: developer
1250 
1251    Notes: Memory must have been obtained with PetscMalloc6()
1252 
1253 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6()
1254 
1255   Concepts: memory allocation
1256 
1257 M*/
1258 #if defined(PETSC_USE_DEBUG)
1259 #define PetscFree6(m1,m2,m3,m4,m5,m6)   (PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1260 #else
1261 #define PetscFree6(m1,m2,m3,m4,m5,m6)   ((m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
1262 #endif
1263 
1264 /*MC
1265    PetscFree7 - Frees 7 chunks of memory obtained with PetscMalloc7()
1266 
1267    Synopsis:
1268     #include "petscsys.h"
1269    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)
1270 
1271    Not Collective
1272 
1273    Input Parameter:
1274 +   m1 - memory to free
1275 .   m2 - 2nd memory to free
1276 .   m3 - 3rd memory to free
1277 .   m4 - 4th memory to free
1278 .   m5 - 5th memory to free
1279 .   m6 - 6th memory to free
1280 -   m7 - 7th memory to free
1281 
1282 
1283    Level: developer
1284 
1285    Notes: Memory must have been obtained with PetscMalloc7()
1286 
1287 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1288           PetscMalloc7()
1289 
1290   Concepts: memory allocation
1291 
1292 M*/
1293 #if defined(PETSC_USE_DEBUG)
1294 #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   (PetscFree(m7) || PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1295 #else
1296 #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   ((m7)=0,(m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
1297 #endif
1298 
1299 PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,int,const char[],const char[],void**);
1300 PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1301 PETSC_EXTERN PetscErrorCode PetscMallocSet(PetscErrorCode (*)(size_t,int,const char[],const char[],void**),PetscErrorCode (*)(void*,int,const char[],const char[]));
1302 PETSC_EXTERN PetscErrorCode PetscMallocClear(void);
1303 
1304 /*
1305     PetscLogDouble variables are used to contain double precision numbers
1306   that are not used in the numerical computations, but rather in logging,
1307   timing etc.
1308 */
1309 typedef double PetscLogDouble;
1310 #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE
1311 
1312 /*
1313    Routines for tracing memory corruption/bleeding with default PETSc  memory allocation
1314 */
1315 PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1316 PETSC_EXTERN PetscErrorCode PetscMallocDumpLog(FILE *);
1317 PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1318 PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1319 PETSC_EXTERN PetscErrorCode PetscMallocDebug(PetscBool);
1320 PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*);
1321 PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1322 PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLog(void);
1323 PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLogThreshold(PetscLogDouble);
1324 PETSC_EXTERN PetscErrorCode PetscMallocGetDumpLog(PetscBool*);
1325 
1326 /*E
1327     PetscDataType - Used for handling different basic data types.
1328 
1329    Level: beginner
1330 
1331    Developer comment: It would be nice if we could always just use MPI Datatypes, why can we not?
1332 
1333 .seealso: PetscBinaryRead(), PetscBinaryWrite(), PetscDataTypeToMPIDataType(),
1334           PetscDataTypeGetSize()
1335 
1336 E*/
1337 typedef enum {PETSC_INT = 0,PETSC_DOUBLE = 1,PETSC_COMPLEX = 2, PETSC_LONG = 3 ,PETSC_SHORT = 4,PETSC_FLOAT = 5,
1338               PETSC_CHAR = 6,PETSC_BIT_LOGICAL = 7,PETSC_ENUM = 8,PETSC_BOOL=9, PETSC___FLOAT128 = 10,PETSC_OBJECT = 11, PETSC_FUNCTION = 12} PetscDataType;
1339 PETSC_EXTERN const char *const PetscDataTypes[];
1340 
1341 #if defined(PETSC_USE_COMPLEX)
1342 #define  PETSC_SCALAR  PETSC_COMPLEX
1343 #else
1344 #if defined(PETSC_USE_REAL_SINGLE)
1345 #define  PETSC_SCALAR  PETSC_FLOAT
1346 #elif defined(PETSC_USE_REAL___FLOAT128)
1347 #define  PETSC_SCALAR  PETSC___FLOAT128
1348 #else
1349 #define  PETSC_SCALAR  PETSC_DOUBLE
1350 #endif
1351 #endif
1352 #if defined(PETSC_USE_REAL_SINGLE)
1353 #define  PETSC_REAL  PETSC_FLOAT
1354 #elif defined(PETSC_USE_REAL___FLOAT128)
1355 #define  PETSC_REAL  PETSC___FLOAT128
1356 #else
1357 #define  PETSC_REAL  PETSC_DOUBLE
1358 #endif
1359 #define  PETSC_FORTRANADDR  PETSC_LONG
1360 
1361 PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1362 PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1363 PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1364 PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);
1365 
1366 /*
1367     Basic memory and string operations. These are usually simple wrappers
1368    around the basic Unix system calls, but a few of them have additional
1369    functionality and/or error checking.
1370 */
1371 PETSC_EXTERN PetscErrorCode PetscBitMemcpy(void*,PetscInt,const void*,PetscInt,PetscInt,PetscDataType);
1372 PETSC_EXTERN PetscErrorCode PetscMemmove(void*,void *,size_t);
1373 PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1374 PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1375 PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***);
1376 PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1377 PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool  *);
1378 PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool  *);
1379 PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1380 PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1381 PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1382 PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1383 PETSC_EXTERN PetscErrorCode PetscStrncat(char[],const char[],size_t);
1384 PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1385 PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1386 PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1387 PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1388 PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1389 PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1390 PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1391 PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1392 PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*);
1393 PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1394 PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1395 PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1396 PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1397 PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);
1398 
1399 PETSC_EXTERN void PetscStrcmpNoError(const char[],const char[],PetscBool  *);
1400 
1401 /*S
1402     PetscToken - 'Token' used for managing tokenizing strings
1403 
1404   Level: intermediate
1405 
1406 .seealso: PetscTokenCreate(), PetscTokenFind(), PetscTokenDestroy()
1407 S*/
1408 typedef struct _p_PetscToken* PetscToken;
1409 
1410 PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1411 PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1412 PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);
1413 
1414 PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1415 PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);
1416 
1417 /*
1418    These are  MPI operations for MPI_Allreduce() etc
1419 */
1420 PETSC_EXTERN MPI_Op PetscMaxSum_Op;
1421 #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128)
1422 PETSC_EXTERN MPI_Op MPIU_SUM;
1423 #else
1424 #define MPIU_SUM MPI_SUM
1425 #endif
1426 #if defined(PETSC_USE_REAL___FLOAT128)
1427 PETSC_EXTERN MPI_Op MPIU_MAX;
1428 PETSC_EXTERN MPI_Op MPIU_MIN;
1429 #else
1430 #define MPIU_MAX MPI_MAX
1431 #define MPIU_MIN MPI_MIN
1432 #endif
1433 PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);
1434 
1435 PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1436 PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1437 
1438 /*S
1439      PetscObject - any PETSc object, PetscViewer, Mat, Vec, KSP etc
1440 
1441    Level: beginner
1442 
1443    Note: This is the base class from which all PETSc objects are derived from.
1444 
1445 .seealso:  PetscObjectDestroy(), PetscObjectView(), PetscObjectGetName(), PetscObjectSetName(), PetscObjectReference(), PetscObjectDereferenc()
1446 S*/
1447 typedef struct _p_PetscObject* PetscObject;
1448 
1449 /*MC
1450     PetscObjectId - unique integer Id for a PetscObject
1451 
1452     Level: developer
1453 
1454     Notes: Unlike pointer values, object ids are never reused.
1455 
1456 .seealso: PetscObjectState, PetscObjectGetId()
1457 M*/
1458 typedef Petsc64bitInt PetscObjectId;
1459 
1460 /*MC
1461     PetscObjectState - integer state for a PetscObject
1462 
1463     Level: developer
1464 
1465     Notes:
1466     Object state is always-increasing and (for objects that track state) can be used to determine if an object has
1467     changed since the last time you interacted with it.  It is 64-bit so that it will not overflow for a very long time.
1468 
1469 .seealso: PetscObjectId, PetscObjectStateGet(), PetscObjectStateIncrease(), PetscObjectStateSet()
1470 M*/
1471 typedef Petsc64bitInt PetscObjectState;
1472 
1473 /*S
1474      PetscFunctionList - Linked list of functions, possibly stored in dynamic libraries, accessed
1475       by string name
1476 
1477    Level: advanced
1478 
1479 .seealso:  PetscFunctionListAdd(), PetscFunctionListDestroy(), PetscOpFlist
1480 S*/
1481 typedef struct _n_PetscFunctionList *PetscFunctionList;
1482 
1483 /*E
1484   PetscFileMode - Access mode for a file.
1485 
1486   Level: beginner
1487 
1488   FILE_MODE_READ - open a file at its beginning for reading
1489 
1490   FILE_MODE_WRITE - open a file at its beginning for writing (will create if the file does not exist)
1491 
1492   FILE_MODE_APPEND - open a file at end for writing
1493 
1494   FILE_MODE_UPDATE - open a file for updating, meaning for reading and writing
1495 
1496   FILE_MODE_APPEND_UPDATE - open a file for updating, meaning for reading and writing, at the end
1497 
1498 .seealso: PetscViewerFileSetMode()
1499 E*/
1500 typedef enum {FILE_MODE_READ, FILE_MODE_WRITE, FILE_MODE_APPEND, FILE_MODE_UPDATE, FILE_MODE_APPEND_UPDATE} PetscFileMode;
1501 extern const char *const PetscFileModes[];
1502 
1503 /*
1504     Defines PETSc error handling.
1505 */
1506 #include <petscerror.h>
1507 
1508 #define PETSC_SMALLEST_CLASSID  1211211
1509 PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1510 PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1511 PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);
1512 
1513 /*
1514    Routines that get memory usage information from the OS
1515 */
1516 PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1517 PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1518 PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1519 
1520 PETSC_EXTERN PetscErrorCode PetscInfoAllow(PetscBool ,const char []);
1521 PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);
1522 
1523 /*
1524    Initialization of PETSc
1525 */
1526 PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1527 PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1528 PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1529 PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1530 PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1531 PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1532 PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1533 PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1534 PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1535 PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);
1536 
1537 PETSC_EXTERN PetscErrorCode PetscEnd(void);
1538 PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);
1539 
1540 PETSC_EXTERN MPI_Comm PETSC_COMM_LOCAL_WORLD;
1541 PETSC_EXTERN PetscErrorCode PetscHMPIMerge(PetscMPIInt,PetscErrorCode (*)(void*),void*);
1542 PETSC_EXTERN PetscErrorCode PetscHMPISpawn(PetscMPIInt);
1543 PETSC_EXTERN PetscErrorCode PetscHMPIFinalize(void);
1544 PETSC_EXTERN PetscErrorCode PetscHMPIRun(MPI_Comm,PetscErrorCode (*)(MPI_Comm,void *),void*);
1545 PETSC_EXTERN PetscErrorCode PetscHMPIRunCtx(MPI_Comm,PetscErrorCode (*)(MPI_Comm,void*,void *),void*);
1546 PETSC_EXTERN PetscErrorCode PetscHMPIFree(MPI_Comm,void*);
1547 PETSC_EXTERN PetscErrorCode PetscHMPIMalloc(MPI_Comm,size_t,void**);
1548 
1549 PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1550 PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1551 PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1552 PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);
1553 
1554 /*
1555      These are so that in extern C code we can caste function pointers to non-extern C
1556    function pointers. Since the regular C++ code expects its function pointers to be C++
1557 */
1558 PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1559 PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1560 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);
1561 
1562 /*
1563     Functions that can act on any PETSc object.
1564 */
1565 PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1566 PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1567 PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1568 PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1569 PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1570 PETSC_EXTERN PetscErrorCode PetscObjectSetPrecision(PetscObject,PetscPrecision);
1571 PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1572 PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1573 PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1574 PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1575 PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1576 PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1577 PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1578 PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1579 PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1580 PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt *);
1581 PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1582 PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1583 PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject *);
1584 PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void));
1585 #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d))
1586 PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1587 PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1588 PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt *);
1589 PETSC_EXTERN PetscErrorCode PetscObjectAddOptionsHandler(PetscObject,PetscErrorCode (*)(PetscObject,void*),PetscErrorCode (*)(PetscObject,void*),void*);
1590 PETSC_EXTERN PetscErrorCode PetscObjectProcessOptionsHandlers(PetscObject);
1591 PETSC_EXTERN PetscErrorCode PetscObjectDestroyOptionsHandlers(PetscObject);
1592 PETSC_EXTERN PetscErrorCode PetscObjectsGetGlobalNumbering(MPI_Comm,PetscInt,PetscObject*,PetscInt*,PetscInt*);
1593 
1594 #include <petscviewertypes.h>
1595 #include <petscoptions.h>
1596 
1597 PETSC_EXTERN PetscErrorCode PetscMemoryShowUsage(PetscViewer,const char[]);
1598 PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1599 PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1600 #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1601 PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1602 PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1603 PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1604 PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1605 PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1606 PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1607 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1608 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1609 PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1610 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1611 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1612 PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1613 PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);
1614 
1615 #if defined(PETSC_HAVE_SAWS)
1616 PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1617 PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1618 PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1619 PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1620 PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1621 PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1622 PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1623 PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1624 PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1625 PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);
1626 
1627 #else
1628 #define PetscSAWsBlock()                        0
1629 #define PetscObjectSAWsViewOff(obj)             0
1630 #define PetscObjectSAWsSetBlock(obj,flg)        0
1631 #define PetscObjectSAWsBlock(obj)               0
1632 #define PetscObjectSAWsGrantAccess(obj)         0
1633 #define PetscObjectSAWsTakeAccess(obj)          0
1634 #define PetscStackViewSAWs()                    0
1635 #define PetscStackSAWsViewOff()                 0
1636 #define PetscStackSAWsTakeAccess()
1637 #define PetscStackSAWsGrantAccess()
1638 
1639 #endif
1640 
1641 typedef void* PetscDLHandle;
1642 typedef enum {PETSC_DL_DECIDE=0,PETSC_DL_NOW=1,PETSC_DL_LOCAL=2} PetscDLMode;
1643 PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1644 PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1645 PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);
1646 
1647 
1648 #if defined(PETSC_USE_DEBUG)
1649 PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1650 #endif
1651 PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);
1652 
1653 /*S
1654      PetscObjectList - Linked list of PETSc objects, each accessable by string name
1655 
1656    Level: developer
1657 
1658    Notes: Used by PetscObjectCompose() and PetscObjectQuery()
1659 
1660 .seealso:  PetscObjectListAdd(), PetscObjectListDestroy(), PetscObjectListFind(), PetscObjectCompose(), PetscObjectQuery(), PetscFunctionList
1661 S*/
1662 typedef struct _n_PetscObjectList *PetscObjectList;
1663 
1664 PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1665 PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1666 PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1667 PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1668 PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1669 PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);
1670 
1671 /*
1672     Dynamic library lists. Lists of names of routines in objects or in dynamic
1673   link libraries that will be loaded as needed.
1674 */
1675 
1676 #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1677 PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1678 PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1679 #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1680 PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1681 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[]);
1682 PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1683 PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1684 PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);
1685 
1686 /*S
1687      PetscDLLibrary - Linked list of dynamics libraries to search for functions
1688 
1689    Level: advanced
1690 
1691 .seealso:  PetscDLLibraryOpen()
1692 S*/
1693 typedef struct _n_PetscDLLibrary *PetscDLLibrary;
1694 PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1695 PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1696 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1697 PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1698 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1699 PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1700 PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1701 PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);
1702 
1703 /*
1704      Useful utility routines
1705 */
1706 PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1707 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1708 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1709 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1710 PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1711 PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);
1712 
1713 /*
1714     PetscNot - negates a logical type value and returns result as a PetscBool
1715 
1716     Notes: This is useful in cases like
1717 $     int        *a;
1718 $     PetscBool  flag = PetscNot(a)
1719      where !a does not return a PetscBool  because we cannot provide a cast from int to PetscBool  in C.
1720 */
1721 #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)
1722 
1723 #if defined(PETSC_HAVE_VALGRIND)
1724 #  include <valgrind/valgrind.h>
1725 #  define PETSC_RUNNING_ON_VALGRIND RUNNING_ON_VALGRIND
1726 #else
1727 #  define PETSC_RUNNING_ON_VALGRIND PETSC_FALSE
1728 #endif
1729 
1730 
1731 /*MC
1732     PetscHelpPrintf - Prints help messages.
1733 
1734    Synopsis:
1735     #include "petscsys.h"
1736      PetscErrorCode (*PetscHelpPrintf)(const char format[],...);
1737 
1738     Not Collective
1739 
1740     Input Parameters:
1741 .   format - the usual printf() format string
1742 
1743    Level: developer
1744 
1745     Fortran Note:
1746     This routine is not supported in Fortran.
1747 
1748     Concepts: help messages^printing
1749     Concepts: printing^help messages
1750 
1751 .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1752 M*/
1753 PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);
1754 
1755 /*
1756      Defines PETSc profiling.
1757 */
1758 #include <petsclog.h>
1759 
1760 
1761 
1762 /*
1763       Simple PETSc parallel IO for ASCII printing
1764 */
1765 PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1766 PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1767 PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1768 PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1769 PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1770 PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1771 PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);
1772 
1773 /* These are used internally by PETSc ASCII IO routines*/
1774 #include <stdarg.h>
1775 PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list);
1776 PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE*,const char[],va_list);
1777 PETSC_EXTERN PetscErrorCode PetscVFPrintfDefault(FILE*,const char[],va_list);
1778 
1779 #if defined(PETSC_HAVE_MATLAB_ENGINE)
1780 PETSC_EXTERN PetscErrorCode PetscVFPrintf_Matlab(FILE*,const char[],va_list);
1781 #endif
1782 
1783 PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1784 PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1785 PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);
1786 
1787 #if defined(PETSC_HAVE_POPEN)
1788 PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1789 PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*,int*);
1790 PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1791 #endif
1792 
1793 PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1794 PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1795 PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*);
1796 PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1797 PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1798 PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1799 PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);
1800 
1801 PETSC_EXTERN PetscErrorCode PetscPopUpSelect(MPI_Comm,const char*,const char*,int,const char**,int*);
1802 
1803 /*S
1804      PetscContainer - Simple PETSc object that contains a pointer to any required data
1805 
1806    Level: advanced
1807 
1808 .seealso:  PetscObject, PetscContainerCreate()
1809 S*/
1810 PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1811 typedef struct _p_PetscContainer*  PetscContainer;
1812 PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void **);
1813 PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void *);
1814 PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1815 PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer *);
1816 PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void*));
1817 
1818 /*
1819    For use in debuggers
1820 */
1821 PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1822 PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1823 PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1824 PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1825 PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);
1826 
1827 #include <stddef.h>
1828 #if defined(PETSC_HAVE_MEMORY_H)
1829 #include <memory.h>
1830 #endif
1831 #if defined(PETSC_HAVE_STDLIB_H)
1832 #include <stdlib.h>
1833 #endif
1834 
1835 #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1836 #include <xmmintrin.h>
1837 #endif
1838 
1839 #undef __FUNCT__
1840 #define __FUNCT__ "PetscMemcpy"
1841 /*@C
1842    PetscMemcpy - Copies n bytes, beginning at location b, to the space
1843    beginning at location a. The two memory regions CANNOT overlap, use
1844    PetscMemmove() in that case.
1845 
1846    Not Collective
1847 
1848    Input Parameters:
1849 +  b - pointer to initial memory space
1850 -  n - length (in bytes) of space to copy
1851 
1852    Output Parameter:
1853 .  a - pointer to copy space
1854 
1855    Level: intermediate
1856 
1857    Compile Option:
1858     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1859                                   for memory copies on double precision values.
1860     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1861                                   for memory copies on double precision values.
1862     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1863                                   for memory copies on double precision values.
1864 
1865    Note:
1866    This routine is analogous to memcpy().
1867 
1868    Developer Note: this is inlined for fastest performance
1869 
1870   Concepts: memory^copying
1871   Concepts: copying^memory
1872 
1873 .seealso: PetscMemmove()
1874 
1875 @*/
1876 PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n)
1877 {
1878 #if defined(PETSC_USE_DEBUG)
1879   unsigned long al = (unsigned long) a,bl = (unsigned long) b;
1880   unsigned long nl = (unsigned long) n;
1881   PetscFunctionBegin;
1882   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1883   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1884 #else
1885   PetscFunctionBegin;
1886 #endif
1887   if (a != b) {
1888 #if defined(PETSC_USE_DEBUG)
1889     if ((al > bl && (al - bl) < nl) || (bl - al) < nl)  SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1890               or make sure your copy regions and lengths are correct. \n\
1891               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1892 #endif
1893 #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1894    if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1895       size_t len = n/sizeof(PetscScalar);
1896 #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1897       PetscBLASInt   one = 1,blen;
1898       PetscErrorCode ierr;
1899       ierr = PetscBLASIntCast(len,&blen);CHKERRQ(ierr);
1900       PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one));
1901 #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1902       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1903 #else
1904       size_t      i;
1905       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1906       for (i=0; i<len; i++) y[i] = x[i];
1907 #endif
1908     } else {
1909       memcpy((char*)(a),(char*)(b),n);
1910     }
1911 #else
1912     memcpy((char*)(a),(char*)(b),n);
1913 #endif
1914   }
1915   PetscFunctionReturn(0);
1916 }
1917 
1918 /*@C
1919    PetscMemzero - Zeros the specified memory.
1920 
1921    Not Collective
1922 
1923    Input Parameters:
1924 +  a - pointer to beginning memory location
1925 -  n - length (in bytes) of memory to initialize
1926 
1927    Level: intermediate
1928 
1929    Compile Option:
1930    PETSC_PREFER_BZERO - on certain machines (the IBM RS6000) the bzero() routine happens
1931   to be faster than the memset() routine. This flag causes the bzero() routine to be used.
1932 
1933    Developer Note: this is inlined for fastest performance
1934 
1935    Concepts: memory^zeroing
1936    Concepts: zeroing^memory
1937 
1938 .seealso: PetscMemcpy()
1939 @*/
1940 PETSC_STATIC_INLINE PetscErrorCode  PetscMemzero(void *a,size_t n)
1941 {
1942   if (n > 0) {
1943 #if defined(PETSC_USE_DEBUG)
1944     if (!a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer");
1945 #endif
1946 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1947     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1948       size_t      i,len = n/sizeof(PetscScalar);
1949       PetscScalar *x = (PetscScalar*)a;
1950       for (i=0; i<len; i++) x[i] = 0.0;
1951     } else {
1952 #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1953     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1954       PetscInt len = n/sizeof(PetscScalar);
1955       fortranzero_(&len,(PetscScalar*)a);
1956     } else {
1957 #endif
1958 #if defined(PETSC_PREFER_BZERO)
1959       bzero((char *)a,n);
1960 #else
1961       memset((char*)a,0,n);
1962 #endif
1963 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1964     }
1965 #endif
1966   }
1967   return 0;
1968 }
1969 
1970 /*MC
1971    PetscPrefetchBlock - Prefetches a block of memory
1972 
1973    Synopsis:
1974     #include "petscsys.h"
1975     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)
1976 
1977    Not Collective
1978 
1979    Input Parameters:
1980 +  a - pointer to first element to fetch (any type but usually PetscInt or PetscScalar)
1981 .  n - number of elements to fetch
1982 .  rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors)
1983 -  t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note
1984 
1985    Level: developer
1986 
1987    Notes:
1988    The last two arguments (rw and t) must be compile-time constants.
1989 
1990    Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality.  Not all architectures offer
1991    equivalent locality hints, but the following macros are always defined to their closest analogue.
1992 +  PETSC_PREFETCH_HINT_NTA - Non-temporal.  Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched).
1993 .  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.
1994 .  PETSC_PREFETCH_HINT_T1 - Fetch to level 2 and higher (not L1).
1995 -  PETSC_PREFETCH_HINT_T2 - Fetch to high-level cache only.  (On many systems, T0 and T1 are equivalent.)
1996 
1997    This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid
1998    address).
1999 
2000    Concepts: memory
2001 M*/
2002 #define PetscPrefetchBlock(a,n,rw,t) do {                               \
2003     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
2004     for ( ; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
2005   } while (0)
2006 
2007 /*
2008       Determine if some of the kernel computation routines use
2009    Fortran (rather than C) for the numerical calculations. On some machines
2010    and compilers (like complex numbers) the Fortran version of the routines
2011    is faster than the C/C++ versions. The flag --with-fortran-kernels
2012    should be used with ./configure to turn these on.
2013 */
2014 #if defined(PETSC_USE_FORTRAN_KERNELS)
2015 
2016 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
2017 #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
2018 #endif
2019 
2020 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
2021 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
2022 #endif
2023 
2024 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
2025 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
2026 #endif
2027 
2028 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
2029 #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
2030 #endif
2031 
2032 #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
2033 #define PETSC_USE_FORTRAN_KERNEL_NORM
2034 #endif
2035 
2036 #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
2037 #define PETSC_USE_FORTRAN_KERNEL_MAXPY
2038 #endif
2039 
2040 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
2041 #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
2042 #endif
2043 
2044 #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
2045 #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
2046 #endif
2047 
2048 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
2049 #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
2050 #endif
2051 
2052 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
2053 #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
2054 #endif
2055 
2056 #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
2057 #define PETSC_USE_FORTRAN_KERNEL_MDOT
2058 #endif
2059 
2060 #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
2061 #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
2062 #endif
2063 
2064 #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
2065 #define PETSC_USE_FORTRAN_KERNEL_AYPX
2066 #endif
2067 
2068 #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
2069 #define PETSC_USE_FORTRAN_KERNEL_WAXPY
2070 #endif
2071 
2072 #endif
2073 
2074 /*
2075     Macros for indicating code that should be compiled with a C interface,
2076    rather than a C++ interface. Any routines that are dynamically loaded
2077    (such as the PCCreate_XXX() routines) must be wrapped so that the name
2078    mangler does not change the functions symbol name. This just hides the
2079    ugly extern "C" {} wrappers.
2080 */
2081 #if defined(__cplusplus)
2082 #define EXTERN_C_BEGIN extern "C" {
2083 #define EXTERN_C_END }
2084 #else
2085 #define EXTERN_C_BEGIN
2086 #define EXTERN_C_END
2087 #endif
2088 
2089 /* --------------------------------------------------------------------*/
2090 
2091 /*MC
2092     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
2093         communication
2094 
2095    Level: beginner
2096 
2097    Note: This manual page is a place-holder because MPICH does not have a manual page for MPI_Comm
2098 
2099 .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
2100 M*/
2101 
2102 /*MC
2103     PetscScalar - PETSc type that represents either a double precision real number, a double precision
2104        complex number, a single precision real number, a long double or an int - if the code is configured
2105        with --with-scalar-type=real,complex --with-precision=single,double,longdouble,int,matsingle
2106 
2107 
2108    Level: beginner
2109 
2110 .seealso: PetscReal, PassiveReal, PassiveScalar, MPIU_SCALAR, PetscInt
2111 M*/
2112 
2113 /*MC
2114     PetscComplex - PETSc type that represents a complex number with precision matching that of PetscReal.
2115 
2116    Synopsis:
2117    #define PETSC_DESIRE_COMPLEX
2118    #include <petscsys.h>
2119    PetscComplex number = 1. + 2.*PETSC_i;
2120 
2121    Level: beginner
2122 
2123    Note:
2124    Complex numbers are automatically available if PETSc was configured --with-scalar-type=complex (in which case
2125    PetscComplex will match PetscScalar), otherwise the macro PETSC_DESIRE_COMPLEX must be defined before including any
2126    PETSc headers. If the compiler supports complex numbers, PetscComplex and associated variables and functions will be
2127    defined and PETSC_HAVE_COMPLEX will be set.
2128 
2129 .seealso: PetscReal, PetscComplex, MPIU_COMPLEX, PetscInt, PETSC_i
2130 M*/
2131 
2132 /*MC
2133     PetscReal - PETSc type that represents a real number version of PetscScalar
2134 
2135    Level: beginner
2136 
2137 .seealso: PetscScalar, PassiveReal, PassiveScalar
2138 M*/
2139 
2140 /*MC
2141     PassiveScalar - PETSc type that represents a PetscScalar
2142    Level: beginner
2143 
2144     This is the same as a PetscScalar except in code that is automatically differentiated it is
2145    treated as a constant (not an indendent or dependent variable)
2146 
2147 .seealso: PetscReal, PassiveReal, PetscScalar
2148 M*/
2149 
2150 /*MC
2151     PassiveReal - PETSc type that represents a PetscReal
2152 
2153    Level: beginner
2154 
2155     This is the same as a PetscReal except in code that is automatically differentiated it is
2156    treated as a constant (not an indendent or dependent variable)
2157 
2158 .seealso: PetscScalar, PetscReal, PassiveScalar
2159 M*/
2160 
2161 /*MC
2162     MPIU_SCALAR - MPI datatype corresponding to PetscScalar
2163 
2164    Level: beginner
2165 
2166     Note: In MPI calls that require an MPI datatype that matches a PetscScalar or array of PetscScalars
2167           pass this value
2168 
2169 .seealso: PetscReal, PassiveReal, PassiveScalar, PetscScalar, MPIU_INT
2170 M*/
2171 
2172 #if defined(PETSC_HAVE_MPIIO)
2173 #if !defined(PETSC_WORDS_BIGENDIAN)
2174 PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2175 PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2176 #else
2177 #define MPIU_File_write_all(a,b,c,d,e) MPI_File_write_all(a,b,c,d,e)
2178 #define MPIU_File_read_all(a,b,c,d,e) MPI_File_read_all(a,b,c,d,e)
2179 #endif
2180 #endif
2181 
2182 /* the following petsc_static_inline require petscerror.h */
2183 
2184 /* Limit MPI to 32-bits */
2185 #define PETSC_MPI_INT_MAX  2147483647
2186 #define PETSC_MPI_INT_MIN -2147483647
2187 /* Limit BLAS to 32-bits */
2188 #define PETSC_BLAS_INT_MAX  2147483647
2189 #define PETSC_BLAS_INT_MIN -2147483647
2190 
2191 #undef __FUNCT__
2192 #define __FUNCT__ "PetscBLASIntCast"
2193 PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
2194 {
2195   PetscFunctionBegin;
2196 #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
2197   if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK");
2198 #endif
2199   *b =  (PetscBLASInt)(a);
2200   PetscFunctionReturn(0);
2201 }
2202 
2203 #undef __FUNCT__
2204 #define __FUNCT__ "PetscMPIIntCast"
2205 PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
2206 {
2207   PetscFunctionBegin;
2208 #if defined(PETSC_USE_64BIT_INDICES)
2209   if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for MPI");
2210 #endif
2211   *b =  (PetscMPIInt)(a);
2212   PetscFunctionReturn(0);
2213 }
2214 
2215 
2216 /*
2217      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2218 */
2219 #if defined(hz)
2220 #undef hz
2221 #endif
2222 
2223 /*  For arrays that contain filenames or paths */
2224 
2225 
2226 #if defined(PETSC_HAVE_LIMITS_H)
2227 #include <limits.h>
2228 #endif
2229 #if defined(PETSC_HAVE_SYS_PARAM_H)
2230 #include <sys/param.h>
2231 #endif
2232 #if defined(PETSC_HAVE_SYS_TYPES_H)
2233 #include <sys/types.h>
2234 #endif
2235 #if defined(MAXPATHLEN)
2236 #  define PETSC_MAX_PATH_LEN     MAXPATHLEN
2237 #elif defined(MAX_PATH)
2238 #  define PETSC_MAX_PATH_LEN     MAX_PATH
2239 #elif defined(_MAX_PATH)
2240 #  define PETSC_MAX_PATH_LEN     _MAX_PATH
2241 #else
2242 #  define PETSC_MAX_PATH_LEN     4096
2243 #endif
2244 
2245 /* Special support for C++ */
2246 #if defined(PETSC_CLANGUAGE_CXX) && defined(__cplusplus)
2247 #include <petscsys.hh>
2248 #endif
2249 
2250 /*MC
2251 
2252     UsingFortran - Fortran can be used with PETSc in four distinct approaches
2253 
2254 $    1) classic Fortran 77 style
2255 $#include "finclude/petscXXX.h" to work with material from the XXX component of PETSc
2256 $       XXX variablename
2257 $      You cannot use this approach if you wish to use the Fortran 90 specific PETSc routines
2258 $      which end in F90; such as VecGetArrayF90()
2259 $
2260 $    2) classic Fortran 90 style
2261 $#include "finclude/petscXXX.h"
2262 $#include "finclude/petscXXX.h90" to work with material from the XXX component of PETSc
2263 $       XXX variablename
2264 $
2265 $    3) Using Fortran modules
2266 $#include "finclude/petscXXXdef.h"
2267 $         use petscXXXX
2268 $       XXX variablename
2269 $
2270 $    4) Use Fortran modules and Fortran data types for PETSc types
2271 $#include "finclude/petscXXXdef.h"
2272 $         use petscXXXX
2273 $       type(XXX) variablename
2274 $      To use this approach you must ./configure PETSc with the additional
2275 $      option --with-fortran-datatypes You cannot use the type(XXX) declaration approach without using Fortran modules
2276 
2277     Finally if you absolutely do not want to use any #include you can use either
2278 
2279 $    3a) skip the #include BUT you cannot use any PETSc data type names like Vec, Mat, PetscInt, PetscErrorCode etc
2280 $        and you must declare the variables as integer, for example
2281 $        integer variablename
2282 $
2283 $    4a) skip the #include, you use the object types like type(Vec) type(Mat) but cannot use the data type
2284 $        names like PetscErrorCode, PetscInt etc. again for those you must use integer
2285 
2286    We recommend either 2 or 3. Approaches 2 and 3 provide type checking for most PETSc function calls; 4 has type checking
2287 for only a few PETSc functions.
2288 
2289    Fortran type checking with interfaces is strick, this means you cannot pass a scalar value when an array value
2290 is expected (even though it is legal Fortran). For example when setting a single value in a matrix with MatSetValues()
2291 you cannot have something like
2292 $      PetscInt row,col
2293 $      PetscScalar val
2294 $        ...
2295 $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)
2296 You must instead have
2297 $      PetscInt row(1),col(1)
2298 $      PetscScalar val(1)
2299 $        ...
2300 $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)
2301 
2302 
2303     See the example src/vec/vec/examples/tutorials/ex20f90.F90 for an example that can use all four approaches
2304 
2305     Developer Notes: The finclude/petscXXXdef.h contain all the #defines (would be typedefs in C code) these
2306      automatically include their predecessors; for example finclude/petscvecdef.h includes finclude/petscisdef.h
2307 
2308      The finclude/petscXXXX.h contain all the parameter statements for that package. These automatically include
2309      their finclude/petscXXXdef.h file but DO NOT automatically include their predecessors;  for example
2310      finclude/petscvec.h does NOT automatically include finclude/petscis.h
2311 
2312      The finclude/ftn-custom/petscXXXdef.h90 are not intended to be used directly in code, they define the
2313      Fortran data type type(XXX) (for example type(Vec)) when PETSc is ./configure with the --with-fortran-datatypes option.
2314 
2315      The finclude/ftn-custom/petscXXX.h90 (not included directly by code) contain interface definitions for
2316      the PETSc Fortran stubs that have different bindings then their C version (for example VecGetArrayF90).
2317 
2318      The finclude/ftn-auto/petscXXX.h90 (not included directly by code) contain interface definitions generated
2319      automatically by "make allfortranstubs".
2320 
2321      The finclude/petscXXX.h90 includes the custom finclude/ftn-custom/petscXXX.h90 and if ./configure
2322      was run with --with-fortran-interfaces it also includes the finclude/ftn-auto/petscXXX.h90 These DO NOT automatically
2323      include their predecessors
2324 
2325     Level: beginner
2326 
2327 M*/
2328 
2329 PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2330 PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2331 PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2332 PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2333 PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2334 PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2335 PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2336 
2337 PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2338 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2339 PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2340 PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2341 PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2342 PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2343 PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt*,PetscInt*,PetscInt*);
2344 PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2345 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2346 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2347 PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2348 PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2349 PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2350 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2351 PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2352 PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2353 PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2354 PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt*,const PetscInt*,PetscInt,const PetscInt*,const PetscInt*,PetscInt*,PetscInt**,PetscInt**);
2355 
2356 PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2357 PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);
2358 
2359 /*J
2360     PetscRandomType - String with the name of a PETSc randomizer
2361 
2362    Level: beginner
2363 
2364    Notes: to use the SPRNG you must have ./configure PETSc
2365    with the option --download-sprng
2366 
2367 .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2368 J*/
2369 typedef const char* PetscRandomType;
2370 #define PETSCRAND       "rand"
2371 #define PETSCRAND48     "rand48"
2372 #define PETSCSPRNG      "sprng"
2373 
2374 /* Logging support */
2375 PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;
2376 
2377 PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);
2378 
2379 /*S
2380      PetscRandom - Abstract PETSc object that manages generating random numbers
2381 
2382    Level: intermediate
2383 
2384   Concepts: random numbers
2385 
2386 .seealso:  PetscRandomCreate(), PetscRandomGetValue(), PetscRandomType
2387 S*/
2388 typedef struct _p_PetscRandom*   PetscRandom;
2389 
2390 /* Dynamic creation and loading functions */
2391 PETSC_EXTERN PetscFunctionList PetscRandomList;
2392 PETSC_EXTERN PetscBool         PetscRandomRegisterAllCalled;
2393 
2394 PETSC_EXTERN PetscErrorCode PetscRandomRegisterAll(void);
2395 PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2396 PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2397 PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2398 PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType*);
2399  PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom,const char[],const char[]);
2400 PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);
2401 
2402 PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2403 PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2404 PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2405 PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2406 PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2407 PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2408 PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2409 PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2410 PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);
2411 
2412 PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2413 PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2414 PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2415 PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2416 PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2417 PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2418 PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);
2419 
2420 PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscDataType);
2421 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscDataType);
2422 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,void*,PetscInt,PetscDataType,PetscBool );
2423 PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,void*,PetscInt,PetscDataType,PetscBool );
2424 PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2425 PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2426 PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2427 PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2428 PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2429 PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2430 PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2431 PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);
2432 PETSC_EXTERN PetscErrorCode PetscWebServe(MPI_Comm,int);
2433 
2434 /*
2435    In binary files variables are stored using the following lengths,
2436   regardless of how they are stored in memory on any one particular
2437   machine. Use these rather then sizeof() in computing sizes for
2438   PetscBinarySeek().
2439 */
2440 #define PETSC_BINARY_INT_SIZE   (32/8)
2441 #define PETSC_BINARY_FLOAT_SIZE  (32/8)
2442 #define PETSC_BINARY_CHAR_SIZE  (8/8)
2443 #define PETSC_BINARY_SHORT_SIZE  (16/8)
2444 #define PETSC_BINARY_DOUBLE_SIZE  (64/8)
2445 #define PETSC_BINARY_SCALAR_SIZE  sizeof(PetscScalar)
2446 
2447 /*E
2448   PetscBinarySeekType - argument to PetscBinarySeek()
2449 
2450   Level: advanced
2451 
2452 .seealso: PetscBinarySeek(), PetscBinarySynchronizedSeek()
2453 E*/
2454 typedef enum {PETSC_BINARY_SEEK_SET = 0,PETSC_BINARY_SEEK_CUR = 1,PETSC_BINARY_SEEK_END = 2} PetscBinarySeekType;
2455 PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2456 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2457 PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);
2458 
2459 PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2460 PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool );
2461 PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2462 PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2463 PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2464 PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2465 
2466 PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2467 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2468 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2469 PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2470 PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2471 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscInt,const PetscMPIInt*,const void*,PetscInt*,PetscMPIInt**,void*) PetscAttrMPIPointerWithType(6,3);
2472 
2473 /*E
2474     PetscBuildTwoSidedType - algorithm for setting up two-sided communication
2475 
2476 $  PETSC_BUILDTWOSIDED_ALLREDUCE - classical algorithm using an MPI_Allreduce with
2477 $      a buffer of length equal to the communicator size. Not memory-scalable due to
2478 $      the large reduction size. Requires only MPI-1.
2479 $  PETSC_BUILDTWOSIDED_IBARRIER - nonblocking algorithm based on MPI_Issend and MPI_Ibarrier.
2480 $      Proved communication-optimal in Hoefler, Siebert, and Lumsdaine (2010). Requires MPI-3.
2481 
2482    Level: developer
2483 
2484 .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedSetType(), PetscCommBuildTwoSidedGetType()
2485 E*/
2486 typedef enum {
2487   PETSC_BUILDTWOSIDED_NOTSET = -1,
2488   PETSC_BUILDTWOSIDED_ALLREDUCE = 0,
2489   PETSC_BUILDTWOSIDED_IBARRIER = 1
2490   /* Updates here must be accompanied by updates in finclude/petscsys.h and the string array in mpits.c */
2491 } PetscBuildTwoSidedType;
2492 PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2493 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2494 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);
2495 
2496 PETSC_EXTERN PetscErrorCode PetscSSEIsEnabled(MPI_Comm,PetscBool  *,PetscBool  *);
2497 
2498 /*E
2499   InsertMode - Whether entries are inserted or added into vectors or matrices
2500 
2501   Level: beginner
2502 
2503 .seealso: VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2504           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(),
2505           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd()
2506 E*/
2507  typedef enum {NOT_SET_VALUES, INSERT_VALUES, ADD_VALUES, MAX_VALUES, INSERT_ALL_VALUES, ADD_ALL_VALUES, INSERT_BC_VALUES, ADD_BC_VALUES} InsertMode;
2508 
2509 /*MC
2510     INSERT_VALUES - Put a value into a vector or matrix, overwrites any previous value
2511 
2512     Level: beginner
2513 
2514 .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2515           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), ADD_VALUES,
2516           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES
2517 
2518 M*/
2519 
2520 /*MC
2521     ADD_VALUES - Adds a value into a vector or matrix, if there previously was no value, just puts the
2522                 value into that location
2523 
2524     Level: beginner
2525 
2526 .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2527           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), INSERT_VALUES,
2528           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES
2529 
2530 M*/
2531 
2532 /*MC
2533     MAX_VALUES - Puts the maximum of the scattered/gathered value and the current value into each location
2534 
2535     Level: beginner
2536 
2537 .seealso: InsertMode, VecScatterBegin(), VecScatterEnd(), ADD_VALUES, INSERT_VALUES
2538 
2539 M*/
2540 
2541 PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);
2542 
2543 typedef enum {PETSC_SUBCOMM_GENERAL=0,PETSC_SUBCOMM_CONTIGUOUS=1,PETSC_SUBCOMM_INTERLACED=2} PetscSubcommType;
2544 PETSC_EXTERN const char *const PetscSubcommTypes[];
2545 
2546 /*S
2547    PetscSubcomm - Context of MPI subcommunicators, used by PCREDUNDANT
2548 
2549    Level: advanced
2550 
2551    Concepts: communicator, create
2552 S*/
2553 typedef struct _n_PetscSubcomm* PetscSubcomm;
2554 
2555 struct _n_PetscSubcomm {
2556   MPI_Comm    parent;           /* parent communicator */
2557   MPI_Comm    dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2558   MPI_Comm    comm;             /* this communicator */
2559   PetscMPIInt n;                /* num of subcommunicators under the parent communicator */
2560   PetscMPIInt color;            /* color of processors belong to this communicator */
2561   PetscMPIInt *subsize;         /* size of subcommunicator[color] */
2562   PetscSubcommType type;
2563 };
2564 
2565 PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2566 PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2567 PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2568 PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2569 PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2570 PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2571 PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2572 
2573 /*S
2574    PetscSegBuffer - a segmented extendable buffer
2575 
2576    Level: developer
2577 
2578 .seealso: PetscSegBufferCreate(), PetscSegBufferGet(), PetscSegBufferExtract(), PetscSegBufferDestroy()
2579 S*/
2580 typedef struct _n_PetscSegBuffer *PetscSegBuffer;
2581 PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2582 PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2583 PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2584 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2585 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2586 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2587 PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2588 PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);
2589 
2590 /* Type-safe wrapper to encourage use of PETSC_RESTRICT. Does not use PetscFunctionBegin because the error handling
2591  * prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as
2592  * possible. */
2593 PETSC_STATIC_INLINE PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg,PetscInt count,PetscInt *PETSC_RESTRICT *slot) {return PetscSegBufferGet(seg,count,(void**)slot);}
2594 
2595 extern PetscSegBuffer PetscCitationsList;
2596 #undef __FUNCT__
2597 #define __FUNCT__ "PetscCitationsRegister"
2598 /*@C
2599       PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code.
2600 
2601      Not Collective - only what is registered on rank 0 of PETSC_COMM_WORLD will be printed
2602 
2603      Input Parameters:
2604 +      cite - the bibtex item, formated to displayed on multiple lines nicely
2605 -      set - a boolean variable initially set to PETSC_FALSE; this is used to insure only a single registration of the citation
2606 
2607    Level: intermediate
2608 
2609      Options Database:
2610 .     -citations [filenmae]   - print out the bibtex entries for the given computation
2611 @*/
2612 PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2613 {
2614   size_t         len;
2615   char           *vstring;
2616   PetscErrorCode ierr;
2617 
2618   PetscFunctionBegin;
2619   if (set && *set) PetscFunctionReturn(0);
2620   ierr = PetscStrlen(cit,&len);CHKERRQ(ierr);
2621   ierr = PetscSegBufferGet(PetscCitationsList,(PetscInt)len,&vstring);CHKERRQ(ierr);
2622   ierr = PetscMemcpy(vstring,cit,len);CHKERRQ(ierr);
2623   if (set) *set = PETSC_TRUE;
2624   PetscFunctionReturn(0);
2625 }
2626 
2627 
2628 /* Reset __FUNCT__ in case the user does not define it themselves */
2629 #undef __FUNCT__
2630 #define __FUNCT__ "User provided function"
2631 
2632 #endif
2633