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