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