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