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