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