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