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