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