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