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