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