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