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