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