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