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