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