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