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