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