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