xref: /petsc/include/petscsys.h (revision 785e854f82a3c614b452fca2cf5ad4f2afe8bdde)
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(), PetscMalloc2()
569 
570   Concepts: memory allocation
571 
572 M*/
573 #define PetscMalloc1(m1,r1) PetscMalloc((m1)*sizeof(**(r1)),r1)
574 
575 /*MC
576    PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN
577 
578    Synopsis:
579     #include "petscsys.h"
580    PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)
581 
582    Not Collective
583 
584    Input Parameter:
585 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
586 -  m2 - number of elements to allocate in 2nd chunk  (may be zero)
587 
588    Output Parameter:
589 +  r1 - memory allocated in first chunk
590 -  r2 - memory allocated in second chunk
591 
592    Level: developer
593 
594 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc1()
595 
596   Concepts: memory allocation
597 
598 M*/
599 #if defined(PETSC_USE_DEBUG)
600 #define PetscMalloc2(m1,r1,m2,r2) (PetscMalloc((m1)*sizeof(**(r1)),r1) || PetscMalloc((m2)*sizeof(**(r2)),r2))
601 #else
602 #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))
603 #endif
604 
605 /*MC
606    PetscMalloc3 - Allocates 3 arrays of memory, all aligned to PETSC_MEMALIGN
607 
608    Synopsis:
609     #include "petscsys.h"
610    PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
611 
612    Not Collective
613 
614    Input Parameter:
615 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
616 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
617 -  m3 - number of elements to allocate in 3rd chunk  (may be zero)
618 
619    Output Parameter:
620 +  r1 - memory allocated in first chunk
621 .  r2 - memory allocated in second chunk
622 -  r3 - memory allocated in third chunk
623 
624    Level: developer
625 
626 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree3()
627 
628   Concepts: memory allocation
629 
630 M*/
631 #if defined(PETSC_USE_DEBUG)
632 #define PetscMalloc3(m1,r1,m2,r2,m3,r3) (PetscMalloc((m1)*sizeof(**(r1)),r1) || PetscMalloc((m2)*sizeof(**(r2)),r2) || PetscMalloc((m3)*sizeof(**(r3)),r3))
633 #else
634 #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)) \
635                                          || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),0))
636 #endif
637 
638 /*MC
639    PetscMalloc4 - Allocates 4 arrays of memory, all aligned to PETSC_MEMALIGN
640 
641    Synopsis:
642     #include "petscsys.h"
643    PetscErrorCode PetscMalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
644 
645    Not Collective
646 
647    Input Parameter:
648 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
649 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
650 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
651 -  m4 - number of elements to allocate in 4th chunk  (may be zero)
652 
653    Output Parameter:
654 +  r1 - memory allocated in first chunk
655 .  r2 - memory allocated in second chunk
656 .  r3 - memory allocated in third chunk
657 -  r4 - memory allocated in fourth chunk
658 
659    Level: developer
660 
661 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree3(), PetscFree4()
662 
663   Concepts: memory allocation
664 
665 M*/
666 #if defined(PETSC_USE_DEBUG)
667 #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))
668 #else
669 #define PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4)                           \
670   ((*(r2) = 0, *(r3) = 0, *(r4) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(m3)*sizeof(**(r3))+(m4)*sizeof(**(r4))+3*(PETSC_MEMALIGN-1),r1)) \
671    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),0))
672 #endif
673 
674 /*MC
675    PetscMalloc5 - Allocates 5 arrays of memory, all aligned to PETSC_MEMALIGN
676 
677    Synopsis:
678     #include "petscsys.h"
679    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)
680 
681    Not Collective
682 
683    Input Parameter:
684 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
685 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
686 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
687 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
688 -  m5 - number of elements to allocate in 5th chunk  (may be zero)
689 
690    Output Parameter:
691 +  r1 - memory allocated in first chunk
692 .  r2 - memory allocated in second chunk
693 .  r3 - memory allocated in third chunk
694 .  r4 - memory allocated in fourth chunk
695 -  r5 - memory allocated in fifth chunk
696 
697    Level: developer
698 
699 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree3(), PetscFree4(), PetscFree5()
700 
701   Concepts: memory allocation
702 
703 M*/
704 #if defined(PETSC_USE_DEBUG)
705 #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))
706 #else
707 #define PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5)      \
708   ((*(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)) \
709    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),*(void**)(r5) = PetscAddrAlign(*(r4)+(m4)),0))
710 #endif
711 
712 
713 /*MC
714    PetscMalloc6 - Allocates 6 arrays of memory, all aligned to PETSC_MEMALIGN
715 
716    Synopsis:
717     #include "petscsys.h"
718    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)
719 
720    Not Collective
721 
722    Input Parameter:
723 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
724 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
725 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
726 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
727 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
728 -  m6 - number of elements to allocate in 6th chunk  (may be zero)
729 
730    Output Parameter:
731 +  r1 - memory allocated in first chunk
732 .  r2 - memory allocated in second chunk
733 .  r3 - memory allocated in third chunk
734 .  r4 - memory allocated in fourth chunk
735 .  r5 - memory allocated in fifth chunk
736 -  r6 - memory allocated in sixth chunk
737 
738    Level: developer
739 
740 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree3(), PetscFree4(), PetscFree5(), PetscFree6()
741 
742   Concepts: memory allocation
743 
744 M*/
745 #if defined(PETSC_USE_DEBUG)
746 #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))
747 #else
748 #define PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) \
749   ((*(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)) \
750    || (*(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))
751 #endif
752 
753 /*MC
754    PetscMalloc7 - Allocates 7 arrays of memory, all aligned to PETSC_MEMALIGN
755 
756    Synopsis:
757     #include "petscsys.h"
758    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)
759 
760    Not Collective
761 
762    Input Parameter:
763 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
764 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
765 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
766 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
767 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
768 .  m6 - number of elements to allocate in 6th chunk  (may be zero)
769 -  m7 - number of elements to allocate in 7th chunk  (may be zero)
770 
771    Output Parameter:
772 +  r1 - memory allocated in first chunk
773 .  r2 - memory allocated in second chunk
774 .  r3 - memory allocated in third chunk
775 .  r4 - memory allocated in fourth chunk
776 .  r5 - memory allocated in fifth chunk
777 .  r6 - memory allocated in sixth chunk
778 -  r7 - memory allocated in seventh chunk
779 
780    Level: developer
781 
782 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree3(), PetscFree4(), PetscFree5(), PetscFree6(), PetscFree7()
783 
784   Concepts: memory allocation
785 
786 M*/
787 #if defined(PETSC_USE_DEBUG)
788 #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))
789 #else
790 #define PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) \
791   ((*(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)) \
792    || (*(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))
793 #endif
794 
795 /*MC
796    PetscNew - Allocates memory of a particular type, zeros the memory! Aligned to PETSC_MEMALIGN
797 
798    Synopsis:
799     #include "petscsys.h"
800    PetscErrorCode PetscNew(struct type,((type *))result)
801 
802    Not Collective
803 
804    Input Parameter:
805 .  type - structure name of space to be allocated. Memory of size sizeof(type) is allocated
806 
807    Output Parameter:
808 .  result - memory allocated
809 
810    Level: beginner
811 
812 .seealso: PetscFree(), PetscMalloc(), PetscNewLog()
813 
814   Concepts: memory allocation
815 
816 M*/
817 #define PetscNew(A,b)      (PetscMalloc(sizeof(A),(b)) || PetscMemzero(*(b),sizeof(A)))
818 
819 /*MC
820    PetscNewLog - Allocates memory of a particular type, zeros the memory! Aligned to PETSC_MEMALIGN. Associates the memory allocated
821          with the given object using PetscLogObjectMemory().
822 
823    Synopsis:
824     #include "petscsys.h"
825    PetscErrorCode PetscNewLog(PetscObject obj,struct type,((type *))result)
826 
827    Not Collective
828 
829    Input Parameter:
830 +  obj - object memory is logged to
831 -  type - structure name of space to be allocated. Memory of size sizeof(type) is allocated
832 
833    Output Parameter:
834 .  result - memory allocated
835 
836    Level: developer
837 
838 .seealso: PetscFree(), PetscMalloc(), PetscNew(), PetscLogObjectMemory()
839 
840   Concepts: memory allocation
841 
842 M*/
843 #define PetscNewLog(o,A,b) (PetscNew(A,b) || ((o) ? PetscLogObjectMemory((PetscObject)o,sizeof(A)) : 0))
844 
845 /*MC
846    PetscFree - Frees memory
847 
848    Synopsis:
849     #include "petscsys.h"
850    PetscErrorCode PetscFree(void *memory)
851 
852    Not Collective
853 
854    Input Parameter:
855 .   memory - memory to free (the pointer is ALWAYS set to 0 upon sucess)
856 
857    Level: beginner
858 
859    Notes: Memory must have been obtained with PetscNew() or PetscMalloc()
860 
861 .seealso: PetscNew(), PetscMalloc(), PetscFreeVoid()
862 
863   Concepts: memory allocation
864 
865 M*/
866 #define PetscFree(a)   ((a) && ((*PetscTrFree)((void*)(a),__LINE__,PETSC_FUNCTION_NAME,__FILE__) || ((a) = 0,0)))
867 
868 /*MC
869    PetscFreeVoid - Frees memory
870 
871    Synopsis:
872     #include "petscsys.h"
873    void PetscFreeVoid(void *memory)
874 
875    Not Collective
876 
877    Input Parameter:
878 .   memory - memory to free
879 
880    Level: beginner
881 
882    Notes: This is different from PetscFree() in that no error code is returned
883 
884 .seealso: PetscFree(), PetscNew(), PetscMalloc()
885 
886   Concepts: memory allocation
887 
888 M*/
889 #define PetscFreeVoid(a) ((*PetscTrFree)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__),(a) = 0)
890 
891 
892 /*MC
893    PetscFree2 - Frees 2 chunks of memory obtained with PetscMalloc2()
894 
895    Synopsis:
896     #include "petscsys.h"
897    PetscErrorCode PetscFree2(void *memory1,void *memory2)
898 
899    Not Collective
900 
901    Input Parameter:
902 +   memory1 - memory to free
903 -   memory2 - 2nd memory to free
904 
905    Level: developer
906 
907    Notes: Memory must have been obtained with PetscMalloc2()
908 
909 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree()
910 
911   Concepts: memory allocation
912 
913 M*/
914 #if defined(PETSC_USE_DEBUG)
915 #define PetscFree2(m1,m2)   (PetscFree(m2) || PetscFree(m1))
916 #else
917 #define PetscFree2(m1,m2)   ((m2)=0, PetscFree(m1))
918 #endif
919 
920 /*MC
921    PetscFree3 - Frees 3 chunks of memory obtained with PetscMalloc3()
922 
923    Synopsis:
924     #include "petscsys.h"
925    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)
926 
927    Not Collective
928 
929    Input Parameter:
930 +   memory1 - memory to free
931 .   memory2 - 2nd memory to free
932 -   memory3 - 3rd memory to free
933 
934    Level: developer
935 
936    Notes: Memory must have been obtained with PetscMalloc3()
937 
938 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3()
939 
940   Concepts: memory allocation
941 
942 M*/
943 #if defined(PETSC_USE_DEBUG)
944 #define PetscFree3(m1,m2,m3)   (PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
945 #else
946 #define PetscFree3(m1,m2,m3)   ((m3)=0,(m2)=0,PetscFree(m1))
947 #endif
948 
949 /*MC
950    PetscFree4 - Frees 4 chunks of memory obtained with PetscMalloc4()
951 
952    Synopsis:
953     #include "petscsys.h"
954    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)
955 
956    Not Collective
957 
958    Input Parameter:
959 +   m1 - memory to free
960 .   m2 - 2nd memory to free
961 .   m3 - 3rd memory to free
962 -   m4 - 4th memory to free
963 
964    Level: developer
965 
966    Notes: Memory must have been obtained with PetscMalloc4()
967 
968 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4()
969 
970   Concepts: memory allocation
971 
972 M*/
973 #if defined(PETSC_USE_DEBUG)
974 #define PetscFree4(m1,m2,m3,m4)   (PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
975 #else
976 #define PetscFree4(m1,m2,m3,m4)   ((m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
977 #endif
978 
979 /*MC
980    PetscFree5 - Frees 5 chunks of memory obtained with PetscMalloc5()
981 
982    Synopsis:
983     #include "petscsys.h"
984    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)
985 
986    Not Collective
987 
988    Input Parameter:
989 +   m1 - memory to free
990 .   m2 - 2nd memory to free
991 .   m3 - 3rd memory to free
992 .   m4 - 4th memory to free
993 -   m5 - 5th memory to free
994 
995    Level: developer
996 
997    Notes: Memory must have been obtained with PetscMalloc5()
998 
999 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5()
1000 
1001   Concepts: memory allocation
1002 
1003 M*/
1004 #if defined(PETSC_USE_DEBUG)
1005 #define PetscFree5(m1,m2,m3,m4,m5)   (PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1006 #else
1007 #define PetscFree5(m1,m2,m3,m4,m5)   ((m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
1008 #endif
1009 
1010 
1011 /*MC
1012    PetscFree6 - Frees 6 chunks of memory obtained with PetscMalloc6()
1013 
1014    Synopsis:
1015     #include "petscsys.h"
1016    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)
1017 
1018    Not Collective
1019 
1020    Input Parameter:
1021 +   m1 - memory to free
1022 .   m2 - 2nd memory to free
1023 .   m3 - 3rd memory to free
1024 .   m4 - 4th memory to free
1025 .   m5 - 5th memory to free
1026 -   m6 - 6th memory to free
1027 
1028 
1029    Level: developer
1030 
1031    Notes: Memory must have been obtained with PetscMalloc6()
1032 
1033 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6()
1034 
1035   Concepts: memory allocation
1036 
1037 M*/
1038 #if defined(PETSC_USE_DEBUG)
1039 #define PetscFree6(m1,m2,m3,m4,m5,m6)   (PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1040 #else
1041 #define PetscFree6(m1,m2,m3,m4,m5,m6)   ((m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
1042 #endif
1043 
1044 /*MC
1045    PetscFree7 - Frees 7 chunks of memory obtained with PetscMalloc7()
1046 
1047    Synopsis:
1048     #include "petscsys.h"
1049    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)
1050 
1051    Not Collective
1052 
1053    Input Parameter:
1054 +   m1 - memory to free
1055 .   m2 - 2nd memory to free
1056 .   m3 - 3rd memory to free
1057 .   m4 - 4th memory to free
1058 .   m5 - 5th memory to free
1059 .   m6 - 6th memory to free
1060 -   m7 - 7th memory to free
1061 
1062 
1063    Level: developer
1064 
1065    Notes: Memory must have been obtained with PetscMalloc7()
1066 
1067 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1068           PetscMalloc7()
1069 
1070   Concepts: memory allocation
1071 
1072 M*/
1073 #if defined(PETSC_USE_DEBUG)
1074 #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   (PetscFree(m7) || PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1075 #else
1076 #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   ((m7)=0,(m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
1077 #endif
1078 
1079 PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,int,const char[],const char[],void**);
1080 PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1081 PETSC_EXTERN PetscErrorCode PetscMallocSet(PetscErrorCode (*)(size_t,int,const char[],const char[],void**),PetscErrorCode (*)(void*,int,const char[],const char[]));
1082 PETSC_EXTERN PetscErrorCode PetscMallocClear(void);
1083 
1084 /*
1085     PetscLogDouble variables are used to contain double precision numbers
1086   that are not used in the numerical computations, but rather in logging,
1087   timing etc.
1088 */
1089 typedef double PetscLogDouble;
1090 #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE
1091 
1092 /*
1093    Routines for tracing memory corruption/bleeding with default PETSc  memory allocation
1094 */
1095 PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1096 PETSC_EXTERN PetscErrorCode PetscMallocDumpLog(FILE *);
1097 PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1098 PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1099 PETSC_EXTERN PetscErrorCode PetscMallocDebug(PetscBool);
1100 PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*);
1101 PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1102 PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLog(void);
1103 PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLogThreshold(PetscLogDouble);
1104 PETSC_EXTERN PetscErrorCode PetscMallocGetDumpLog(PetscBool*);
1105 
1106 /*E
1107     PetscDataType - Used for handling different basic data types.
1108 
1109    Level: beginner
1110 
1111    Developer comment: It would be nice if we could always just use MPI Datatypes, why can we not?
1112 
1113 .seealso: PetscBinaryRead(), PetscBinaryWrite(), PetscDataTypeToMPIDataType(),
1114           PetscDataTypeGetSize()
1115 
1116 E*/
1117 typedef enum {PETSC_INT = 0,PETSC_DOUBLE = 1,PETSC_COMPLEX = 2, PETSC_LONG = 3 ,PETSC_SHORT = 4,PETSC_FLOAT = 5,
1118               PETSC_CHAR = 6,PETSC_BIT_LOGICAL = 7,PETSC_ENUM = 8,PETSC_BOOL=9, PETSC___FLOAT128 = 10,PETSC_OBJECT = 11, PETSC_FUNCTION = 12} PetscDataType;
1119 PETSC_EXTERN const char *const PetscDataTypes[];
1120 
1121 #if defined(PETSC_USE_COMPLEX)
1122 #define  PETSC_SCALAR  PETSC_COMPLEX
1123 #else
1124 #if defined(PETSC_USE_REAL_SINGLE)
1125 #define  PETSC_SCALAR  PETSC_FLOAT
1126 #elif defined(PETSC_USE_REAL___FLOAT128)
1127 #define  PETSC_SCALAR  PETSC___FLOAT128
1128 #else
1129 #define  PETSC_SCALAR  PETSC_DOUBLE
1130 #endif
1131 #endif
1132 #if defined(PETSC_USE_REAL_SINGLE)
1133 #define  PETSC_REAL  PETSC_FLOAT
1134 #elif defined(PETSC_USE_REAL___FLOAT128)
1135 #define  PETSC_REAL  PETSC___FLOAT128
1136 #else
1137 #define  PETSC_REAL  PETSC_DOUBLE
1138 #endif
1139 #define  PETSC_FORTRANADDR  PETSC_LONG
1140 
1141 PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1142 PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1143 PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1144 PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);
1145 
1146 /*
1147     Basic memory and string operations. These are usually simple wrappers
1148    around the basic Unix system calls, but a few of them have additional
1149    functionality and/or error checking.
1150 */
1151 PETSC_EXTERN PetscErrorCode PetscBitMemcpy(void*,PetscInt,const void*,PetscInt,PetscInt,PetscDataType);
1152 PETSC_EXTERN PetscErrorCode PetscMemmove(void*,void *,size_t);
1153 PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1154 PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1155 PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***);
1156 PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1157 PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool  *);
1158 PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool  *);
1159 PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1160 PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1161 PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1162 PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1163 PETSC_EXTERN PetscErrorCode PetscStrncat(char[],const char[],size_t);
1164 PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1165 PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1166 PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1167 PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1168 PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1169 PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1170 PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1171 PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1172 PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*);
1173 PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1174 PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1175 PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1176 PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1177 PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);
1178 
1179 PETSC_EXTERN void PetscStrcmpNoError(const char[],const char[],PetscBool  *);
1180 
1181 /*S
1182     PetscToken - 'Token' used for managing tokenizing strings
1183 
1184   Level: intermediate
1185 
1186 .seealso: PetscTokenCreate(), PetscTokenFind(), PetscTokenDestroy()
1187 S*/
1188 typedef struct _p_PetscToken* PetscToken;
1189 
1190 PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1191 PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1192 PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);
1193 
1194 PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1195 PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);
1196 
1197 /*
1198    These are  MPI operations for MPI_Allreduce() etc
1199 */
1200 PETSC_EXTERN MPI_Op PetscMaxSum_Op;
1201 #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128)
1202 PETSC_EXTERN MPI_Op MPIU_SUM;
1203 #else
1204 #define MPIU_SUM MPI_SUM
1205 #endif
1206 #if defined(PETSC_USE_REAL___FLOAT128)
1207 PETSC_EXTERN MPI_Op MPIU_MAX;
1208 PETSC_EXTERN MPI_Op MPIU_MIN;
1209 #else
1210 #define MPIU_MAX MPI_MAX
1211 #define MPIU_MIN MPI_MIN
1212 #endif
1213 PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);
1214 
1215 PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1216 PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1217 
1218 /*S
1219      PetscObject - any PETSc object, PetscViewer, Mat, Vec, KSP etc
1220 
1221    Level: beginner
1222 
1223    Note: This is the base class from which all PETSc objects are derived from.
1224 
1225 .seealso:  PetscObjectDestroy(), PetscObjectView(), PetscObjectGetName(), PetscObjectSetName(), PetscObjectReference(), PetscObjectDereferenc()
1226 S*/
1227 typedef struct _p_PetscObject* PetscObject;
1228 
1229 /*MC
1230     PetscObjectId - unique integer Id for a PetscObject
1231 
1232     Level: developer
1233 
1234     Notes: Unlike pointer values, object ids are never reused.
1235 
1236 .seealso: PetscObjectState, PetscObjectGetId()
1237 M*/
1238 typedef Petsc64bitInt PetscObjectId;
1239 
1240 /*MC
1241     PetscObjectState - integer state for a PetscObject
1242 
1243     Level: developer
1244 
1245     Notes:
1246     Object state is always-increasing and (for objects that track state) can be used to determine if an object has
1247     changed since the last time you interacted with it.  It is 64-bit so that it will not overflow for a very long time.
1248 
1249 .seealso: PetscObjectId, PetscObjectStateGet(), PetscObjectStateIncrease(), PetscObjectStateSet()
1250 M*/
1251 typedef Petsc64bitInt PetscObjectState;
1252 
1253 /*S
1254      PetscFunctionList - Linked list of functions, possibly stored in dynamic libraries, accessed
1255       by string name
1256 
1257    Level: advanced
1258 
1259 .seealso:  PetscFunctionListAdd(), PetscFunctionListDestroy(), PetscOpFlist
1260 S*/
1261 typedef struct _n_PetscFunctionList *PetscFunctionList;
1262 
1263 /*E
1264   PetscFileMode - Access mode for a file.
1265 
1266   Level: beginner
1267 
1268   FILE_MODE_READ - open a file at its beginning for reading
1269 
1270   FILE_MODE_WRITE - open a file at its beginning for writing (will create if the file does not exist)
1271 
1272   FILE_MODE_APPEND - open a file at end for writing
1273 
1274   FILE_MODE_UPDATE - open a file for updating, meaning for reading and writing
1275 
1276   FILE_MODE_APPEND_UPDATE - open a file for updating, meaning for reading and writing, at the end
1277 
1278 .seealso: PetscViewerFileSetMode()
1279 E*/
1280 typedef enum {FILE_MODE_READ, FILE_MODE_WRITE, FILE_MODE_APPEND, FILE_MODE_UPDATE, FILE_MODE_APPEND_UPDATE} PetscFileMode;
1281 extern const char *const PetscFileModes[];
1282 
1283 /*
1284     Defines PETSc error handling.
1285 */
1286 #include <petscerror.h>
1287 
1288 #define PETSC_SMALLEST_CLASSID  1211211
1289 PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1290 PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1291 PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);
1292 
1293 /*
1294    Routines that get memory usage information from the OS
1295 */
1296 PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1297 PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1298 PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1299 
1300 PETSC_EXTERN PetscErrorCode PetscInfoAllow(PetscBool ,const char []);
1301 PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);
1302 
1303 /*
1304    Initialization of PETSc
1305 */
1306 PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1307 PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1308 PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1309 PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1310 PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1311 PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1312 PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1313 PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1314 PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1315 PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);
1316 
1317 PETSC_EXTERN PetscErrorCode PetscEnd(void);
1318 PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);
1319 
1320 PETSC_EXTERN MPI_Comm PETSC_COMM_LOCAL_WORLD;
1321 PETSC_EXTERN PetscErrorCode PetscHMPIMerge(PetscMPIInt,PetscErrorCode (*)(void*),void*);
1322 PETSC_EXTERN PetscErrorCode PetscHMPISpawn(PetscMPIInt);
1323 PETSC_EXTERN PetscErrorCode PetscHMPIFinalize(void);
1324 PETSC_EXTERN PetscErrorCode PetscHMPIRun(MPI_Comm,PetscErrorCode (*)(MPI_Comm,void *),void*);
1325 PETSC_EXTERN PetscErrorCode PetscHMPIRunCtx(MPI_Comm,PetscErrorCode (*)(MPI_Comm,void*,void *),void*);
1326 PETSC_EXTERN PetscErrorCode PetscHMPIFree(MPI_Comm,void*);
1327 PETSC_EXTERN PetscErrorCode PetscHMPIMalloc(MPI_Comm,size_t,void**);
1328 
1329 PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1330 PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1331 PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1332 PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);
1333 
1334 /*
1335      These are so that in extern C code we can caste function pointers to non-extern C
1336    function pointers. Since the regular C++ code expects its function pointers to be C++
1337 */
1338 PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1339 PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1340 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);
1341 
1342 /*
1343     Functions that can act on any PETSc object.
1344 */
1345 PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1346 PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1347 PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1348 PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1349 PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1350 PETSC_EXTERN PetscErrorCode PetscObjectSetPrecision(PetscObject,PetscPrecision);
1351 PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1352 PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1353 PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1354 PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1355 PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1356 PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1357 PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1358 PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1359 PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1360 PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt *);
1361 PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1362 PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1363 PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject *);
1364 PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void));
1365 #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d))
1366 PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1367 PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1368 PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt *);
1369 PETSC_EXTERN PetscErrorCode PetscObjectAddOptionsHandler(PetscObject,PetscErrorCode (*)(PetscObject,void*),PetscErrorCode (*)(PetscObject,void*),void*);
1370 PETSC_EXTERN PetscErrorCode PetscObjectProcessOptionsHandlers(PetscObject);
1371 PETSC_EXTERN PetscErrorCode PetscObjectDestroyOptionsHandlers(PetscObject);
1372 PETSC_EXTERN PetscErrorCode PetscObjectsGetGlobalNumbering(MPI_Comm,PetscInt,PetscObject*,PetscInt*,PetscInt*);
1373 
1374 #include <petscviewertypes.h>
1375 #include <petscoptions.h>
1376 
1377 PETSC_EXTERN PetscErrorCode PetscMemoryShowUsage(PetscViewer,const char[]);
1378 PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1379 PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1380 #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1381 PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1382 PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1383 PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1384 PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1385 PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1386 PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1387 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1388 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1389 PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1390 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1391 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1392 PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1393 PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);
1394 
1395 #if defined(PETSC_HAVE_SAWS)
1396 PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1397 PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1398 PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1399 PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1400 PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1401 PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1402 PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1403 PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1404 PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1405 PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);
1406 
1407 #else
1408 #define PetscSAWsBlock()                        0
1409 #define PetscObjectSAWsViewOff(obj)             0
1410 #define PetscObjectSAWsSetBlock(obj,flg)        0
1411 #define PetscObjectSAWsBlock(obj)               0
1412 #define PetscObjectSAWsGrantAccess(obj)         0
1413 #define PetscObjectSAWsTakeAccess(obj)          0
1414 #define PetscStackViewSAWs()                    0
1415 #define PetscStackSAWsViewOff()                 0
1416 #define PetscStackSAWsTakeAccess()
1417 #define PetscStackSAWsGrantAccess()
1418 
1419 #endif
1420 
1421 typedef void* PetscDLHandle;
1422 typedef enum {PETSC_DL_DECIDE=0,PETSC_DL_NOW=1,PETSC_DL_LOCAL=2} PetscDLMode;
1423 PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1424 PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1425 PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);
1426 
1427 
1428 #if defined(PETSC_USE_DEBUG)
1429 PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1430 #endif
1431 PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);
1432 
1433 /*S
1434      PetscObjectList - Linked list of PETSc objects, each accessable by string name
1435 
1436    Level: developer
1437 
1438    Notes: Used by PetscObjectCompose() and PetscObjectQuery()
1439 
1440 .seealso:  PetscObjectListAdd(), PetscObjectListDestroy(), PetscObjectListFind(), PetscObjectCompose(), PetscObjectQuery(), PetscFunctionList
1441 S*/
1442 typedef struct _n_PetscObjectList *PetscObjectList;
1443 
1444 PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1445 PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1446 PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1447 PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1448 PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1449 PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);
1450 
1451 /*
1452     Dynamic library lists. Lists of names of routines in objects or in dynamic
1453   link libraries that will be loaded as needed.
1454 */
1455 
1456 #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1457 PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1458 PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1459 #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1460 PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1461 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[]);
1462 PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1463 PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1464 PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);
1465 
1466 /*S
1467      PetscDLLibrary - Linked list of dynamics libraries to search for functions
1468 
1469    Level: advanced
1470 
1471 .seealso:  PetscDLLibraryOpen()
1472 S*/
1473 typedef struct _n_PetscDLLibrary *PetscDLLibrary;
1474 PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1475 PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1476 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1477 PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1478 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1479 PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1480 PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1481 PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);
1482 
1483 /*
1484      Useful utility routines
1485 */
1486 PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1487 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1488 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1489 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1490 PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1491 PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);
1492 
1493 /*
1494     PetscNot - negates a logical type value and returns result as a PetscBool
1495 
1496     Notes: This is useful in cases like
1497 $     int        *a;
1498 $     PetscBool  flag = PetscNot(a)
1499      where !a does not return a PetscBool  because we cannot provide a cast from int to PetscBool  in C.
1500 */
1501 #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)
1502 
1503 #if defined(PETSC_HAVE_VALGRIND)
1504 #  include <valgrind/valgrind.h>
1505 #  define PETSC_RUNNING_ON_VALGRIND RUNNING_ON_VALGRIND
1506 #else
1507 #  define PETSC_RUNNING_ON_VALGRIND PETSC_FALSE
1508 #endif
1509 
1510 
1511 /*MC
1512     PetscHelpPrintf - Prints help messages.
1513 
1514    Synopsis:
1515     #include "petscsys.h"
1516      PetscErrorCode (*PetscHelpPrintf)(const char format[],...);
1517 
1518     Not Collective
1519 
1520     Input Parameters:
1521 .   format - the usual printf() format string
1522 
1523    Level: developer
1524 
1525     Fortran Note:
1526     This routine is not supported in Fortran.
1527 
1528     Concepts: help messages^printing
1529     Concepts: printing^help messages
1530 
1531 .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1532 M*/
1533 PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);
1534 
1535 /*
1536      Defines PETSc profiling.
1537 */
1538 #include <petsclog.h>
1539 
1540 
1541 
1542 /*
1543       Simple PETSc parallel IO for ASCII printing
1544 */
1545 PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1546 PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1547 PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1548 PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1549 PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1550 PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1551 PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);
1552 
1553 /* These are used internally by PETSc ASCII IO routines*/
1554 #include <stdarg.h>
1555 PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list);
1556 PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE*,const char[],va_list);
1557 PETSC_EXTERN PetscErrorCode PetscVFPrintfDefault(FILE*,const char[],va_list);
1558 
1559 #if defined(PETSC_HAVE_MATLAB_ENGINE)
1560 PETSC_EXTERN PetscErrorCode PetscVFPrintf_Matlab(FILE*,const char[],va_list);
1561 #endif
1562 
1563 PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1564 PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1565 PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);
1566 
1567 #if defined(PETSC_HAVE_POPEN)
1568 PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1569 PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*,int*);
1570 PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1571 #endif
1572 
1573 PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1574 PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1575 PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*);
1576 PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1577 PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1578 PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1579 PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);
1580 
1581 PETSC_EXTERN PetscErrorCode PetscPopUpSelect(MPI_Comm,const char*,const char*,int,const char**,int*);
1582 
1583 /*S
1584      PetscContainer - Simple PETSc object that contains a pointer to any required data
1585 
1586    Level: advanced
1587 
1588 .seealso:  PetscObject, PetscContainerCreate()
1589 S*/
1590 PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1591 typedef struct _p_PetscContainer*  PetscContainer;
1592 PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void **);
1593 PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void *);
1594 PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1595 PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer *);
1596 PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void*));
1597 
1598 /*
1599    For use in debuggers
1600 */
1601 PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1602 PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1603 PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1604 PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1605 PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);
1606 
1607 #include <stddef.h>
1608 #if defined(PETSC_HAVE_MEMORY_H)
1609 #include <memory.h>
1610 #endif
1611 #if defined(PETSC_HAVE_STDLIB_H)
1612 #include <stdlib.h>
1613 #endif
1614 
1615 #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1616 #include <xmmintrin.h>
1617 #endif
1618 
1619 #undef __FUNCT__
1620 #define __FUNCT__ "PetscMemcpy"
1621 /*@C
1622    PetscMemcpy - Copies n bytes, beginning at location b, to the space
1623    beginning at location a. The two memory regions CANNOT overlap, use
1624    PetscMemmove() in that case.
1625 
1626    Not Collective
1627 
1628    Input Parameters:
1629 +  b - pointer to initial memory space
1630 -  n - length (in bytes) of space to copy
1631 
1632    Output Parameter:
1633 .  a - pointer to copy space
1634 
1635    Level: intermediate
1636 
1637    Compile Option:
1638     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1639                                   for memory copies on double precision values.
1640     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1641                                   for memory copies on double precision values.
1642     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1643                                   for memory copies on double precision values.
1644 
1645    Note:
1646    This routine is analogous to memcpy().
1647 
1648    Developer Note: this is inlined for fastest performance
1649 
1650   Concepts: memory^copying
1651   Concepts: copying^memory
1652 
1653 .seealso: PetscMemmove()
1654 
1655 @*/
1656 PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n)
1657 {
1658 #if defined(PETSC_USE_DEBUG)
1659   unsigned long al = (unsigned long) a,bl = (unsigned long) b;
1660   unsigned long nl = (unsigned long) n;
1661   PetscFunctionBegin;
1662   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1663   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1664 #else
1665   PetscFunctionBegin;
1666 #endif
1667   if (a != b) {
1668 #if defined(PETSC_USE_DEBUG)
1669     if ((al > bl && (al - bl) < nl) || (bl - al) < nl)  SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1670               or make sure your copy regions and lengths are correct. \n\
1671               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1672 #endif
1673 #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1674    if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1675       size_t len = n/sizeof(PetscScalar);
1676 #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1677       PetscBLASInt   one = 1,blen;
1678       PetscErrorCode ierr;
1679       ierr = PetscBLASIntCast(len,&blen);CHKERRQ(ierr);
1680       PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one));
1681 #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1682       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1683 #else
1684       size_t      i;
1685       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1686       for (i=0; i<len; i++) y[i] = x[i];
1687 #endif
1688     } else {
1689       memcpy((char*)(a),(char*)(b),n);
1690     }
1691 #else
1692     memcpy((char*)(a),(char*)(b),n);
1693 #endif
1694   }
1695   PetscFunctionReturn(0);
1696 }
1697 
1698 /*@C
1699    PetscMemzero - Zeros the specified memory.
1700 
1701    Not Collective
1702 
1703    Input Parameters:
1704 +  a - pointer to beginning memory location
1705 -  n - length (in bytes) of memory to initialize
1706 
1707    Level: intermediate
1708 
1709    Compile Option:
1710    PETSC_PREFER_BZERO - on certain machines (the IBM RS6000) the bzero() routine happens
1711   to be faster than the memset() routine. This flag causes the bzero() routine to be used.
1712 
1713    Developer Note: this is inlined for fastest performance
1714 
1715    Concepts: memory^zeroing
1716    Concepts: zeroing^memory
1717 
1718 .seealso: PetscMemcpy()
1719 @*/
1720 PETSC_STATIC_INLINE PetscErrorCode  PetscMemzero(void *a,size_t n)
1721 {
1722   if (n > 0) {
1723 #if defined(PETSC_USE_DEBUG)
1724     if (!a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer");
1725 #endif
1726 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1727     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1728       size_t      i,len = n/sizeof(PetscScalar);
1729       PetscScalar *x = (PetscScalar*)a;
1730       for (i=0; i<len; i++) x[i] = 0.0;
1731     } else {
1732 #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1733     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1734       PetscInt len = n/sizeof(PetscScalar);
1735       fortranzero_(&len,(PetscScalar*)a);
1736     } else {
1737 #endif
1738 #if defined(PETSC_PREFER_BZERO)
1739       bzero((char *)a,n);
1740 #else
1741       memset((char*)a,0,n);
1742 #endif
1743 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1744     }
1745 #endif
1746   }
1747   return 0;
1748 }
1749 
1750 /*MC
1751    PetscPrefetchBlock - Prefetches a block of memory
1752 
1753    Synopsis:
1754     #include "petscsys.h"
1755     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)
1756 
1757    Not Collective
1758 
1759    Input Parameters:
1760 +  a - pointer to first element to fetch (any type but usually PetscInt or PetscScalar)
1761 .  n - number of elements to fetch
1762 .  rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors)
1763 -  t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note
1764 
1765    Level: developer
1766 
1767    Notes:
1768    The last two arguments (rw and t) must be compile-time constants.
1769 
1770    Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality.  Not all architectures offer
1771    equivalent locality hints, but the following macros are always defined to their closest analogue.
1772 +  PETSC_PREFETCH_HINT_NTA - Non-temporal.  Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched).
1773 .  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.
1774 .  PETSC_PREFETCH_HINT_T1 - Fetch to level 2 and higher (not L1).
1775 -  PETSC_PREFETCH_HINT_T2 - Fetch to high-level cache only.  (On many systems, T0 and T1 are equivalent.)
1776 
1777    This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid
1778    address).
1779 
1780    Concepts: memory
1781 M*/
1782 #define PetscPrefetchBlock(a,n,rw,t) do {                               \
1783     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
1784     for ( ; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
1785   } while (0)
1786 
1787 /*
1788       Determine if some of the kernel computation routines use
1789    Fortran (rather than C) for the numerical calculations. On some machines
1790    and compilers (like complex numbers) the Fortran version of the routines
1791    is faster than the C/C++ versions. The flag --with-fortran-kernels
1792    should be used with ./configure to turn these on.
1793 */
1794 #if defined(PETSC_USE_FORTRAN_KERNELS)
1795 
1796 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1797 #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1798 #endif
1799 
1800 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
1801 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
1802 #endif
1803 
1804 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1805 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
1806 #endif
1807 
1808 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1809 #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
1810 #endif
1811 
1812 #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
1813 #define PETSC_USE_FORTRAN_KERNEL_NORM
1814 #endif
1815 
1816 #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
1817 #define PETSC_USE_FORTRAN_KERNEL_MAXPY
1818 #endif
1819 
1820 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1821 #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
1822 #endif
1823 
1824 #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1825 #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
1826 #endif
1827 
1828 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
1829 #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
1830 #endif
1831 
1832 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1833 #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
1834 #endif
1835 
1836 #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
1837 #define PETSC_USE_FORTRAN_KERNEL_MDOT
1838 #endif
1839 
1840 #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
1841 #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
1842 #endif
1843 
1844 #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
1845 #define PETSC_USE_FORTRAN_KERNEL_AYPX
1846 #endif
1847 
1848 #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
1849 #define PETSC_USE_FORTRAN_KERNEL_WAXPY
1850 #endif
1851 
1852 #endif
1853 
1854 /*
1855     Macros for indicating code that should be compiled with a C interface,
1856    rather than a C++ interface. Any routines that are dynamically loaded
1857    (such as the PCCreate_XXX() routines) must be wrapped so that the name
1858    mangler does not change the functions symbol name. This just hides the
1859    ugly extern "C" {} wrappers.
1860 */
1861 #if defined(__cplusplus)
1862 #define EXTERN_C_BEGIN extern "C" {
1863 #define EXTERN_C_END }
1864 #else
1865 #define EXTERN_C_BEGIN
1866 #define EXTERN_C_END
1867 #endif
1868 
1869 /* --------------------------------------------------------------------*/
1870 
1871 /*MC
1872     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
1873         communication
1874 
1875    Level: beginner
1876 
1877    Note: This manual page is a place-holder because MPICH does not have a manual page for MPI_Comm
1878 
1879 .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
1880 M*/
1881 
1882 /*MC
1883     PetscScalar - PETSc type that represents either a double precision real number, a double precision
1884        complex number, a single precision real number, a long double or an int - if the code is configured
1885        with --with-scalar-type=real,complex --with-precision=single,double,longdouble,int,matsingle
1886 
1887 
1888    Level: beginner
1889 
1890 .seealso: PetscReal, PassiveReal, PassiveScalar, MPIU_SCALAR, PetscInt
1891 M*/
1892 
1893 /*MC
1894     PetscComplex - PETSc type that represents a complex number with precision matching that of PetscReal.
1895 
1896    Synopsis:
1897    #define PETSC_DESIRE_COMPLEX
1898    #include <petscsys.h>
1899    PetscComplex number = 1. + 2.*PETSC_i;
1900 
1901    Level: beginner
1902 
1903    Note:
1904    Complex numbers are automatically available if PETSc was configured --with-scalar-type=complex (in which case
1905    PetscComplex will match PetscScalar), otherwise the macro PETSC_DESIRE_COMPLEX must be defined before including any
1906    PETSc headers. If the compiler supports complex numbers, PetscComplex and associated variables and functions will be
1907    defined and PETSC_HAVE_COMPLEX will be set.
1908 
1909 .seealso: PetscReal, PetscComplex, MPIU_COMPLEX, PetscInt, PETSC_i
1910 M*/
1911 
1912 /*MC
1913     PetscReal - PETSc type that represents a real number version of PetscScalar
1914 
1915    Level: beginner
1916 
1917 .seealso: PetscScalar, PassiveReal, PassiveScalar
1918 M*/
1919 
1920 /*MC
1921     PassiveScalar - PETSc type that represents a PetscScalar
1922    Level: beginner
1923 
1924     This is the same as a PetscScalar except in code that is automatically differentiated it is
1925    treated as a constant (not an indendent or dependent variable)
1926 
1927 .seealso: PetscReal, PassiveReal, PetscScalar
1928 M*/
1929 
1930 /*MC
1931     PassiveReal - PETSc type that represents a PetscReal
1932 
1933    Level: beginner
1934 
1935     This is the same as a PetscReal except in code that is automatically differentiated it is
1936    treated as a constant (not an indendent or dependent variable)
1937 
1938 .seealso: PetscScalar, PetscReal, PassiveScalar
1939 M*/
1940 
1941 /*MC
1942     MPIU_SCALAR - MPI datatype corresponding to PetscScalar
1943 
1944    Level: beginner
1945 
1946     Note: In MPI calls that require an MPI datatype that matches a PetscScalar or array of PetscScalars
1947           pass this value
1948 
1949 .seealso: PetscReal, PassiveReal, PassiveScalar, PetscScalar, MPIU_INT
1950 M*/
1951 
1952 #if defined(PETSC_HAVE_MPIIO)
1953 #if !defined(PETSC_WORDS_BIGENDIAN)
1954 PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1955 PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1956 #else
1957 #define MPIU_File_write_all(a,b,c,d,e) MPI_File_write_all(a,b,c,d,e)
1958 #define MPIU_File_read_all(a,b,c,d,e) MPI_File_read_all(a,b,c,d,e)
1959 #endif
1960 #endif
1961 
1962 /* the following petsc_static_inline require petscerror.h */
1963 
1964 /* Limit MPI to 32-bits */
1965 #define PETSC_MPI_INT_MAX  2147483647
1966 #define PETSC_MPI_INT_MIN -2147483647
1967 /* Limit BLAS to 32-bits */
1968 #define PETSC_BLAS_INT_MAX  2147483647
1969 #define PETSC_BLAS_INT_MIN -2147483647
1970 
1971 #undef __FUNCT__
1972 #define __FUNCT__ "PetscBLASIntCast"
1973 PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
1974 {
1975   PetscFunctionBegin;
1976 #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
1977   if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK");
1978 #endif
1979   *b =  (PetscBLASInt)(a);
1980   PetscFunctionReturn(0);
1981 }
1982 
1983 #undef __FUNCT__
1984 #define __FUNCT__ "PetscMPIIntCast"
1985 PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
1986 {
1987   PetscFunctionBegin;
1988 #if defined(PETSC_USE_64BIT_INDICES)
1989   if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for MPI");
1990 #endif
1991   *b =  (PetscMPIInt)(a);
1992   PetscFunctionReturn(0);
1993 }
1994 
1995 
1996 /*
1997      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
1998 */
1999 #if defined(hz)
2000 #undef hz
2001 #endif
2002 
2003 /*  For arrays that contain filenames or paths */
2004 
2005 
2006 #if defined(PETSC_HAVE_LIMITS_H)
2007 #include <limits.h>
2008 #endif
2009 #if defined(PETSC_HAVE_SYS_PARAM_H)
2010 #include <sys/param.h>
2011 #endif
2012 #if defined(PETSC_HAVE_SYS_TYPES_H)
2013 #include <sys/types.h>
2014 #endif
2015 #if defined(MAXPATHLEN)
2016 #  define PETSC_MAX_PATH_LEN     MAXPATHLEN
2017 #elif defined(MAX_PATH)
2018 #  define PETSC_MAX_PATH_LEN     MAX_PATH
2019 #elif defined(_MAX_PATH)
2020 #  define PETSC_MAX_PATH_LEN     _MAX_PATH
2021 #else
2022 #  define PETSC_MAX_PATH_LEN     4096
2023 #endif
2024 
2025 /* Special support for C++ */
2026 #if defined(PETSC_CLANGUAGE_CXX) && defined(__cplusplus)
2027 #include <petscsys.hh>
2028 #endif
2029 
2030 /*MC
2031 
2032     UsingFortran - Fortran can be used with PETSc in four distinct approaches
2033 
2034 $    1) classic Fortran 77 style
2035 $#include "finclude/petscXXX.h" to work with material from the XXX component of PETSc
2036 $       XXX variablename
2037 $      You cannot use this approach if you wish to use the Fortran 90 specific PETSc routines
2038 $      which end in F90; such as VecGetArrayF90()
2039 $
2040 $    2) classic Fortran 90 style
2041 $#include "finclude/petscXXX.h"
2042 $#include "finclude/petscXXX.h90" to work with material from the XXX component of PETSc
2043 $       XXX variablename
2044 $
2045 $    3) Using Fortran modules
2046 $#include "finclude/petscXXXdef.h"
2047 $         use petscXXXX
2048 $       XXX variablename
2049 $
2050 $    4) Use Fortran modules and Fortran data types for PETSc types
2051 $#include "finclude/petscXXXdef.h"
2052 $         use petscXXXX
2053 $       type(XXX) variablename
2054 $      To use this approach you must ./configure PETSc with the additional
2055 $      option --with-fortran-datatypes You cannot use the type(XXX) declaration approach without using Fortran modules
2056 
2057     Finally if you absolutely do not want to use any #include you can use either
2058 
2059 $    3a) skip the #include BUT you cannot use any PETSc data type names like Vec, Mat, PetscInt, PetscErrorCode etc
2060 $        and you must declare the variables as integer, for example
2061 $        integer variablename
2062 $
2063 $    4a) skip the #include, you use the object types like type(Vec) type(Mat) but cannot use the data type
2064 $        names like PetscErrorCode, PetscInt etc. again for those you must use integer
2065 
2066    We recommend either 2 or 3. Approaches 2 and 3 provide type checking for most PETSc function calls; 4 has type checking
2067 for only a few PETSc functions.
2068 
2069    Fortran type checking with interfaces is strick, this means you cannot pass a scalar value when an array value
2070 is expected (even though it is legal Fortran). For example when setting a single value in a matrix with MatSetValues()
2071 you cannot have something like
2072 $      PetscInt row,col
2073 $      PetscScalar val
2074 $        ...
2075 $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)
2076 You must instead have
2077 $      PetscInt row(1),col(1)
2078 $      PetscScalar val(1)
2079 $        ...
2080 $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)
2081 
2082 
2083     See the example src/vec/vec/examples/tutorials/ex20f90.F90 for an example that can use all four approaches
2084 
2085     Developer Notes: The finclude/petscXXXdef.h contain all the #defines (would be typedefs in C code) these
2086      automatically include their predecessors; for example finclude/petscvecdef.h includes finclude/petscisdef.h
2087 
2088      The finclude/petscXXXX.h contain all the parameter statements for that package. These automatically include
2089      their finclude/petscXXXdef.h file but DO NOT automatically include their predecessors;  for example
2090      finclude/petscvec.h does NOT automatically include finclude/petscis.h
2091 
2092      The finclude/ftn-custom/petscXXXdef.h90 are not intended to be used directly in code, they define the
2093      Fortran data type type(XXX) (for example type(Vec)) when PETSc is ./configure with the --with-fortran-datatypes option.
2094 
2095      The finclude/ftn-custom/petscXXX.h90 (not included directly by code) contain interface definitions for
2096      the PETSc Fortran stubs that have different bindings then their C version (for example VecGetArrayF90).
2097 
2098      The finclude/ftn-auto/petscXXX.h90 (not included directly by code) contain interface definitions generated
2099      automatically by "make allfortranstubs".
2100 
2101      The finclude/petscXXX.h90 includes the custom finclude/ftn-custom/petscXXX.h90 and if ./configure
2102      was run with --with-fortran-interfaces it also includes the finclude/ftn-auto/petscXXX.h90 These DO NOT automatically
2103      include their predecessors
2104 
2105     Level: beginner
2106 
2107 M*/
2108 
2109 PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2110 PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2111 PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2112 PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2113 PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2114 PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2115 PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2116 
2117 PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2118 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2119 PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2120 PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2121 PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2122 PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2123 PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt*,PetscInt*,PetscInt*);
2124 PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2125 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2126 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2127 PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2128 PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2129 PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2130 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2131 PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2132 PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2133 PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2134 PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt*,const PetscInt*,PetscInt,const PetscInt*,const PetscInt*,PetscInt*,PetscInt**,PetscInt**);
2135 
2136 PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2137 PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);
2138 
2139 /*J
2140     PetscRandomType - String with the name of a PETSc randomizer
2141 
2142    Level: beginner
2143 
2144    Notes: to use the SPRNG you must have ./configure PETSc
2145    with the option --download-sprng
2146 
2147 .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2148 J*/
2149 typedef const char* PetscRandomType;
2150 #define PETSCRAND       "rand"
2151 #define PETSCRAND48     "rand48"
2152 #define PETSCSPRNG      "sprng"
2153 
2154 /* Logging support */
2155 PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;
2156 
2157 PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);
2158 
2159 /*S
2160      PetscRandom - Abstract PETSc object that manages generating random numbers
2161 
2162    Level: intermediate
2163 
2164   Concepts: random numbers
2165 
2166 .seealso:  PetscRandomCreate(), PetscRandomGetValue(), PetscRandomType
2167 S*/
2168 typedef struct _p_PetscRandom*   PetscRandom;
2169 
2170 /* Dynamic creation and loading functions */
2171 PETSC_EXTERN PetscFunctionList PetscRandomList;
2172 PETSC_EXTERN PetscBool         PetscRandomRegisterAllCalled;
2173 
2174 PETSC_EXTERN PetscErrorCode PetscRandomRegisterAll(void);
2175 PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2176 PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2177 PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2178 PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType*);
2179  PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom,const char[],const char[]);
2180 PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);
2181 
2182 PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2183 PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2184 PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2185 PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2186 PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2187 PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2188 PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2189 PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2190 PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);
2191 
2192 PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2193 PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2194 PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2195 PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2196 PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2197 PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2198 PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);
2199 
2200 PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscDataType);
2201 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscDataType);
2202 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,void*,PetscInt,PetscDataType,PetscBool );
2203 PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,void*,PetscInt,PetscDataType,PetscBool );
2204 PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2205 PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2206 PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2207 PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2208 PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2209 PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2210 PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2211 PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);
2212 PETSC_EXTERN PetscErrorCode PetscWebServe(MPI_Comm,int);
2213 
2214 /*
2215    In binary files variables are stored using the following lengths,
2216   regardless of how they are stored in memory on any one particular
2217   machine. Use these rather then sizeof() in computing sizes for
2218   PetscBinarySeek().
2219 */
2220 #define PETSC_BINARY_INT_SIZE   (32/8)
2221 #define PETSC_BINARY_FLOAT_SIZE  (32/8)
2222 #define PETSC_BINARY_CHAR_SIZE  (8/8)
2223 #define PETSC_BINARY_SHORT_SIZE  (16/8)
2224 #define PETSC_BINARY_DOUBLE_SIZE  (64/8)
2225 #define PETSC_BINARY_SCALAR_SIZE  sizeof(PetscScalar)
2226 
2227 /*E
2228   PetscBinarySeekType - argument to PetscBinarySeek()
2229 
2230   Level: advanced
2231 
2232 .seealso: PetscBinarySeek(), PetscBinarySynchronizedSeek()
2233 E*/
2234 typedef enum {PETSC_BINARY_SEEK_SET = 0,PETSC_BINARY_SEEK_CUR = 1,PETSC_BINARY_SEEK_END = 2} PetscBinarySeekType;
2235 PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2236 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2237 PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);
2238 
2239 PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2240 PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool );
2241 PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2242 PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2243 PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2244 PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2245 
2246 PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2247 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2248 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2249 PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2250 PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2251 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscInt,const PetscMPIInt*,const void*,PetscInt*,PetscMPIInt**,void*) PetscAttrMPIPointerWithType(6,3);
2252 
2253 /*E
2254     PetscBuildTwoSidedType - algorithm for setting up two-sided communication
2255 
2256 $  PETSC_BUILDTWOSIDED_ALLREDUCE - classical algorithm using an MPI_Allreduce with
2257 $      a buffer of length equal to the communicator size. Not memory-scalable due to
2258 $      the large reduction size. Requires only MPI-1.
2259 $  PETSC_BUILDTWOSIDED_IBARRIER - nonblocking algorithm based on MPI_Issend and MPI_Ibarrier.
2260 $      Proved communication-optimal in Hoefler, Siebert, and Lumsdaine (2010). Requires MPI-3.
2261 
2262    Level: developer
2263 
2264 .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedSetType(), PetscCommBuildTwoSidedGetType()
2265 E*/
2266 typedef enum {
2267   PETSC_BUILDTWOSIDED_NOTSET = -1,
2268   PETSC_BUILDTWOSIDED_ALLREDUCE = 0,
2269   PETSC_BUILDTWOSIDED_IBARRIER = 1
2270   /* Updates here must be accompanied by updates in finclude/petscsys.h and the string array in mpits.c */
2271 } PetscBuildTwoSidedType;
2272 PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2273 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2274 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);
2275 
2276 PETSC_EXTERN PetscErrorCode PetscSSEIsEnabled(MPI_Comm,PetscBool  *,PetscBool  *);
2277 
2278 /*E
2279   InsertMode - Whether entries are inserted or added into vectors or matrices
2280 
2281   Level: beginner
2282 
2283 .seealso: VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2284           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(),
2285           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd()
2286 E*/
2287  typedef enum {NOT_SET_VALUES, INSERT_VALUES, ADD_VALUES, MAX_VALUES, INSERT_ALL_VALUES, ADD_ALL_VALUES, INSERT_BC_VALUES, ADD_BC_VALUES} InsertMode;
2288 
2289 /*MC
2290     INSERT_VALUES - Put a value into a vector or matrix, overwrites any previous value
2291 
2292     Level: beginner
2293 
2294 .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2295           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), ADD_VALUES,
2296           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES
2297 
2298 M*/
2299 
2300 /*MC
2301     ADD_VALUES - Adds a value into a vector or matrix, if there previously was no value, just puts the
2302                 value into that location
2303 
2304     Level: beginner
2305 
2306 .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2307           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), INSERT_VALUES,
2308           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES
2309 
2310 M*/
2311 
2312 /*MC
2313     MAX_VALUES - Puts the maximum of the scattered/gathered value and the current value into each location
2314 
2315     Level: beginner
2316 
2317 .seealso: InsertMode, VecScatterBegin(), VecScatterEnd(), ADD_VALUES, INSERT_VALUES
2318 
2319 M*/
2320 
2321 PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);
2322 
2323 typedef enum {PETSC_SUBCOMM_GENERAL=0,PETSC_SUBCOMM_CONTIGUOUS=1,PETSC_SUBCOMM_INTERLACED=2} PetscSubcommType;
2324 PETSC_EXTERN const char *const PetscSubcommTypes[];
2325 
2326 /*S
2327    PetscSubcomm - Context of MPI subcommunicators, used by PCREDUNDANT
2328 
2329    Level: advanced
2330 
2331    Concepts: communicator, create
2332 S*/
2333 typedef struct _n_PetscSubcomm* PetscSubcomm;
2334 
2335 struct _n_PetscSubcomm {
2336   MPI_Comm    parent;           /* parent communicator */
2337   MPI_Comm    dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2338   MPI_Comm    comm;             /* this communicator */
2339   PetscMPIInt n;                /* num of subcommunicators under the parent communicator */
2340   PetscMPIInt color;            /* color of processors belong to this communicator */
2341   PetscMPIInt *subsize;         /* size of subcommunicator[color] */
2342   PetscSubcommType type;
2343 };
2344 
2345 PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2346 PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2347 PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2348 PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2349 PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2350 PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2351 PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2352 
2353 /*S
2354    PetscSegBuffer - a segmented extendable buffer
2355 
2356    Level: developer
2357 
2358 .seealso: PetscSegBufferCreate(), PetscSegBufferGet(), PetscSegBufferExtract(), PetscSegBufferDestroy()
2359 S*/
2360 typedef struct _n_PetscSegBuffer *PetscSegBuffer;
2361 PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2362 PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2363 PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2364 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2365 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2366 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2367 PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2368 PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);
2369 
2370 /* Type-safe wrapper to encourage use of PETSC_RESTRICT. Does not use PetscFunctionBegin because the error handling
2371  * prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as
2372  * possible. */
2373 PETSC_STATIC_INLINE PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg,PetscInt count,PetscInt *PETSC_RESTRICT *slot) {return PetscSegBufferGet(seg,count,(void**)slot);}
2374 
2375 extern PetscSegBuffer PetscCitationsList;
2376 #undef __FUNCT__
2377 #define __FUNCT__ "PetscCitationsRegister"
2378 /*@C
2379       PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code.
2380 
2381      Not Collective - only what is registered on rank 0 of PETSC_COMM_WORLD will be printed
2382 
2383      Input Parameters:
2384 +      cite - the bibtex item, formated to displayed on multiple lines nicely
2385 -      set - a boolean variable initially set to PETSC_FALSE; this is used to insure only a single registration of the citation
2386 
2387    Level: intermediate
2388 
2389      Options Database:
2390 .     -citations [filenmae]   - print out the bibtex entries for the given computation
2391 @*/
2392 PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2393 {
2394   size_t         len;
2395   char           *vstring;
2396   PetscErrorCode ierr;
2397 
2398   PetscFunctionBegin;
2399   if (set && *set) PetscFunctionReturn(0);
2400   ierr = PetscStrlen(cit,&len);CHKERRQ(ierr);
2401   ierr = PetscSegBufferGet(PetscCitationsList,(PetscInt)len,&vstring);CHKERRQ(ierr);
2402   ierr = PetscMemcpy(vstring,cit,len);CHKERRQ(ierr);
2403   if (set) *set = PETSC_TRUE;
2404   PetscFunctionReturn(0);
2405 }
2406 
2407 
2408 /* Reset __FUNCT__ in case the user does not define it themselves */
2409 #undef __FUNCT__
2410 #define __FUNCT__ "User provided function"
2411 
2412 #endif
2413