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