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