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