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