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