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