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