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