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