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