xref: /petsc/include/petscsys.h (revision d6cbdb99fb73f50d903551d9fc83deb7cc6e3d53)
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 
1137 
1138 /*E
1139     PetscDataType - Used for handling different basic data types.
1140 
1141    Level: beginner
1142 
1143    Developer comment: It would be nice if we could always just use MPI Datatypes, why can we not?
1144 
1145 .seealso: PetscBinaryRead(), PetscBinaryWrite(), PetscDataTypeToMPIDataType(),
1146           PetscDataTypeGetSize()
1147 
1148 E*/
1149 typedef enum {PETSC_INT = 0,PETSC_DOUBLE = 1,PETSC_COMPLEX = 2, PETSC_LONG = 3 ,PETSC_SHORT = 4,PETSC_FLOAT = 5,
1150               PETSC_CHAR = 6,PETSC_BIT_LOGICAL = 7,PETSC_ENUM = 8,PETSC_BOOL=9, PETSC___FLOAT128 = 10} PetscDataType;
1151 extern const char *PetscDataTypes[];
1152 
1153 #if defined(PETSC_USE_COMPLEX)
1154 #define  PETSC_SCALAR  PETSC_COMPLEX
1155 #else
1156 #if defined(PETSC_USE_REAL_SINGLE)
1157 #define  PETSC_SCALAR  PETSC_FLOAT
1158 #elif defined(PETSC_USE_REAL___FLOAT128)
1159 #define  PETSC_SCALAR  PETSC___FLOAT128
1160 #else
1161 #define  PETSC_SCALAR  PETSC_DOUBLE
1162 #endif
1163 #endif
1164 #if defined(PETSC_USE_REAL_SINGLE)
1165 #define  PETSC_REAL  PETSC_FLOAT
1166 #elif defined(PETSC_USE_REAL___FLOAT128)
1167 #define  PETSC_REAL  PETSC___FLOAT128
1168 #else
1169 #define  PETSC_REAL  PETSC_DOUBLE
1170 #endif
1171 #define  PETSC_FORTRANADDR  PETSC_LONG
1172 
1173 extern PetscErrorCode  PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1174 extern PetscErrorCode  PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1175 extern PetscErrorCode  PetscDataTypeGetSize(PetscDataType,size_t*);
1176 
1177 /*
1178     Basic memory and string operations. These are usually simple wrappers
1179    around the basic Unix system calls, but a few of them have additional
1180    functionality and/or error checking.
1181 */
1182 extern PetscErrorCode    PetscBitMemcpy(void*,PetscInt,const void*,PetscInt,PetscInt,PetscDataType);
1183 extern PetscErrorCode    PetscMemmove(void*,void *,size_t);
1184 extern PetscErrorCode    PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1185 extern PetscErrorCode    PetscStrlen(const char[],size_t*);
1186 extern PetscErrorCode    PetscStrToArray(const char[],int*,char ***);
1187 extern PetscErrorCode    PetscStrToArrayDestroy(int,char **);
1188 extern PetscErrorCode    PetscStrcmp(const char[],const char[],PetscBool  *);
1189 extern PetscErrorCode    PetscStrgrt(const char[],const char[],PetscBool  *);
1190 extern PetscErrorCode    PetscStrcasecmp(const char[],const char[],PetscBool *);
1191 extern PetscErrorCode    PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1192 extern PetscErrorCode    PetscStrcpy(char[],const char[]);
1193 extern PetscErrorCode    PetscStrcat(char[],const char[]);
1194 extern PetscErrorCode    PetscStrncat(char[],const char[],size_t);
1195 extern PetscErrorCode    PetscStrncpy(char[],const char[],size_t);
1196 extern PetscErrorCode    PetscStrchr(const char[],char,char *[]);
1197 extern PetscErrorCode    PetscStrtolower(char[]);
1198 extern PetscErrorCode    PetscStrrchr(const char[],char,char *[]);
1199 extern PetscErrorCode    PetscStrstr(const char[],const char[],char *[]);
1200 extern PetscErrorCode    PetscStrrstr(const char[],const char[],char *[]);
1201 extern PetscErrorCode    PetscStrendswith(const char[],const char[],PetscBool*);
1202 extern PetscErrorCode    PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1203 extern PetscErrorCode    PetscStrallocpy(const char[],char *[]);
1204 extern PetscErrorCode    PetscStrreplace(MPI_Comm,const char[],char[],size_t);
1205 
1206 /*S
1207     PetscToken - 'Token' used for managing tokenizing strings
1208 
1209   Level: intermediate
1210 
1211 .seealso: PetscTokenCreate(), PetscTokenFind(), PetscTokenDestroy()
1212 S*/
1213 typedef struct _p_PetscToken* PetscToken;
1214 
1215 extern PetscErrorCode    PetscTokenCreate(const char[],const char,PetscToken*);
1216 extern PetscErrorCode    PetscTokenFind(PetscToken,char *[]);
1217 extern PetscErrorCode    PetscTokenDestroy(PetscToken*);
1218 
1219 /*
1220    These are  MPI operations for MPI_Allreduce() etc
1221 */
1222 extern  MPI_Op PetscMaxSum_Op;
1223 #if (defined(PETSC_USE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128)
1224 extern  MPI_Op MPIU_SUM;
1225 #else
1226 #define MPIU_SUM MPI_SUM
1227 #endif
1228 #if defined(PETSC_USE_REAL___FLOAT128)
1229 extern  MPI_Op MPIU_MAX;
1230 extern  MPI_Op MPIU_MIN;
1231 #else
1232 #define MPIU_MAX MPI_MAX
1233 #define MPIU_MIN MPI_MIN
1234 #endif
1235 extern PetscErrorCode  PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);
1236 
1237 extern PetscErrorCode MPILong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1238 extern PetscErrorCode MPILong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1239 
1240 /*S
1241      PetscObject - any PETSc object, PetscViewer, Mat, Vec, KSP etc
1242 
1243    Level: beginner
1244 
1245    Note: This is the base class from which all objects appear.
1246 
1247 .seealso:  PetscObjectDestroy(), PetscObjectView(), PetscObjectGetName(), PetscObjectSetName(), PetscObjectReference(), PetscObjectDereferenc()
1248 S*/
1249 typedef struct _p_PetscObject* PetscObject;
1250 
1251 /*S
1252      PetscFList - Linked list of functions, possibly stored in dynamic libraries, accessed
1253       by string name
1254 
1255    Level: advanced
1256 
1257 .seealso:  PetscFListAdd(), PetscFListDestroy()
1258 S*/
1259 typedef struct _n_PetscFList *PetscFList;
1260 
1261 /*E
1262   PetscFileMode - Access mode for a file.
1263 
1264   Level: beginner
1265 
1266   FILE_MODE_READ - open a file at its beginning for reading
1267 
1268   FILE_MODE_WRITE - open a file at its beginning for writing (will create if the file does not exist)
1269 
1270   FILE_MODE_APPEND - open a file at end for writing
1271 
1272   FILE_MODE_UPDATE - open a file for updating, meaning for reading and writing
1273 
1274   FILE_MODE_APPEND_UPDATE - open a file for updating, meaning for reading and writing, at the end
1275 
1276 .seealso: PetscViewerFileSetMode()
1277 E*/
1278 typedef enum {FILE_MODE_READ, FILE_MODE_WRITE, FILE_MODE_APPEND, FILE_MODE_UPDATE, FILE_MODE_APPEND_UPDATE} PetscFileMode;
1279 
1280 #include "petscviewer.h"
1281 #include "petscoptions.h"
1282 
1283 #define PETSC_SMALLEST_CLASSID  1211211
1284 extern PetscClassId   PETSC_LARGEST_CLASSID;
1285 extern PetscClassId   PETSC_OBJECT_CLASSID;
1286 extern PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);
1287 
1288 /*
1289    Routines that get memory usage information from the OS
1290 */
1291 extern PetscErrorCode  PetscMemoryGetCurrentUsage(PetscLogDouble *);
1292 extern PetscErrorCode  PetscMemoryGetMaximumUsage(PetscLogDouble *);
1293 extern PetscErrorCode  PetscMemorySetGetMaximumUsage(void);
1294 extern PetscErrorCode  PetscMemoryShowUsage(PetscViewer,const char[]);
1295 
1296 extern PetscErrorCode  PetscInfoAllow(PetscBool ,const char []);
1297 extern PetscErrorCode  PetscGetTime(PetscLogDouble*);
1298 extern PetscErrorCode  PetscGetCPUTime(PetscLogDouble*);
1299 extern PetscErrorCode  PetscSleep(PetscReal);
1300 
1301 /*
1302    Initialization of PETSc
1303 */
1304 extern PetscErrorCode  PetscInitialize(int*,char***,const char[],const char[]);
1305 PetscPolymorphicSubroutine(PetscInitialize,(int *argc,char ***args),(argc,args,PETSC_NULL,PETSC_NULL))
1306 extern PetscErrorCode  PetscInitializeNoArguments(void);
1307 extern PetscErrorCode  PetscInitialized(PetscBool  *);
1308 extern PetscErrorCode  PetscFinalized(PetscBool  *);
1309 extern PetscErrorCode  PetscFinalize(void);
1310 extern PetscErrorCode  PetscInitializeFortran(void);
1311 extern PetscErrorCode  PetscGetArgs(int*,char ***);
1312 extern PetscErrorCode  PetscGetArguments(char ***);
1313 extern PetscErrorCode  PetscFreeArguments(char **);
1314 
1315 extern PetscErrorCode  PetscEnd(void);
1316 extern PetscErrorCode  PetscSysInitializePackage(const char[]);
1317 
1318 extern MPI_Comm PETSC_COMM_LOCAL_WORLD;
1319 extern PetscErrorCode  PetscHMPIMerge(PetscMPIInt,PetscErrorCode (*)(void*),void*);
1320 extern PetscErrorCode  PetscHMPISpawn(PetscMPIInt);
1321 extern PetscErrorCode  PetscHMPIFinalize(void);
1322 extern PetscErrorCode  PetscHMPIRun(MPI_Comm,PetscErrorCode (*)(MPI_Comm,void *),void*);
1323 extern PetscErrorCode  PetscHMPIRunCtx(MPI_Comm,PetscErrorCode (*)(MPI_Comm,void*,void *),void*);
1324 extern PetscErrorCode  PetscHMPIFree(MPI_Comm,void*);
1325 extern PetscErrorCode  PetscHMPIMalloc(MPI_Comm,size_t,void**);
1326 
1327 extern PetscErrorCode  PetscPythonInitialize(const char[],const char[]);
1328 extern PetscErrorCode  PetscPythonFinalize(void);
1329 extern PetscErrorCode  PetscPythonPrintError(void);
1330 extern PetscErrorCode  PetscPythonMonitorSet(PetscObject,const char[]);
1331 
1332 /*
1333      These are so that in extern C code we can caste function pointers to non-extern C
1334    function pointers. Since the regular C++ code expects its function pointers to be
1335    C++.
1336 */
1337 typedef void (**PetscVoidStarFunction)(void);
1338 typedef void (*PetscVoidFunction)(void);
1339 typedef PetscErrorCode (*PetscErrorCodeFunction)(void);
1340 
1341 /*
1342    PetscTryMethod - Queries an object for a method, if it exists then calls it.
1343               These are intended to be used only inside PETSc functions.
1344 
1345    Level: developer
1346 
1347 .seealso: PetscUseMethod()
1348 */
1349 #define  PetscTryMethod(obj,A,B,C) \
1350   0;{ PetscErrorCode (*f)B, __ierr; \
1351     __ierr = PetscObjectQueryFunction((PetscObject)obj,A,(PetscVoidStarFunction)&f);CHKERRQ(__ierr); \
1352     if (f) {__ierr = (*f)C;CHKERRQ(__ierr);}\
1353   }
1354 
1355 /*
1356    PetscUseMethod - Queries an object for a method, if it exists then calls it, otherwise generates an error.
1357               These are intended to be used only inside PETSc functions.
1358 
1359    Level: developer
1360 
1361 .seealso: PetscTryMethod()
1362 */
1363 #define  PetscUseMethod(obj,A,B,C) \
1364   0;{ PetscErrorCode (*f)B, __ierr; \
1365     __ierr = PetscObjectQueryFunction((PetscObject)obj,A,(PetscVoidStarFunction)&f);CHKERRQ(__ierr); \
1366     if (f) {__ierr = (*f)C;CHKERRQ(__ierr);}\
1367     else SETERRQ1(((PetscObject)obj)->comm,PETSC_ERR_SUP,"Cannot locate function %s in object",A); \
1368   }
1369 
1370 /*
1371     Functions that can act on any PETSc object.
1372 */
1373 extern PetscErrorCode  PetscObjectDestroy(PetscObject*);
1374 extern PetscErrorCode  PetscObjectGetComm(PetscObject,MPI_Comm *);
1375 extern PetscErrorCode  PetscObjectGetClassId(PetscObject,PetscClassId *);
1376 extern PetscErrorCode  PetscObjectGetClassName(PetscObject,const char *[]);
1377 extern PetscErrorCode  PetscObjectSetType(PetscObject,const char []);
1378 extern PetscErrorCode  PetscObjectSetPrecision(PetscObject,PetscPrecision);
1379 extern PetscErrorCode  PetscObjectGetType(PetscObject,const char *[]);
1380 extern PetscErrorCode  PetscObjectSetName(PetscObject,const char[]);
1381 extern PetscErrorCode  PetscObjectGetName(PetscObject,const char*[]);
1382 extern PetscErrorCode  PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer,const char[]);
1383 extern PetscErrorCode  PetscObjectSetTabLevel(PetscObject,PetscInt);
1384 extern PetscErrorCode  PetscObjectGetTabLevel(PetscObject,PetscInt*);
1385 extern PetscErrorCode  PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1386 extern PetscErrorCode  PetscObjectReference(PetscObject);
1387 extern PetscErrorCode  PetscObjectGetReference(PetscObject,PetscInt*);
1388 extern PetscErrorCode  PetscObjectDereference(PetscObject);
1389 extern PetscErrorCode  PetscObjectGetNewTag(PetscObject,PetscMPIInt *);
1390 extern PetscErrorCode  PetscObjectView(PetscObject,PetscViewer);
1391 extern PetscErrorCode  PetscObjectCompose(PetscObject,const char[],PetscObject);
1392 extern PetscErrorCode  PetscObjectRemoveReference(PetscObject,const char[]);
1393 extern PetscErrorCode  PetscObjectQuery(PetscObject,const char[],PetscObject *);
1394 extern PetscErrorCode  PetscObjectComposeFunction(PetscObject,const char[],const char[],void (*)(void));
1395 extern PetscErrorCode  PetscObjectSetFromOptions(PetscObject);
1396 extern PetscErrorCode  PetscObjectSetUp(PetscObject);
1397 extern PetscErrorCode  PetscCommGetNewTag(MPI_Comm,PetscMPIInt *);
1398 extern PetscErrorCode  PetscObjectAddOptionsHandler(PetscObject,PetscErrorCode (*)(PetscObject,void*),PetscErrorCode (*)(PetscObject,void*),void*);
1399 extern PetscErrorCode  PetscObjectProcessOptionsHandlers(PetscObject);
1400 extern PetscErrorCode  PetscObjectDestroyOptionsHandlers(PetscObject);
1401 
1402 /*MC
1403    PetscObjectComposeFunctionDynamic - Associates a function with a given PETSc object.
1404 
1405     Synopsis:
1406     PetscErrorCode PetscObjectComposeFunctionDynamic(PetscObject obj,const char name[],const char fname[],void *ptr)
1407 
1408    Logically Collective on PetscObject
1409 
1410    Input Parameters:
1411 +  obj - the PETSc object; this must be cast with a (PetscObject), for example,
1412          PetscObjectCompose((PetscObject)mat,...);
1413 .  name - name associated with the child function
1414 .  fname - name of the function
1415 -  ptr - function pointer (or PETSC_NULL if using dynamic libraries)
1416 
1417    Level: advanced
1418 
1419 
1420    Notes:
1421    To remove a registered routine, pass in a PETSC_NULL rname and fnc().
1422 
1423    PetscObjectComposeFunctionDynamic() can be used with any PETSc object (such as
1424    Mat, Vec, KSP, SNES, etc.) or any user-provided object.
1425 
1426    The composed function must be wrapped in a EXTERN_C_BEGIN/END for this to
1427    work in C++/complex with dynamic link libraries (./configure options --with-shared-libraries --with-dynamic-loading)
1428    enabled.
1429 
1430    Concepts: objects^composing functions
1431    Concepts: composing functions
1432    Concepts: functions^querying
1433    Concepts: objects^querying
1434    Concepts: querying objects
1435 
1436 .seealso: PetscObjectQueryFunction()
1437 M*/
1438 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1439 #define PetscObjectComposeFunctionDynamic(a,b,c,d) PetscObjectComposeFunction(a,b,c,0)
1440 #else
1441 #define PetscObjectComposeFunctionDynamic(a,b,c,d) PetscObjectComposeFunction(a,b,c,(PetscVoidFunction)(d))
1442 #endif
1443 
1444 extern PetscErrorCode  PetscObjectQueryFunction(PetscObject,const char[],void (**)(void));
1445 extern PetscErrorCode  PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1446 extern PetscErrorCode  PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1447 extern PetscErrorCode  PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1448 extern PetscErrorCode  PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1449 extern PetscErrorCode  PetscObjectAMSPublish(PetscObject);
1450 extern PetscErrorCode  PetscObjectUnPublish(PetscObject);
1451 extern PetscErrorCode  PetscObjectChangeTypeName(PetscObject,const char[]);
1452 extern PetscErrorCode  PetscObjectRegisterDestroy(PetscObject);
1453 extern PetscErrorCode  PetscObjectRegisterDestroyAll(void);
1454 extern PetscErrorCode  PetscObjectName(PetscObject);
1455 extern PetscErrorCode  PetscTypeCompare(PetscObject,const char[],PetscBool *);
1456 extern PetscErrorCode  PetscTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1457 extern PetscErrorCode  PetscRegisterFinalize(PetscErrorCode (*)(void));
1458 extern PetscErrorCode  PetscRegisterFinalizeAll(void);
1459 
1460 /*
1461     Defines PETSc error handling.
1462 */
1463 #include "petscerror.h"
1464 
1465 /*S
1466      PetscOList - Linked list of PETSc objects, each accessable by string name
1467 
1468    Level: developer
1469 
1470    Notes: Used by PetscObjectCompose() and PetscObjectQuery()
1471 
1472 .seealso:  PetscOListAdd(), PetscOListDestroy(), PetscOListFind(), PetscObjectCompose(), PetscObjectQuery()
1473 S*/
1474 typedef struct _n_PetscOList *PetscOList;
1475 
1476 extern PetscErrorCode  PetscOListDestroy(PetscOList*);
1477 extern PetscErrorCode  PetscOListFind(PetscOList,const char[],PetscObject*);
1478 extern PetscErrorCode  PetscOListReverseFind(PetscOList,PetscObject,char**,PetscBool*);
1479 extern PetscErrorCode  PetscOListAdd(PetscOList *,const char[],PetscObject);
1480 extern PetscErrorCode  PetscOListRemoveReference(PetscOList *,const char[]);
1481 extern PetscErrorCode  PetscOListDuplicate(PetscOList,PetscOList *);
1482 
1483 /*
1484     Dynamic library lists. Lists of names of routines in objects or in dynamic
1485   link libraries that will be loaded as needed.
1486 */
1487 extern PetscErrorCode  PetscFListAdd(PetscFList*,const char[],const char[],void (*)(void));
1488 extern PetscErrorCode  PetscFListDestroy(PetscFList*);
1489 extern PetscErrorCode  PetscFListFind(PetscFList,MPI_Comm,const char[],PetscBool,void (**)(void));
1490 extern PetscErrorCode  PetscFListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFList,const char[]);
1491 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1492 #define    PetscFListAddDynamic(a,b,p,c) PetscFListAdd(a,b,p,0)
1493 #else
1494 #define    PetscFListAddDynamic(a,b,p,c) PetscFListAdd(a,b,p,(void (*)(void))c)
1495 #endif
1496 extern PetscErrorCode  PetscFListDuplicate(PetscFList,PetscFList *);
1497 extern PetscErrorCode  PetscFListView(PetscFList,PetscViewer);
1498 extern PetscErrorCode  PetscFListConcat(const char [],const char [],char []);
1499 extern PetscErrorCode  PetscFListGet(PetscFList,char ***,int*);
1500 
1501 /*S
1502      PetscDLLibrary - Linked list of dynamics libraries to search for functions
1503 
1504    Level: advanced
1505 
1506    --with-shared-libraries --with-dynamic-loading must be used with ./configure to use dynamic libraries
1507 
1508 .seealso:  PetscDLLibraryOpen()
1509 S*/
1510 typedef struct _n_PetscDLLibrary *PetscDLLibrary;
1511 extern  PetscDLLibrary DLLibrariesLoaded;
1512 extern PetscErrorCode  PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1513 extern PetscErrorCode  PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1514 extern PetscErrorCode  PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1515 extern PetscErrorCode  PetscDLLibraryPrintPath(PetscDLLibrary);
1516 extern PetscErrorCode  PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1517 extern PetscErrorCode  PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1518 extern PetscErrorCode  PetscDLLibraryClose(PetscDLLibrary);
1519 extern PetscErrorCode  PetscDLLibraryCCAAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1520 
1521 /*
1522   PetscFwk support.  Needs to be documented.
1523   Logically it is an extension of PetscDLLXXX, PetscObjectCompose, etc.
1524 */
1525 #include "petscfwk.h"
1526 
1527 /*
1528      Useful utility routines
1529 */
1530 extern PetscErrorCode  PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1531 extern PetscErrorCode  PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1532 extern PetscErrorCode  PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1533 PetscPolymorphicSubroutine(PetscSequentialPhaseBegin,(MPI_Comm comm),(comm,1))
1534 PetscPolymorphicSubroutine(PetscSequentialPhaseBegin,(void),(PETSC_COMM_WORLD,1))
1535 extern PetscErrorCode  PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1536 PetscPolymorphicSubroutine(PetscSequentialPhaseEnd,(MPI_Comm comm),(comm,1))
1537 PetscPolymorphicSubroutine(PetscSequentialPhaseEnd,(void),(PETSC_COMM_WORLD,1))
1538 extern PetscErrorCode  PetscBarrier(PetscObject);
1539 extern PetscErrorCode  PetscMPIDump(FILE*);
1540 
1541 /*
1542     PetscNot - negates a logical type value and returns result as a PetscBool
1543 
1544     Notes: This is useful in cases like
1545 $     int        *a;
1546 $     PetscBool  flag = PetscNot(a)
1547      where !a does not return a PetscBool  because we cannot provide a cast from int to PetscBool  in C.
1548 */
1549  #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)
1550 
1551 /*
1552     Defines basic graphics available from PETSc.
1553 */
1554 #include "petscdraw.h"
1555 
1556 /*
1557     Defines the base data structures for all PETSc objects
1558 */
1559 #include "private/petscimpl.h"
1560 
1561 
1562 /*MC
1563     PetscErrorPrintf - Prints error messages.
1564 
1565    Synopsis:
1566      PetscErrorCode (*PetscErrorPrintf)(const char format[],...);
1567 
1568     Not Collective
1569 
1570     Input Parameters:
1571 .   format - the usual printf() format string
1572 
1573    Options Database Keys:
1574 +    -error_output_stdout - cause error messages to be printed to stdout instead of the
1575          (default) stderr
1576 -    -error_output_none to turn off all printing of error messages (does not change the way the
1577           error is handled.)
1578 
1579    Notes: Use
1580 $     PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the
1581 $                        error is handled.) and
1582 $     PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on
1583 $        of you can use your own function
1584 
1585           Use
1586      PETSC_STDERR = FILE* obtained from a file open etc. to have stderr printed to the file.
1587      PETSC_STDOUT = FILE* obtained from a file open etc. to have stdout printed to the file.
1588 
1589           Use
1590       PetscPushErrorHandler() to provide your own error handler that determines what kind of messages to print
1591 
1592    Level: developer
1593 
1594     Fortran Note:
1595     This routine is not supported in Fortran.
1596 
1597     Concepts: error messages^printing
1598     Concepts: printing^error messages
1599 
1600 .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscHelpPrintf(), PetscPrintf(), PetscErrorHandlerPush(), PetscVFPrintf(), PetscHelpPrintf()
1601 M*/
1602 extern  PetscErrorCode (*PetscErrorPrintf)(const char[],...);
1603 
1604 /*MC
1605     PetscHelpPrintf - Prints help messages.
1606 
1607    Synopsis:
1608      PetscErrorCode (*PetscHelpPrintf)(const char format[],...);
1609 
1610     Not Collective
1611 
1612     Input Parameters:
1613 .   format - the usual printf() format string
1614 
1615    Level: developer
1616 
1617     Fortran Note:
1618     This routine is not supported in Fortran.
1619 
1620     Concepts: help messages^printing
1621     Concepts: printing^help messages
1622 
1623 .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1624 M*/
1625 extern  PetscErrorCode  (*PetscHelpPrintf)(MPI_Comm,const char[],...);
1626 
1627 /*
1628      Defines PETSc profiling.
1629 */
1630 #include "petsclog.h"
1631 
1632 /*
1633           For locking, unlocking and destroying AMS memories associated with  PETSc objects. ams.h is included in petscviewer.h
1634 */
1635 #if defined(PETSC_HAVE_AMS)
1636 extern PetscBool  PetscAMSPublishAll;
1637 #define PetscObjectTakeAccess(obj)  ((((PetscObject)(obj))->amem == -1) ? 0 : AMS_Memory_take_access(((PetscObject)(obj))->amem))
1638 #define PetscObjectGrantAccess(obj) ((((PetscObject)(obj))->amem == -1) ? 0 : AMS_Memory_grant_access(((PetscObject)(obj))->amem))
1639 #define PetscObjectDepublish(obj)   ((((PetscObject)(obj))->amem == -1) ? 0 : AMS_Memory_destroy(((PetscObject)(obj))->amem));((PetscObject)(obj))->amem = -1;
1640 #else
1641 #define PetscObjectTakeAccess(obj)   0
1642 #define PetscObjectGrantAccess(obj)  0
1643 #define PetscObjectDepublish(obj)      0
1644 #endif
1645 
1646 /*
1647       Simple PETSc parallel IO for ASCII printing
1648 */
1649 extern PetscErrorCode   PetscFixFilename(const char[],char[]);
1650 extern PetscErrorCode   PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1651 extern PetscErrorCode   PetscFClose(MPI_Comm,FILE*);
1652 extern PetscErrorCode   PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1653 extern PetscErrorCode   PetscPrintf(MPI_Comm,const char[],...);
1654 extern PetscErrorCode   PetscSNPrintf(char*,size_t,const char [],...);
1655 extern PetscErrorCode   PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);
1656 
1657 
1658 
1659 /* These are used internally by PETSc ASCII IO routines*/
1660 #include <stdarg.h>
1661 extern PetscErrorCode   PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list);
1662 extern PetscErrorCode   (*PetscVFPrintf)(FILE*,const char[],va_list);
1663 extern PetscErrorCode   PetscVFPrintfDefault(FILE*,const char[],va_list);
1664 
1665 #if defined(PETSC_HAVE_MATLAB_ENGINE)
1666 extern PetscErrorCode  PetscVFPrintf_Matlab(FILE*,const char[],va_list);
1667 #endif
1668 
1669 extern PetscErrorCode  PetscErrorPrintfDefault(const char [],...);
1670 extern PetscErrorCode  PetscErrorPrintfNone(const char [],...);
1671 extern PetscErrorCode  PetscHelpPrintfDefault(MPI_Comm,const char [],...);
1672 
1673 #if defined(PETSC_HAVE_POPEN)
1674 extern PetscErrorCode   PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1675 extern PetscErrorCode   PetscPClose(MPI_Comm,FILE*);
1676 #endif
1677 
1678 extern PetscErrorCode   PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1679 extern PetscErrorCode   PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1680 extern PetscErrorCode   PetscSynchronizedFlush(MPI_Comm);
1681 extern PetscErrorCode   PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1682 extern PetscErrorCode   PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1683 extern PetscErrorCode   PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1684 extern PetscErrorCode   PetscGetPetscDir(const char*[]);
1685 
1686 extern PetscErrorCode   PetscPopUpSelect(MPI_Comm,const char*,const char*,int,const char**,int*);
1687 
1688 /*S
1689      PetscContainer - Simple PETSc object that contains a pointer to any required data
1690 
1691    Level: advanced
1692 
1693 .seealso:  PetscObject, PetscContainerCreate()
1694 S*/
1695 extern PetscClassId  PETSC_CONTAINER_CLASSID;
1696 typedef struct _p_PetscContainer*  PetscContainer;
1697 extern PetscErrorCode  PetscContainerGetPointer(PetscContainer,void **);
1698 extern PetscErrorCode  PetscContainerSetPointer(PetscContainer,void *);
1699 extern PetscErrorCode  PetscContainerDestroy(PetscContainer*);
1700 extern PetscErrorCode  PetscContainerCreate(MPI_Comm,PetscContainer *);
1701 extern PetscErrorCode  PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void*));
1702 
1703 /*
1704    For use in debuggers
1705 */
1706 extern  PetscMPIInt PetscGlobalRank;
1707 extern  PetscMPIInt PetscGlobalSize;
1708 extern PetscErrorCode  PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1709 extern PetscErrorCode  PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1710 extern PetscErrorCode  PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);
1711 
1712 #if defined(PETSC_HAVE_MEMORY_H)
1713 #include <memory.h>
1714 #endif
1715 #if defined(PETSC_HAVE_STDLIB_H)
1716 #include <stdlib.h>
1717 #endif
1718 #if defined(PETSC_HAVE_STRINGS_H)
1719 #include <strings.h>
1720 #endif
1721 #if defined(PETSC_HAVE_STRING_H)
1722 #include <string.h>
1723 #endif
1724 
1725 
1726 #if defined(PETSC_HAVE_XMMINTRIN_H)
1727 #include <xmmintrin.h>
1728 #endif
1729 #if defined(PETSC_HAVE_STDINT_H)
1730 #include <stdint.h>
1731 #endif
1732 
1733 #undef __FUNCT__
1734 #define __FUNCT__ "PetscMemcpy"
1735 /*@C
1736    PetscMemcpy - Copies n bytes, beginning at location b, to the space
1737    beginning at location a. The two memory regions CANNOT overlap, use
1738    PetscMemmove() in that case.
1739 
1740    Not Collective
1741 
1742    Input Parameters:
1743 +  b - pointer to initial memory space
1744 -  n - length (in bytes) of space to copy
1745 
1746    Output Parameter:
1747 .  a - pointer to copy space
1748 
1749    Level: intermediate
1750 
1751    Compile Option:
1752     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1753                                   for memory copies on double precision values.
1754     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1755                                   for memory copies on double precision values.
1756     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1757                                   for memory copies on double precision values.
1758 
1759    Note:
1760    This routine is analogous to memcpy().
1761 
1762    Developer Note: this is inlined for fastest performance
1763 
1764   Concepts: memory^copying
1765   Concepts: copying^memory
1766 
1767 .seealso: PetscMemmove()
1768 
1769 @*/
1770 PETSC_STATIC_INLINE PetscErrorCode  PetscMemcpy(void *a,const void *b,size_t n)
1771 {
1772 #if defined(PETSC_USE_DEBUG)
1773   unsigned long al = (unsigned long) a,bl = (unsigned long) b;
1774   unsigned long nl = (unsigned long) n;
1775   PetscFunctionBegin;
1776   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1777   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1778 #else
1779   PetscFunctionBegin;
1780 #endif
1781   if (a != b) {
1782 #if defined(PETSC_USE_DEBUG)
1783     if ((al > bl && (al - bl) < nl) || (bl - al) < nl) {
1784       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1785               or make sure your copy regions and lengths are correct. \n\
1786               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1787     }
1788 #endif
1789 #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1790    if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1791       size_t len = n/sizeof(PetscScalar);
1792 #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1793       PetscBLASInt one = 1,blen = PetscBLASIntCast(len);
1794       BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one);
1795 #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1796       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1797 #else
1798       size_t      i;
1799       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1800       for (i=0; i<len; i++) y[i] = x[i];
1801 #endif
1802     } else {
1803       memcpy((char*)(a),(char*)(b),n);
1804     }
1805 #else
1806     memcpy((char*)(a),(char*)(b),n);
1807 #endif
1808   }
1809   PetscFunctionReturn(0);
1810 }
1811 
1812 /*@C
1813    PetscMemzero - Zeros the specified memory.
1814 
1815    Not Collective
1816 
1817    Input Parameters:
1818 +  a - pointer to beginning memory location
1819 -  n - length (in bytes) of memory to initialize
1820 
1821    Level: intermediate
1822 
1823    Compile Option:
1824    PETSC_PREFER_BZERO - on certain machines (the IBM RS6000) the bzero() routine happens
1825   to be faster than the memset() routine. This flag causes the bzero() routine to be used.
1826 
1827    Developer Note: this is inlined for fastest performance
1828 
1829    Concepts: memory^zeroing
1830    Concepts: zeroing^memory
1831 
1832 .seealso: PetscMemcpy()
1833 @*/
1834 PETSC_STATIC_INLINE PetscErrorCode  PetscMemzero(void *a,size_t n)
1835 {
1836   if (n > 0) {
1837 #if defined(PETSC_USE_DEBUG)
1838     if (!a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer");
1839 #endif
1840 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1841     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1842       size_t      i,len = n/sizeof(PetscScalar);
1843       PetscScalar *x = (PetscScalar*)a;
1844       for (i=0; i<len; i++) x[i] = 0.0;
1845     } else {
1846 #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1847     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1848       PetscInt len = n/sizeof(PetscScalar);
1849       fortranzero_(&len,(PetscScalar*)a);
1850     } else {
1851 #endif
1852 #if defined(PETSC_PREFER_BZERO)
1853       bzero((char *)a,n);
1854 #else
1855       memset((char*)a,0,n);
1856 #endif
1857 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1858     }
1859 #endif
1860   }
1861   return 0;
1862 }
1863 
1864 /*MC
1865    PetscPrefetchBlock - Prefetches a block of memory
1866 
1867    Synopsis:
1868     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)
1869 
1870    Not Collective
1871 
1872    Input Parameters:
1873 +  a - pointer to first element to fetch (any type but usually PetscInt or PetscScalar)
1874 .  n - number of elements to fetch
1875 .  rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors)
1876 -  t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note
1877 
1878    Level: developer
1879 
1880    Notes:
1881    The last two arguments (rw and t) must be compile-time constants.
1882 
1883    Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality.  Not all architectures offer
1884    equivalent locality hints, but the following macros are always defined to their closest analogue.
1885 +  PETSC_PREFETCH_HINT_NTA - Non-temporal.  Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched).
1886 .  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.
1887 .  PETSC_PREFETCH_HINT_T1 - Fetch to level 2 and higher (not L1).
1888 -  PETSC_PREFETCH_HINT_T2 - Fetch to high-level cache only.  (On many systems, T0 and T1 are equivalent.)
1889 
1890    This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid
1891    address).
1892 
1893    Concepts: memory
1894 M*/
1895 #define PetscPrefetchBlock(a,n,rw,t) do {                               \
1896     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
1897     for ( ; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
1898   } while (0)
1899 
1900 /*
1901     Allows accessing MATLAB Engine
1902 */
1903 #include "petscmatlab.h"
1904 
1905 /*
1906       Determine if some of the kernel computation routines use
1907    Fortran (rather than C) for the numerical calculations. On some machines
1908    and compilers (like complex numbers) the Fortran version of the routines
1909    is faster than the C/C++ versions. The flag --with-fortran-kernels
1910    should be used with ./configure to turn these on.
1911 */
1912 #if defined(PETSC_USE_FORTRAN_KERNELS)
1913 
1914 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1915 #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1916 #endif
1917 
1918 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
1919 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
1920 #endif
1921 
1922 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1923 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
1924 #endif
1925 
1926 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1927 #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
1928 #endif
1929 
1930 #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
1931 #define PETSC_USE_FORTRAN_KERNEL_NORM
1932 #endif
1933 
1934 #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
1935 #define PETSC_USE_FORTRAN_KERNEL_MAXPY
1936 #endif
1937 
1938 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1939 #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
1940 #endif
1941 
1942 #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1943 #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
1944 #endif
1945 
1946 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
1947 #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
1948 #endif
1949 
1950 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1951 #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
1952 #endif
1953 
1954 #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
1955 #define PETSC_USE_FORTRAN_KERNEL_MDOT
1956 #endif
1957 
1958 #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
1959 #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
1960 #endif
1961 
1962 #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
1963 #define PETSC_USE_FORTRAN_KERNEL_AYPX
1964 #endif
1965 
1966 #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
1967 #define PETSC_USE_FORTRAN_KERNEL_WAXPY
1968 #endif
1969 
1970 #endif
1971 
1972 /*
1973     Macros for indicating code that should be compiled with a C interface,
1974    rather than a C++ interface. Any routines that are dynamically loaded
1975    (such as the PCCreate_XXX() routines) must be wrapped so that the name
1976    mangler does not change the functions symbol name. This just hides the
1977    ugly extern "C" {} wrappers.
1978 */
1979 #if defined(__cplusplus)
1980 #define EXTERN_C_BEGIN extern "C" {
1981 #define EXTERN_C_END }
1982 #else
1983 #define EXTERN_C_BEGIN
1984 #define EXTERN_C_END
1985 #endif
1986 
1987 /* --------------------------------------------------------------------*/
1988 
1989 /*MC
1990     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
1991         communication
1992 
1993    Level: beginner
1994 
1995    Note: This manual page is a place-holder because MPICH does not have a manual page for MPI_Comm
1996 
1997 .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
1998 M*/
1999 
2000 /*MC
2001     PetscScalar - PETSc type that represents either a double precision real number, a double precision
2002        complex number, a single precision real number, a long double or an int - if the code is configured
2003        with --with-scalar-type=real,complex --with-precision=single,double,longdouble,int,matsingle
2004 
2005 
2006    Level: beginner
2007 
2008 .seealso: PetscReal, PassiveReal, PassiveScalar, MPIU_SCALAR, PetscInt
2009 M*/
2010 
2011 /*MC
2012     PetscReal - PETSc type that represents a real number version of PetscScalar
2013 
2014    Level: beginner
2015 
2016 .seealso: PetscScalar, PassiveReal, PassiveScalar
2017 M*/
2018 
2019 /*MC
2020     PassiveScalar - PETSc type that represents a PetscScalar
2021    Level: beginner
2022 
2023     This is the same as a PetscScalar except in code that is automatically differentiated it is
2024    treated as a constant (not an indendent or dependent variable)
2025 
2026 .seealso: PetscReal, PassiveReal, PetscScalar
2027 M*/
2028 
2029 /*MC
2030     PassiveReal - PETSc type that represents a PetscReal
2031 
2032    Level: beginner
2033 
2034     This is the same as a PetscReal except in code that is automatically differentiated it is
2035    treated as a constant (not an indendent or dependent variable)
2036 
2037 .seealso: PetscScalar, PetscReal, PassiveScalar
2038 M*/
2039 
2040 /*MC
2041     MPIU_SCALAR - MPI datatype corresponding to PetscScalar
2042 
2043    Level: beginner
2044 
2045     Note: In MPI calls that require an MPI datatype that matches a PetscScalar or array of PetscScalars
2046           pass this value
2047 
2048 .seealso: PetscReal, PassiveReal, PassiveScalar, PetscScalar, MPIU_INT
2049 M*/
2050 
2051 #if defined(PETSC_HAVE_MPIIO)
2052 #if !defined(PETSC_WORDS_BIGENDIAN)
2053 extern PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2054 extern PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2055 #else
2056 #define MPIU_File_write_all(a,b,c,d,e) MPI_File_write_all(a,b,c,d,e)
2057 #define MPIU_File_read_all(a,b,c,d,e) MPI_File_read_all(a,b,c,d,e)
2058 #endif
2059 #endif
2060 
2061 /* the following petsc_static_inline require petscerror.h */
2062 
2063 /* Limit MPI to 32-bits */
2064 #define PETSC_MPI_INT_MAX  2147483647
2065 #define PETSC_MPI_INT_MIN -2147483647
2066 /* Limit BLAS to 32-bits */
2067 #define PETSC_BLAS_INT_MAX  2147483647
2068 #define PETSC_BLAS_INT_MIN -2147483647
2069 /* On 32 bit systems HDF5 is limited by size of integer, because hsize_t is defined as size_t */
2070 #define PETSC_HDF5_INT_MAX  2147483647
2071 #define PETSC_HDF5_INT_MIN -2147483647
2072 
2073 #if defined(PETSC_USE_64BIT_INDICES)
2074 #define PetscMPIIntCheck(a)  if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Message too long for MPI")
2075 #define PetscBLASIntCheck(a)  if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK")
2076 #define PetscMPIIntCast(a) (PetscMPIInt)(a);PetscMPIIntCheck(a)
2077 #define PetscBLASIntCast(a) (PetscBLASInt)(a);PetscBLASIntCheck(a)
2078 
2079 #if (PETSC_SIZEOF_SIZE_T == 4)
2080 #define PetscHDF5IntCheck(a)  if ((a) > PETSC_HDF5_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for HDF5")
2081 #define PetscHDF5IntCast(a) (hsize_t)(a);PetscHDF5IntCheck(a)
2082 #else
2083 #define PetscHDF5IntCheck(a)
2084 #define PetscHDF5IntCast(a) a
2085 #endif
2086 
2087 #else
2088 #define PetscMPIIntCheck(a)
2089 #define PetscBLASIntCheck(a)
2090 #define PetscHDF5IntCheck(a)
2091 #define PetscMPIIntCast(a) a
2092 #define PetscBLASIntCast(a) a
2093 #define PetscHDF5IntCast(a) a
2094 #endif
2095 
2096 
2097 /*
2098      The IBM include files define hz, here we hide it so that it may be used
2099    as a regular user variable.
2100 */
2101 #if defined(hz)
2102 #undef hz
2103 #endif
2104 
2105 /*  For arrays that contain filenames or paths */
2106 
2107 
2108 #if defined(PETSC_HAVE_LIMITS_H)
2109 #include <limits.h>
2110 #endif
2111 #if defined(PETSC_HAVE_SYS_PARAM_H)
2112 #include <sys/param.h>
2113 #endif
2114 #if defined(PETSC_HAVE_SYS_TYPES_H)
2115 #include <sys/types.h>
2116 #endif
2117 #if defined(MAXPATHLEN)
2118 #  define PETSC_MAX_PATH_LEN     MAXPATHLEN
2119 #elif defined(MAX_PATH)
2120 #  define PETSC_MAX_PATH_LEN     MAX_PATH
2121 #elif defined(_MAX_PATH)
2122 #  define PETSC_MAX_PATH_LEN     _MAX_PATH
2123 #else
2124 #  define PETSC_MAX_PATH_LEN     4096
2125 #endif
2126 
2127 /* Special support for C++ */
2128 #include "petscsys.hh"
2129 
2130 
2131 /*MC
2132 
2133     UsingFortran - Fortran can be used with PETSc in four distinct approaches
2134 
2135 $    1) classic Fortran 77 style
2136 $#include "finclude/petscXXX.h" to work with material from the XXX component of PETSc
2137 $       XXX variablename
2138 $      You cannot use this approach if you wish to use the Fortran 90 specific PETSc routines
2139 $      which end in F90; such as VecGetArrayF90()
2140 $
2141 $    2) classic Fortran 90 style
2142 $#include "finclude/petscXXX.h"
2143 $#include "finclude/petscXXX.h90" to work with material from the XXX component of PETSc
2144 $       XXX variablename
2145 $
2146 $    3) Using Fortran modules
2147 $#include "finclude/petscXXXdef.h"
2148 $         use petscXXXX
2149 $       XXX variablename
2150 $
2151 $    4) Use Fortran modules and Fortran data types for PETSc types
2152 $#include "finclude/petscXXXdef.h"
2153 $         use petscXXXX
2154 $       type(XXX) variablename
2155 $      To use this approach you must ./configure PETSc with the additional
2156 $      option --with-fortran-datatypes You cannot use the type(XXX) declaration approach without using Fortran modules
2157 
2158     Finally if you absolutely do not want to use any #include you can use either
2159 
2160 $    3a) skip the #include BUT you cannot use any PETSc data type names like Vec, Mat, PetscInt, PetscErrorCode etc
2161 $        and you must declare the variables as integer, for example
2162 $        integer variablename
2163 $
2164 $    4a) skip the #include, you use the object types like type(Vec) type(Mat) but cannot use the data type
2165 $        names like PetscErrorCode, PetscInt etc. again for those you must use integer
2166 
2167    We recommend either 2 or 3. Approaches 2 and 3 provide type checking for most PETSc function calls; 4 has type checking
2168 for only a few PETSc functions.
2169 
2170    Fortran type checking with interfaces is strick, this means you cannot pass a scalar value when an array value
2171 is expected (even though it is legal Fortran). For example when setting a single value in a matrix with MatSetValues()
2172 you cannot have something like
2173 $      PetscInt row,col
2174 $      PetscScalar val
2175 $        ...
2176 $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)
2177 You must instead have
2178 $      PetscInt row(1),col(1)
2179 $      PetscScalar val(1)
2180 $        ...
2181 $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)
2182 
2183 
2184     See the example src/vec/vec/examples/tutorials/ex20f90.F90 for an example that can use all four approaches
2185 
2186     Developer Notes: The finclude/petscXXXdef.h contain all the #defines (would be typedefs in C code) these
2187      automatically include their predecessors; for example finclude/petscvecdef.h includes finclude/petscisdef.h
2188 
2189      The finclude/petscXXXX.h contain all the parameter statements for that package. These automatically include
2190      their finclude/petscXXXdef.h file but DO NOT automatically include their predecessors;  for example
2191      finclude/petscvec.h does NOT automatically include finclude/petscis.h
2192 
2193      The finclude/ftn-custom/petscXXXdef.h90 are not intended to be used directly in code, they define the
2194      Fortran data type type(XXX) (for example type(Vec)) when PETSc is ./configure with the --with-fortran-datatypes option.
2195 
2196      The finclude/ftn-custom/petscXXX.h90 (not included directly by code) contain interface definitions for
2197      the PETSc Fortran stubs that have different bindings then their C version (for example VecGetArrayF90).
2198 
2199      The finclude/ftn-auto/petscXXX.h90 (not included directly by code) contain interface definitions generated
2200      automatically by "make allfortranstubs".
2201 
2202      The finclude/petscXXX.h90 includes the custom finclude/ftn-custom/petscXXX.h90 and if ./configure
2203      was run with --with-fortran-interfaces it also includes the finclude/ftn-auto/petscXXX.h90 These DO NOT automatically
2204      include their predecessors
2205 
2206     Level: beginner
2207 
2208 M*/
2209 
2210 extern PetscErrorCode  PetscGetArchType(char[],size_t);
2211 extern PetscErrorCode  PetscGetHostName(char[],size_t);
2212 extern PetscErrorCode  PetscGetUserName(char[],size_t);
2213 extern PetscErrorCode  PetscGetProgramName(char[],size_t);
2214 extern PetscErrorCode  PetscSetProgramName(const char[]);
2215 extern PetscErrorCode  PetscGetDate(char[],size_t);
2216 
2217 extern PetscErrorCode  PetscSortInt(PetscInt,PetscInt[]);
2218 extern PetscErrorCode  PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2219 extern PetscErrorCode  PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2220 extern PetscErrorCode  PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2221 extern PetscErrorCode  PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2222 extern PetscErrorCode  PetscSortIntWithArrayPair(PetscInt,PetscInt*,PetscInt*,PetscInt*);
2223 extern PetscErrorCode  PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2224 extern PetscErrorCode  PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2225 extern PetscErrorCode  PetscSortReal(PetscInt,PetscReal[]);
2226 extern PetscErrorCode  PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2227 extern PetscErrorCode  PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2228 extern PetscErrorCode  PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2229 extern PetscErrorCode  PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2230 extern PetscErrorCode  PetscMergeIntArrayPair(PetscInt,const PetscInt*,const PetscInt*,PetscInt,const PetscInt*,const PetscInt*,PetscInt*,PetscInt**,PetscInt**);
2231 
2232 extern PetscErrorCode  PetscSetDisplay(void);
2233 extern PetscErrorCode  PetscGetDisplay(char[],size_t);
2234 
2235 /*J
2236     PetscRandomType - String with the name of a PETSc randomizer
2237        with an optional dynamic library name, for example
2238        http://www.mcs.anl.gov/petsc/lib.a:myrandcreate()
2239 
2240    Level: beginner
2241 
2242    Notes: to use the SPRNG you must have ./configure PETSc
2243    with the option --download-sprng
2244 
2245 .seealso: PetscRandomSetType(), PetscRandom
2246 J*/
2247 #define PetscRandomType char*
2248 #define PETSCRAND       "rand"
2249 #define PETSCRAND48     "rand48"
2250 #define PETSCSPRNG      "sprng"
2251 
2252 /* Logging support */
2253 extern  PetscClassId PETSC_RANDOM_CLASSID;
2254 
2255 extern PetscErrorCode  PetscRandomInitializePackage(const char[]);
2256 
2257 /*S
2258      PetscRandom - Abstract PETSc object that manages generating random numbers
2259 
2260    Level: intermediate
2261 
2262   Concepts: random numbers
2263 
2264 .seealso:  PetscRandomCreate(), PetscRandomGetValue(), PetscRandomType
2265 S*/
2266 typedef struct _p_PetscRandom*   PetscRandom;
2267 
2268 /* Dynamic creation and loading functions */
2269 extern PetscFList PetscRandomList;
2270 extern PetscBool  PetscRandomRegisterAllCalled;
2271 
2272 extern PetscErrorCode  PetscRandomRegisterAll(const char []);
2273 extern PetscErrorCode  PetscRandomRegister(const char[],const char[],const char[],PetscErrorCode (*)(PetscRandom));
2274 extern PetscErrorCode  PetscRandomRegisterDestroy(void);
2275 extern PetscErrorCode  PetscRandomSetType(PetscRandom, const PetscRandomType);
2276 extern PetscErrorCode  PetscRandomSetFromOptions(PetscRandom);
2277 extern PetscErrorCode  PetscRandomGetType(PetscRandom, const PetscRandomType*);
2278 extern PetscErrorCode  PetscRandomViewFromOptions(PetscRandom,char*);
2279 extern PetscErrorCode  PetscRandomView(PetscRandom,PetscViewer);
2280 
2281 /*MC
2282   PetscRandomRegisterDynamic - Adds a new PetscRandom component implementation
2283 
2284   Synopsis:
2285   PetscErrorCode PetscRandomRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(PetscRandom))
2286 
2287   Not Collective
2288 
2289   Input Parameters:
2290 + name        - The name of a new user-defined creation routine
2291 . path        - The path (either absolute or relative) of the library containing this routine
2292 . func_name   - The name of routine to create method context
2293 - create_func - The creation routine itself
2294 
2295   Notes:
2296   PetscRandomRegisterDynamic() may be called multiple times to add several user-defined randome number generators
2297 
2298   If dynamic libraries are used, then the fourth input argument (routine_create) is ignored.
2299 
2300   Sample usage:
2301 .vb
2302     PetscRandomRegisterDynamic("my_rand","/home/username/my_lib/lib/libO/solaris/libmy.a", "MyPetscRandomtorCreate", MyPetscRandomtorCreate);
2303 .ve
2304 
2305   Then, your random type can be chosen with the procedural interface via
2306 .vb
2307     PetscRandomCreate(MPI_Comm, PetscRandom *);
2308     PetscRandomSetType(PetscRandom,"my_random_name");
2309 .ve
2310    or at runtime via the option
2311 .vb
2312     -random_type my_random_name
2313 .ve
2314 
2315   Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
2316 
2317          For an example of the code needed to interface your own random number generator see
2318          src/sys/random/impls/rand/rand.c
2319 
2320   Level: advanced
2321 
2322 .keywords: PetscRandom, register
2323 .seealso: PetscRandomRegisterAll(), PetscRandomRegisterDestroy(), PetscRandomRegister()
2324 M*/
2325 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
2326 #define PetscRandomRegisterDynamic(a,b,c,d) PetscRandomRegister(a,b,c,0)
2327 #else
2328 #define PetscRandomRegisterDynamic(a,b,c,d) PetscRandomRegister(a,b,c,d)
2329 #endif
2330 
2331 extern PetscErrorCode  PetscRandomCreate(MPI_Comm,PetscRandom*);
2332 extern PetscErrorCode  PetscRandomGetValue(PetscRandom,PetscScalar*);
2333 extern PetscErrorCode  PetscRandomGetValueReal(PetscRandom,PetscReal*);
2334 extern PetscErrorCode  PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2335 extern PetscErrorCode  PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2336 extern PetscErrorCode  PetscRandomSetSeed(PetscRandom,unsigned long);
2337 extern PetscErrorCode  PetscRandomGetSeed(PetscRandom,unsigned long *);
2338 extern PetscErrorCode  PetscRandomSeed(PetscRandom);
2339 extern PetscErrorCode  PetscRandomDestroy(PetscRandom*);
2340 
2341 extern PetscErrorCode  PetscGetFullPath(const char[],char[],size_t);
2342 extern PetscErrorCode  PetscGetRelativePath(const char[],char[],size_t);
2343 extern PetscErrorCode  PetscGetWorkingDirectory(char[],size_t);
2344 extern PetscErrorCode  PetscGetRealPath(const char[],char[]);
2345 extern PetscErrorCode  PetscGetHomeDirectory(char[],size_t);
2346 extern PetscErrorCode  PetscTestFile(const char[],char,PetscBool *);
2347 extern PetscErrorCode  PetscTestDirectory(const char[],char,PetscBool *);
2348 
2349 extern PetscErrorCode  PetscBinaryRead(int,void*,PetscInt,PetscDataType);
2350 extern PetscErrorCode  PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscDataType);
2351 extern PetscErrorCode  PetscBinarySynchronizedWrite(MPI_Comm,int,void*,PetscInt,PetscDataType,PetscBool );
2352 extern PetscErrorCode  PetscBinaryWrite(int,void*,PetscInt,PetscDataType,PetscBool );
2353 extern PetscErrorCode  PetscBinaryOpen(const char[],PetscFileMode,int *);
2354 extern PetscErrorCode  PetscBinaryClose(int);
2355 extern PetscErrorCode  PetscSharedTmp(MPI_Comm,PetscBool  *);
2356 extern PetscErrorCode  PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2357 extern PetscErrorCode  PetscGetTmp(MPI_Comm,char[],size_t);
2358 extern PetscErrorCode  PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2359 extern PetscErrorCode  PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2360 extern PetscErrorCode  PetscOpenSocket(char*,int,int*);
2361 extern PetscErrorCode  PetscWebServe(MPI_Comm,int);
2362 
2363 /*
2364    In binary files variables are stored using the following lengths,
2365   regardless of how they are stored in memory on any one particular
2366   machine. Use these rather then sizeof() in computing sizes for
2367   PetscBinarySeek().
2368 */
2369 #define PETSC_BINARY_INT_SIZE   (32/8)
2370 #define PETSC_BINARY_FLOAT_SIZE  (32/8)
2371 #define PETSC_BINARY_CHAR_SIZE  (8/8)
2372 #define PETSC_BINARY_SHORT_SIZE  (16/8)
2373 #define PETSC_BINARY_DOUBLE_SIZE  (64/8)
2374 #define PETSC_BINARY_SCALAR_SIZE  sizeof(PetscScalar)
2375 
2376 /*E
2377   PetscBinarySeekType - argument to PetscBinarySeek()
2378 
2379   Level: advanced
2380 
2381 .seealso: PetscBinarySeek(), PetscBinarySynchronizedSeek()
2382 E*/
2383 typedef enum {PETSC_BINARY_SEEK_SET = 0,PETSC_BINARY_SEEK_CUR = 1,PETSC_BINARY_SEEK_END = 2} PetscBinarySeekType;
2384 extern PetscErrorCode  PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2385 extern PetscErrorCode  PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2386 extern PetscErrorCode  PetscByteSwap(void *,PetscDataType,PetscInt);
2387 
2388 extern PetscErrorCode  PetscSetDebugTerminal(const char[]);
2389 extern PetscErrorCode  PetscSetDebugger(const char[],PetscBool );
2390 extern PetscErrorCode  PetscSetDefaultDebugger(void);
2391 extern PetscErrorCode  PetscSetDebuggerFromString(char*);
2392 extern PetscErrorCode  PetscAttachDebugger(void);
2393 extern PetscErrorCode  PetscStopForDebugger(void);
2394 
2395 extern PetscErrorCode  PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2396 extern PetscErrorCode  PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2397 extern PetscErrorCode  PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2398 extern PetscErrorCode  PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2399 extern PetscErrorCode  PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2400 
2401 extern PetscErrorCode  PetscSSEIsEnabled(MPI_Comm,PetscBool  *,PetscBool  *);
2402 
2403 /*E
2404   InsertMode - Whether entries are inserted or added into vectors or matrices
2405 
2406   Level: beginner
2407 
2408 .seealso: VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2409           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(),
2410           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd()
2411 E*/
2412  typedef enum {NOT_SET_VALUES, INSERT_VALUES, ADD_VALUES, MAX_VALUES, INSERT_ALL_VALUES, ADD_ALL_VALUES} InsertMode;
2413 
2414 /*MC
2415     INSERT_VALUES - Put a value into a vector or matrix, overwrites any previous value
2416 
2417     Level: beginner
2418 
2419 .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2420           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), ADD_VALUES,
2421           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES
2422 
2423 M*/
2424 
2425 /*MC
2426     ADD_VALUES - Adds a value into a vector or matrix, if there previously was no value, just puts the
2427                 value into that location
2428 
2429     Level: beginner
2430 
2431 .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2432           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), INSERT_VALUES,
2433           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES
2434 
2435 M*/
2436 
2437 /*MC
2438     MAX_VALUES - Puts the maximum of the scattered/gathered value and the current value into each location
2439 
2440     Level: beginner
2441 
2442 .seealso: InsertMode, VecScatterBegin(), VecScatterEnd(), ADD_VALUES, INSERT_VALUES
2443 
2444 M*/
2445 
2446 /*S
2447    PetscSubcomm - Context of MPI subcommunicators, used by PCREDUNDANT
2448 
2449    Level: advanced
2450 
2451    Concepts: communicator, create
2452 S*/
2453 typedef struct _n_PetscSubcomm* PetscSubcomm;
2454 
2455 struct _n_PetscSubcomm {
2456   MPI_Comm   parent;      /* parent communicator */
2457   MPI_Comm   dupparent;   /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2458   MPI_Comm   comm;        /* this communicator */
2459   PetscInt   n;           /* num of subcommunicators under the parent communicator */
2460   PetscInt   color;       /* color of processors belong to this communicator */
2461 };
2462 
2463 typedef enum {PETSC_SUBCOMM_GENERAL=0,PETSC_SUBCOMM_CONTIGUOUS=1,PETSC_SUBCOMM_INTERLACED=2} PetscSubcommType;
2464 extern const char *PetscSubcommTypes[];
2465 
2466 extern PetscErrorCode  PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2467 extern PetscErrorCode  PetscSubcommDestroy(PetscSubcomm*);
2468 extern PetscErrorCode  PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2469 extern PetscErrorCode  PetscSubcommSetType(PetscSubcomm,const PetscSubcommType);
2470 extern PetscErrorCode  PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt,PetscMPIInt);
2471 
2472 #include <petscctable.h>
2473 
2474 PETSC_EXTERN_CXX_END
2475 
2476 /* Reset __FUNCT__ in case the user does not define it themselves */
2477 #undef __FUNCT__
2478 #define __FUNCT__ "User provided function"
2479 
2480 #endif
2481