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