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