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