Actual source code: petscsys.h

petsc-3.3-p7 2013-05-11
  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: */
  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. For --prefix installs the ${PETSC_ARCH}/
 12:    does not exist and petscconf.h is in the same directory as the other PETSc include files.
 13: */
 14: #include "petscconf.h"
 15: #include "petscfix.h"

 17: #if defined(PETSC_DESIRE_FEATURE_TEST_MACROS)
 18: /*
 19:    Feature test macros must be included before headers defined by IEEE Std 1003.1-2001
 20:    We only turn these in PETSc source files that require them by setting PETSC_DESIRE_FEATURE_TEST_MACROS
 21: */
 22: #if defined(PETSC__POSIX_C_SOURCE_200112L) && !defined(_POSIX_C_SOURCE)
 23: #define _POSIX_C_SOURCE 200112L
 24: #endif
 25: #if defined(PETSC__BSD_SOURCE) && !defined(_BSD_SOURCE)
 26: #define _BSD_SOURCE
 27: #endif
 28: #if defined(PETSC__GNU_SOURCE) && !defined(_GNU_SOURCE)
 29: #define _GNU_SOURCE
 30: #endif
 31: #endif

 33: /* ========================================================================== */
 34: /* 
 35:    This facilitates using the C version of PETSc from C++ and the
 36:    C++ version from C. Use --with-c-support --with-clanguage=c++ with ./configure for the latter)
 37: */
 38: #if defined(PETSC_CLANGUAGE_CXX) && !defined(PETSC_USE_EXTERN_CXX) && !defined(__cplusplus)
 39: #error "PETSc configured with --with-clanguage=c++ and NOT --with-c-support - it can be used only with a C++ compiler"
 40: #endif

 42: #if defined(__cplusplus)
 43: #  define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_CXX
 44: #else
 45: #  define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_C
 46: #endif

 48: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES) /* For Win32 shared libraries */
 49: #  define  __declspec(dllexport)
 50: #  define PETSC_DLLIMPORT __declspec(dllimport)
 51: #elif defined(PETSC_USE_VISIBILITY)
 52: #  define  __attribute__((visibility ("default")))
 53: #  define PETSC_DLLIMPORT __attribute__((visibility ("default")))
 54: #else
 55: #  define 
 56: #  define PETSC_DLLIMPORT
 57: #endif

 59: #if defined(petsc_EXPORTS)      /* CMake defines this when building the shared library */
 60: #  define PETSC_VISIBILITY_PUBLIC 
 61: #else  /* Win32 users need this to import symbols from petsc.dll */
 62: #  define PETSC_VISIBILITY_PUBLIC PETSC_DLLIMPORT
 63: #endif

 65: #if defined(PETSC_USE_EXTERN_CXX) && defined(__cplusplus)
 66: #define PETSC_EXTERN extern "C" PETSC_VISIBILITY_PUBLIC
 67: #define PETSC_EXTERN_TYPEDEF extern "C"
 68: #else
 69: #define PETSC_EXTERN extern PETSC_VISIBILITY_PUBLIC
 70: #define PETSC_EXTERN_TYPEDEF
 71: #endif
 72: /* ========================================================================== */
 73: /* 
 74:    Current PETSc version number and release date. Also listed in
 75:     Web page
 76:     src/docs/tex/manual/intro.tex,
 77:     src/docs/tex/manual/manual.tex.
 78:     src/docs/website/index.html.
 79: */
 80:  #include petscversion.h
 81: #define PETSC_AUTHOR_INFO        "       The PETSc Team\n    petsc-maint@mcs.anl.gov\n http://www.mcs.anl.gov/petsc/\n"
 82: #if (PETSC_VERSION_RELEASE == 1)
 83: #define PetscGetVersion(version,len) PetscSNPrintf(version,len,"Petsc Release Version %d.%d.%d, Patch %d, %s ", \
 84:                                          PETSC_VERSION_MAJOR,PETSC_VERSION_MINOR, PETSC_VERSION_SUBMINOR, \
 85:                                          PETSC_VERSION_PATCH,PETSC_VERSION_PATCH_DATE)
 86: #else
 87: #define PetscGetVersion(version,len) PetscSNPrintf(version,len,"Petsc Development HG revision: %s  HG Date: %s", \
 88:                                         PETSC_VERSION_HG, PETSC_VERSION_DATE_HG)
 89: #endif

 91: /*MC
 92:     PetscGetVersion - Gets the PETSc version information in a string.

 94:     Input Parameter:
 95: .   len - length of the string

 97:     Output Parameter:
 98: .   version - version string

100:     Level: developer

102:     Usage:
103:     char version[256];
104:     PetscGetVersion(version,256);CHKERRQ(ierr)

106:     Fortran Note:
107:     This routine is not supported in Fortran.

109: .seealso: PetscGetProgramName()

111: M*/

113: /* ========================================================================== */

115: /*
116:     Defines the interface to MPI allowing the use of all MPI functions.

118:     PETSc does not use the C++ binding of MPI at ALL. The following flag
119:     makes sure the C++ bindings are not included. The C++ bindings REQUIRE
120:     putting mpi.h before ANY C++ include files, we cannot control this
121:     with all PETSc users. Users who want to use the MPI C++ bindings can include 
122:     mpicxx.h directly in their code
123: */
124: #define MPICH_SKIP_MPICXX 1
125: #define OMPI_SKIP_MPICXX 1
126: #include "mpi.h"

128: /*
129:     Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler 
130:     see the top of mpicxx.h in the MPICH2 distribution.
131: */
132: #include <stdio.h>

134: /* MSMPI on 32bit windows requires this yukky hack - that breaks MPI standard compliance */
135: #if !defined(MPIAPI)
136: #define MPIAPI
137: #endif


140: /*MC
141:     PetscErrorCode - datatype used for return error code from all PETSc functions

143:     Level: beginner

145: .seealso: CHKERRQ, SETERRQ
146: M*/
147: typedef int PetscErrorCode;

149: /*MC

151:     PetscClassId - A unique id used to identify each PETSc class.
152:          (internal integer in the data structure used for error
153:          checking). These are all computed by an offset from the lowest
154:          one, PETSC_SMALLEST_CLASSID. 

156:     Level: advanced

158: .seealso: PetscClassIdRegister(), PetscLogEventRegister(), PetscHeaderCreate()
159: M*/
160: typedef int PetscClassId;


163: /*MC
164:     PetscMPIInt - datatype used to represent 'int' parameters to MPI functions.

166:     Level: intermediate

168:     Notes: usually this is the same as PetscInt, but if PETSc was built with --with-64-bit-indices but 
169:            standard C/Fortran integers are 32 bit then this is NOT the same as PetscInt it remains 32 bit

171:     PetscMPIIntCheck(a) checks if the given PetscInt a will fit in a PetscMPIInt, if not it generates a 
172:       PETSC_ERR_ARG_OUTOFRANGE error.

174:     PetscMPIInt b = PetscMPIIntCast(a) checks if the given PetscInt a will fit in a PetscMPIInt, if not it 
175:       generates a PETSC_ERR_ARG_OUTOFRANGE error.

177: .seealso: PetscBLASInt, PetscInt

179: M*/
180: typedef int PetscMPIInt;

182: /*MC
183:     PetscEnum - datatype used to pass enum types within PETSc functions.

185:     Level: intermediate

187: .seealso: PetscOptionsGetEnum(), PetscOptionsEnum(), PetscBagRegisterEnum()
188: M*/
189: typedef enum { ENUM_DUMMY } PetscEnum;

191: /*MC
192:     PetscInt - PETSc type that represents integer - used primarily to
193:       represent size of arrays and indexing into arrays. Its size can be configured with the option
194:       --with-64-bit-indices - to be either 32bit or 64bit [default 32 bit ints]

196:    Level: intermediate

198: .seealso: PetscScalar, PetscBLASInt, PetscMPIInt
199: M*/
200: #if (PETSC_SIZEOF_LONG_LONG == 8)
201: typedef long long Petsc64bitInt;
202: #elif defined(PETSC_HAVE___INT64)
203: typedef __int64 Petsc64bitInt;
204: #else
205: typedef unknown64bit Petsc64bitInt
206: #endif
207: #if defined(PETSC_USE_64BIT_INDICES)
208: typedef Petsc64bitInt PetscInt;
209: #define MPIU_INT MPI_LONG_LONG_INT
210: #else
211: typedef int PetscInt;
212: #define MPIU_INT MPI_INT
213: #endif


216: /*MC
217:     PetscBLASInt - datatype used to represent 'int' parameters to BLAS/LAPACK functions.

219:     Level: intermediate

221:     Notes: usually this is the same as PetscInt, but if PETSc was built with --with-64-bit-indices but
222:            standard C/Fortran integers are 32 bit then this is NOT the same as PetscInt it remains 32 bit
223:            (except on very rare BLAS/LAPACK implementations that support 64 bit integers see the note below).

225:     PetscBLASIntCheck(a) checks if the given PetscInt a will fit in a PetscBLASInt, if not it generates a
226:       PETSC_ERR_ARG_OUTOFRANGE error.

228:     PetscBLASInt b = PetscBLASIntCast(a) checks if the given PetscInt a will fit in a PetscBLASInt, if not it
229:       generates a PETSC_ERR_ARG_OUTOFRANGE error

231:     Installation Notes: The 64bit versions of MATLAB ship with BLAS and LAPACK that use 64 bit integers for sizes etc,
232:      if you run ./configure with the option
233:      --with-blas-lapack-lib=[/Applications/MATLAB_R2010b.app/bin/maci64/libmwblas.dylib,/Applications/MATLAB_R2010b.app/bin/maci64/libmwlapack.dylib]
234:      but you need to also use --known-64-bit-blas-indices. 

236:         MKL also ships with 64 bit integer versions of the BLAS and LAPACK, if you select those you must also ./configure with --known-64-bit-blas-indices

238:      Developer Notes: Eventually ./configure should automatically determine the size of the integers used by BLAS/LAPACK.

240:      External packages such as hypre, ML, SuperLU etc do not provide any support for passing 64 bit integers to BLAS/LAPACK so cannot
241:      be used with PETSc if you have set PetscBLASInt to long int.

243: .seealso: PetscMPIInt, PetscInt

245: M*/
246: #if defined(PETSC_HAVE_64BIT_BLAS_INDICES)
247: typedef Petsc64bitInt PetscBLASInt;
248: #else
249: typedef int PetscBLASInt;
250: #endif

252: /*EC

254:     PetscPrecision - indicates what precision the object is using. This is currently not used.

256:     Level: advanced

258: .seealso: PetscObjectSetPrecision()
259: E*/
260: typedef enum { PETSC_PRECISION_SINGLE=4,PETSC_PRECISION_DOUBLE=8 } PetscPrecision;
261: PETSC_EXTERN const char *PetscPrecisions[];

263: /* 
264:     For the rare cases when one needs to send a size_t object with MPI
265: */
266: #if (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_INT)
267: #define MPIU_SIZE_T MPI_INT
268: #elif  (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_LONG)
269: #define MPIU_SIZE_T MPI_LONG
270: #elif  (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_LONG_LONG)
271: #define MPIU_SIZE_T MPI_LONG_LONG_INT
272: #else
273: #error "Unknown size for size_t! Send us a bugreport at petsc-maint@mcs.anl.gov"
274: #endif


277: /*
278:       You can use PETSC_STDOUT as a replacement of stdout. You can also change
279:     the value of PETSC_STDOUT to redirect all standard output elsewhere
280: */
281: PETSC_EXTERN FILE* PETSC_STDOUT;

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: PETSC_EXTERN FILE* PETSC_STDERR;

289: /*MC
290:     PetscUnlikely - hints the compiler that the given condition is usually FALSE

292:     Synopsis:
293:     PetscBool  PetscUnlikely(PetscBool  cond)

295:     Not Collective

297:     Input Parameters:
298: .   cond - condition or expression

300:     Note: This returns the same truth value, it is only a hint to compilers that the resulting
301:     branch is unlikely.

303:     Level: advanced

305: .seealso: PetscLikely(), CHKERRQ
306: M*/

308: /*MC
309:     PetscLikely - hints the compiler that the given condition is usually TRUE

311:     Synopsis:
312:     PetscBool  PetscUnlikely(PetscBool  cond)

314:     Not Collective

316:     Input Parameters:
317: .   cond - condition or expression

319:     Note: This returns the same truth value, it is only a hint to compilers that the resulting
320:     branch is likely.

322:     Level: advanced

324: .seealso: PetscUnlikely()
325: M*/
326: #if defined(PETSC_HAVE_BUILTIN_EXPECT)
327: #  define PetscUnlikely(cond)   __builtin_expect(!!(cond),0)
328: #  define PetscLikely(cond)     __builtin_expect(!!(cond),1)
329: #else
330: #  define PetscUnlikely(cond)   (cond)
331: #  define PetscLikely(cond)     (cond)
332: #endif

334: /*
335:     Defines some elementary mathematics functions and constants.
336: */
337:  #include petscmath.h

339: /*
340:     Declare extern C stuff after including external header files
341: */


344: /*
345:        Basic PETSc constants
346: */

348: /*E
349:     PetscBool  - Logical variable. Actually an int in C and a logical in Fortran. 

351:    Level: beginner

353:    Developer Note: Why have PetscBool , why not use bool in C? The problem is that K and R C, C99 and C++ all have different mechanisms for 
354:       boolean values. It is not easy to have a simple macro that that will work properly in all circumstances with all three mechanisms.

356: E*/
357: typedef enum { PETSC_FALSE,PETSC_TRUE } PetscBool;
358: PETSC_EXTERN const char *PetscBools[];

360: /*E
361:     PetscCopyMode  - Determines how an array passed to certain functions is copied or retained

363:    Level: beginner

365: $   PETSC_COPY_VALUES - the array values are copied into new space, the user is free to reuse or delete the passed in array
366: $   PETSC_OWN_POINTER - the array values are NOT copied, the object takes ownership of the array and will free it later, the user cannot change or 
367: $                       delete the array. The array MUST have been obtained with PetscMalloc(). Hence this mode cannot be used in Fortran.
368: $   PETSC_USE_POINTER - the array values are NOT copied, the object uses the array but does NOT take ownership of the array. The user cannot use 
369:                         the array but the user must delete the array after the object is destroyed.

371: E*/
372: typedef enum { PETSC_COPY_VALUES, PETSC_OWN_POINTER, PETSC_USE_POINTER} PetscCopyMode;
373: PETSC_EXTERN const char *PetscCopyModes[];

375: /*MC
376:     PETSC_FALSE - False value of PetscBool 

378:     Level: beginner

380:     Note: Zero integer

382: .seealso: PetscBool , PETSC_TRUE
383: M*/

385: /*MC
386:     PETSC_TRUE - True value of PetscBool 

388:     Level: beginner

390:     Note: Nonzero integer

392: .seealso: PetscBool , PETSC_FALSE
393: M*/

395: /*MC
396:     PETSC_NULL - standard way of passing in a null or array or pointer

398:    Level: beginner

400:    Notes: accepted by many PETSc functions to not set a parameter and instead use
401:           some default

403:           This macro does not exist in Fortran; you must use PETSC_NULL_INTEGER, 
404:           PETSC_NULL_DOUBLE_PRECISION, PETSC_NULL_FUNCTION, PETSC_NULL_OBJECT etc

406:   Developer Note: Why have PETSC_NULL, why not just use NULL? The problem is that NULL is defined in different include files under 
407:       different versions of Unix. It is tricky to insure the correct include file is always included.

409: .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_IGNORE, PETSC_DETERMINE

411: M*/
412: #define PETSC_NULL           0

414: /*MC
415:     PETSC_IGNORE - same as PETSC_NULL, means PETSc will ignore this argument

417:    Level: beginner

419:    Note: accepted by many PETSc functions to not set a parameter and instead use
420:           some default

422:    Fortran Notes: This macro does not exist in Fortran; you must use PETSC_NULL_INTEGER, 
423:           PETSC_NULL_DOUBLE_PRECISION etc

425: .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_NULL, PETSC_DETERMINE

427: M*/
428: #define PETSC_IGNORE         PETSC_NULL

430: /*MC
431:     PETSC_DECIDE - standard way of passing in integer or floating point parameter
432:        where you wish PETSc to use the default.

434:    Level: beginner

436: .seealso: PETSC_NULL, PETSC_DEFAULT, PETSC_IGNORE, PETSC_DETERMINE

438: M*/
439: #define PETSC_DECIDE  -1

441: /*MC
442:     PETSC_DETERMINE - standard way of passing in integer or floating point parameter
443:        where you wish PETSc to compute the required value.

445:    Level: beginner


448:    Developer Note: I would like to use const PetscInt PETSC_DETERMINE = PETSC_DECIDE; but for
449:      some reason this is not allowed by the standard even though PETSC_DECIDE is a constant value.

451: .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_IGNORE, PETSC_NULL, VecSetSizes()

453: M*/
454: #define PETSC_DETERMINE PETSC_DECIDE

456: /*MC
457:     PETSC_DEFAULT - standard way of passing in integer or floating point parameter
458:        where you wish PETSc to use the default.

460:    Level: beginner

462:    Fortran Notes: You need to use PETSC_DEFAULT_INTEGER or PETSC_DEFAULT_DOUBLE_PRECISION.

464: .seealso: PETSC_DECIDE, PETSC_NULL, PETSC_IGNORE, PETSC_DETERMINE

466: M*/
467: #define PETSC_DEFAULT  -2

469: /*MC
470:     PETSC_COMM_WORLD - the equivalent of the MPI_COMM_WORLD communicator which represents
471:            all the processs that PETSc knows about. 

473:    Level: beginner

475:    Notes: By default PETSC_COMM_WORLD and MPI_COMM_WORLD are identical unless you wish to 
476:           run PETSc on ONLY a subset of MPI_COMM_WORLD. In that case create your new (smaller)
477:           communicator, call it, say comm, and set PETSC_COMM_WORLD = comm BEFORE calling
478:           PetscInitialize()

480: .seealso: PETSC_COMM_SELF

482: M*/
483: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;

485: /*MC
486:     PETSC_COMM_SELF - This is always MPI_COMM_SELF

488:    Level: beginner

490: .seealso: PETSC_COMM_WORLD

492: M*/
493: #define PETSC_COMM_SELF MPI_COMM_SELF

495: PETSC_EXTERN PetscBool PetscInitializeCalled;
496: PETSC_EXTERN PetscBool PetscFinalizeCalled;
497: PETSC_EXTERN PetscBool PetscCUSPSynchronize;

499: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm));
500: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*);
501: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm*);

503: /*MC
504:    PetscMalloc - Allocates memory

506:    Synopsis:
507:    PetscErrorCode PetscMalloc(size_t m,void **result)

509:    Not Collective

511:    Input Parameter:
512: .  m - number of bytes to allocate

514:    Output Parameter:
515: .  result - memory allocated

517:    Level: beginner

519:    Notes: Memory is always allocated at least double aligned

521:           If you request memory of zero size it will allocate no space and assign the pointer to 0; PetscFree() will 
522:           properly handle not freeing the null pointer.

524: .seealso: PetscFree(), PetscNew()

526:   Concepts: memory allocation

528: M*/
529: #define PetscMalloc(a,b)  ((a != 0) ? (*PetscTrMalloc)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__,__SDIR__,(void**)(b)) : (*(b) = 0,0) )

531: /*MC
532:    PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment

534:    Synopsis:
535:    void *PetscAddrAlign(void *addr)

537:    Not Collective

539:    Input Parameters:
540: .  addr - address to align (any pointer type)

542:    Level: developer

544: .seealso: PetscMallocAlign()

546:   Concepts: memory allocation
547: M*/
548: #define PetscAddrAlign(a) (void*)((((PETSC_UINTPTR_T)(a))+(PETSC_MEMALIGN-1)) & ~(PETSC_MEMALIGN-1))

550: /*MC
551:    PetscMalloc2 - Allocates 2 chunks of  memory both aligned to PETSC_MEMALIGN

553:    Synopsis:
554:    PetscErrorCode PetscMalloc2(size_t m1,type, t1,void **r1,size_t m2,type t2,void **r2)

556:    Not Collective

558:    Input Parameter:
559: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
560: .  t1 - type of first memory elements 
561: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
562: -  t2 - type of second memory elements

564:    Output Parameter:
565: +  r1 - memory allocated in first chunk
566: -  r2 - memory allocated in second chunk

568:    Level: developer

570: .seealso: PetscFree(), PetscNew(), PetscMalloc()

572:   Concepts: memory allocation

574: M*/
575: #if defined(PETSC_USE_DEBUG)
576: #define PetscMalloc2(m1,t1,r1,m2,t2,r2) (PetscMalloc((m1)*sizeof(t1),r1) || PetscMalloc((m2)*sizeof(t2),r2))
577: #else
578: #define PetscMalloc2(m1,t1,r1,m2,t2,r2) ((*(r2) = 0,PetscMalloc((m1)*sizeof(t1)+(m2)*sizeof(t2)+(PETSC_MEMALIGN-1),r1)) || (*(r2) = (t2*)PetscAddrAlign(*(r1)+m1),0))
579: #endif

581: /*MC
582:    PetscMalloc3 - Allocates 3 chunks of  memory  all aligned to PETSC_MEMALIGN

584:    Synopsis:
585:    PetscErrorCode PetscMalloc3(size_t m1,type, t1,void **r1,size_t m2,type t2,void **r2,size_t m3,type t3,void **r3)

587:    Not Collective

589:    Input Parameter:
590: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
591: .  t1 - type of first memory elements 
592: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
593: .  t2 - type of second memory elements
594: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
595: -  t3 - type of third memory elements

597:    Output Parameter:
598: +  r1 - memory allocated in first chunk
599: .  r2 - memory allocated in second chunk
600: -  r3 - memory allocated in third chunk

602:    Level: developer

604: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree3()

606:   Concepts: memory allocation

608: M*/
609: #if defined(PETSC_USE_DEBUG)
610: #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))
611: #else
612: #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)) \
613:                                                   || (*(r2) = (t2*)PetscAddrAlign(*(r1)+m1),*(r3) = (t3*)PetscAddrAlign(*(r2)+m2),0))
614: #endif

616: /*MC
617:    PetscMalloc4 - Allocates 4 chunks of  memory  all aligned to PETSC_MEMALIGN

619:    Synopsis:
620:    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)

622:    Not Collective

624:    Input Parameter:
625: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
626: .  t1 - type of first memory elements 
627: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
628: .  t2 - type of second memory elements
629: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
630: .  t3 - type of third memory elements
631: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
632: -  t4 - type of fourth memory elements

634:    Output Parameter:
635: +  r1 - memory allocated in first chunk
636: .  r2 - memory allocated in second chunk
637: .  r3 - memory allocated in third chunk
638: -  r4 - memory allocated in fourth chunk

640:    Level: developer

642: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree3(), PetscFree4()

644:   Concepts: memory allocation

646: M*/
647: #if defined(PETSC_USE_DEBUG)
648: #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))
649: #else
650: #define PetscMalloc4(m1,t1,r1,m2,t2,r2,m3,t3,r3,m4,t4,r4)               \
651:   ((*(r2) = 0, *(r3) = 0, *(r4) = 0,PetscMalloc((m1)*sizeof(t1)+(m2)*sizeof(t2)+(m3)*sizeof(t3)+(m4)*sizeof(t4)+3*(PETSC_MEMALIGN-1),r1)) \
652:    || (*(r2) = (t2*)PetscAddrAlign(*(r1)+m1),*(r3) = (t3*)PetscAddrAlign(*(r2)+m2),*(r4) = (t4*)PetscAddrAlign(*(r3)+m3),0))
653: #endif

655: /*MC
656:    PetscMalloc5 - Allocates 5 chunks of  memory all aligned to PETSC_MEMALIGN

658:    Synopsis:
659:    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)

661:    Not Collective

663:    Input Parameter:
664: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
665: .  t1 - type of first memory elements 
666: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
667: .  t2 - type of second memory elements
668: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
669: .  t3 - type of third memory elements
670: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
671: .  t4 - type of fourth memory elements
672: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
673: -  t5 - type of fifth memory elements

675:    Output Parameter:
676: +  r1 - memory allocated in first chunk
677: .  r2 - memory allocated in second chunk
678: .  r3 - memory allocated in third chunk
679: .  r4 - memory allocated in fourth chunk
680: -  r5 - memory allocated in fifth chunk

682:    Level: developer

684: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree3(), PetscFree4(), PetscFree5()

686:   Concepts: memory allocation

688: M*/
689: #if defined(PETSC_USE_DEBUG)
690: #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))
691: #else
692: #define PetscMalloc5(m1,t1,r1,m2,t2,r2,m3,t3,r3,m4,t4,r4,m5,t5,r5)      \
693:   ((*(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)) \
694:    || (*(r2) = (t2*)PetscAddrAlign(*(r1)+m1),*(r3) = (t3*)PetscAddrAlign(*(r2)+m2),*(r4) = (t4*)PetscAddrAlign(*(r3)+m3),*(r5) = (t5*)PetscAddrAlign(*(r4)+m4),0))
695: #endif


698: /*MC
699:    PetscMalloc6 - Allocates 6 chunks of  memory all aligned to PETSC_MEMALIGN

701:    Synopsis:
702:    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)

704:    Not Collective

706:    Input Parameter:
707: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
708: .  t1 - type of first memory elements 
709: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
710: .  t2 - type of second memory elements
711: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
712: .  t3 - type of third memory elements
713: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
714: .  t4 - type of fourth memory elements
715: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
716: .  t5 - type of fifth memory elements
717: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
718: -  t6 - type of sixth memory elements

720:    Output Parameter:
721: +  r1 - memory allocated in first chunk
722: .  r2 - memory allocated in second chunk
723: .  r3 - memory allocated in third chunk
724: .  r4 - memory allocated in fourth chunk
725: .  r5 - memory allocated in fifth chunk
726: -  r6 - memory allocated in sixth chunk

728:    Level: developer

730: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree3(), PetscFree4(), PetscFree5(), PetscFree6()

732:   Concepts: memory allocation

734: M*/
735: #if defined(PETSC_USE_DEBUG)
736: #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))
737: #else
738: #define PetscMalloc6(m1,t1,r1,m2,t2,r2,m3,t3,r3,m4,t4,r4,m5,t5,r5,m6,t6,r6) \
739:   ((*(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)) \
740:    || (*(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))
741: #endif

743: /*MC
744:    PetscMalloc7 - Allocates 7 chunks of  memory all aligned to PETSC_MEMALIGN

746:    Synopsis:
747:    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)

749:    Not Collective

751:    Input Parameter:
752: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
753: .  t1 - type of first memory elements 
754: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
755: .  t2 - type of second memory elements
756: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
757: .  t3 - type of third memory elements
758: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
759: .  t4 - type of fourth memory elements
760: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
761: .  t5 - type of fifth memory elements
762: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
763: .  t6 - type of sixth memory elements
764: .  m7 - number of elements to allocate in 7th chunk  (may be zero)
765: -  t7 - type of sixth memory elements

767:    Output Parameter:
768: +  r1 - memory allocated in first chunk
769: .  r2 - memory allocated in second chunk
770: .  r3 - memory allocated in third chunk
771: .  r4 - memory allocated in fourth chunk
772: .  r5 - memory allocated in fifth chunk
773: .  r6 - memory allocated in sixth chunk
774: -  r7 - memory allocated in seventh chunk

776:    Level: developer

778: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree3(), PetscFree4(), PetscFree5(), PetscFree6(), PetscFree7()

780:   Concepts: memory allocation

782: M*/
783: #if defined(PETSC_USE_DEBUG)
784: #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))
785: #else
786: #define PetscMalloc7(m1,t1,r1,m2,t2,r2,m3,t3,r3,m4,t4,r4,m5,t5,r5,m6,t6,r6,m7,t7,r7) \
787:   ((*(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)) \
788:    || (*(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))
789: #endif

791: /*MC
792:    PetscNew - Allocates memory of a particular type, zeros the memory! Aligned to PETSC_MEMALIGN

794:    Synopsis:
795:    PetscErrorCode PetscNew(struct type,((type *))result)

797:    Not Collective

799:    Input Parameter:
800: .  type - structure name of space to be allocated. Memory of size sizeof(type) is allocated

802:    Output Parameter:
803: .  result - memory allocated

805:    Level: beginner

807: .seealso: PetscFree(), PetscMalloc(), PetscNewLog()

809:   Concepts: memory allocation

811: M*/
812: #define PetscNew(A,b)      (PetscMalloc(sizeof(A),(b)) || PetscMemzero(*(b),sizeof(A)))

814: /*MC
815:    PetscNewLog - Allocates memory of a particular type, zeros the memory! Aligned to PETSC_MEMALIGN. Associates the memory allocated 
816:          with the given object using PetscLogObjectMemory().

818:    Synopsis:
819:    PetscErrorCode PetscNewLog(PetscObject obj,struct type,((type *))result)

821:    Not Collective

823:    Input Parameter:
824: +  obj - object memory is logged to
825: -  type - structure name of space to be allocated. Memory of size sizeof(type) is allocated

827:    Output Parameter:
828: .  result - memory allocated

830:    Level: developer

832: .seealso: PetscFree(), PetscMalloc(), PetscNew(), PetscLogObjectMemory()

834:   Concepts: memory allocation

836: M*/
837: #define PetscNewLog(o,A,b) (PetscNew(A,b) || ((o) ? PetscLogObjectMemory(o,sizeof(A)) : 0))

839: /*MC
840:    PetscFree - Frees memory

842:    Synopsis:
843:    PetscErrorCode PetscFree(void *memory)

845:    Not Collective

847:    Input Parameter:
848: .   memory - memory to free (the pointer is ALWAYS set to 0 upon sucess)

850:    Level: beginner

852:    Notes: Memory must have been obtained with PetscNew() or PetscMalloc()

854: .seealso: PetscNew(), PetscMalloc(), PetscFreeVoid()

856:   Concepts: memory allocation

858: M*/
859: #define PetscFree(a)   ((a) && ((*PetscTrFree)((void*)(a),__LINE__,PETSC_FUNCTION_NAME,__FILE__,__SDIR__) || ((a) = 0,0)))

861: /*MC
862:    PetscFreeVoid - Frees memory

864:    Synopsis:
865:    void PetscFreeVoid(void *memory)

867:    Not Collective

869:    Input Parameter:
870: .   memory - memory to free

872:    Level: beginner

874:    Notes: This is different from PetscFree() in that no error code is returned

876: .seealso: PetscFree(), PetscNew(), PetscMalloc()

878:   Concepts: memory allocation

880: M*/
881: #define PetscFreeVoid(a) ((*PetscTrFree)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__,__SDIR__),(a) = 0)


884: /*MC
885:    PetscFree2 - Frees 2 chunks of memory obtained with PetscMalloc2()

887:    Synopsis:
888:    PetscErrorCode PetscFree2(void *memory1,void *memory2)

890:    Not Collective

892:    Input Parameter:
893: +   memory1 - memory to free
894: -   memory2 - 2nd memory to free

896:    Level: developer

898:    Notes: Memory must have been obtained with PetscMalloc2()

900: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree()

902:   Concepts: memory allocation

904: M*/
905: #if defined(PETSC_USE_DEBUG)
906: #define PetscFree2(m1,m2)   (PetscFree(m2) || PetscFree(m1))
907: #else
908: #define PetscFree2(m1,m2)   ((m2)=0, PetscFree(m1))
909: #endif

911: /*MC
912:    PetscFree3 - Frees 3 chunks of memory obtained with PetscMalloc3()

914:    Synopsis:
915:    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)

917:    Not Collective

919:    Input Parameter:
920: +   memory1 - memory to free
921: .   memory2 - 2nd memory to free
922: -   memory3 - 3rd memory to free

924:    Level: developer

926:    Notes: Memory must have been obtained with PetscMalloc3()

928: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3()

930:   Concepts: memory allocation

932: M*/
933: #if defined(PETSC_USE_DEBUG)
934: #define PetscFree3(m1,m2,m3)   (PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
935: #else
936: #define PetscFree3(m1,m2,m3)   ((m3)=0,(m2)=0,PetscFree(m1))
937: #endif

939: /*MC
940:    PetscFree4 - Frees 4 chunks of memory obtained with PetscMalloc4()

942:    Synopsis:
943:    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)

945:    Not Collective

947:    Input Parameter:
948: +   m1 - memory to free
949: .   m2 - 2nd memory to free
950: .   m3 - 3rd memory to free
951: -   m4 - 4th memory to free

953:    Level: developer

955:    Notes: Memory must have been obtained with PetscMalloc4()

957: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4()

959:   Concepts: memory allocation

961: M*/
962: #if defined(PETSC_USE_DEBUG)
963: #define PetscFree4(m1,m2,m3,m4)   (PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
964: #else
965: #define PetscFree4(m1,m2,m3,m4)   ((m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
966: #endif

968: /*MC
969:    PetscFree5 - Frees 5 chunks of memory obtained with PetscMalloc5()

971:    Synopsis:
972:    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)

974:    Not Collective

976:    Input Parameter:
977: +   m1 - memory to free
978: .   m2 - 2nd memory to free
979: .   m3 - 3rd memory to free
980: .   m4 - 4th memory to free
981: -   m5 - 5th memory to free

983:    Level: developer

985:    Notes: Memory must have been obtained with PetscMalloc5()

987: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5()

989:   Concepts: memory allocation

991: M*/
992: #if defined(PETSC_USE_DEBUG)
993: #define PetscFree5(m1,m2,m3,m4,m5)   (PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
994: #else
995: #define PetscFree5(m1,m2,m3,m4,m5)   ((m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
996: #endif


999: /*MC
1000:    PetscFree6 - Frees 6 chunks of memory obtained with PetscMalloc6()

1002:    Synopsis:
1003:    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)

1005:    Not Collective

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: .   m5 - 5th memory to free
1013: -   m6 - 6th memory to free


1016:    Level: developer

1018:    Notes: Memory must have been obtained with PetscMalloc6()

1020: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6()

1022:   Concepts: memory allocation

1024: M*/
1025: #if defined(PETSC_USE_DEBUG)
1026: #define PetscFree6(m1,m2,m3,m4,m5,m6)   (PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1027: #else
1028: #define PetscFree6(m1,m2,m3,m4,m5,m6)   ((m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
1029: #endif

1031: /*MC
1032:    PetscFree7 - Frees 7 chunks of memory obtained with PetscMalloc7()

1034:    Synopsis:
1035:    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)

1037:    Not Collective

1039:    Input Parameter:
1040: +   m1 - memory to free
1041: .   m2 - 2nd memory to free
1042: .   m3 - 3rd memory to free
1043: .   m4 - 4th memory to free
1044: .   m5 - 5th memory to free
1045: .   m6 - 6th memory to free
1046: -   m7 - 7th memory to free


1049:    Level: developer

1051:    Notes: Memory must have been obtained with PetscMalloc7()

1053: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1054:           PetscMalloc7()

1056:   Concepts: memory allocation

1058: M*/
1059: #if defined(PETSC_USE_DEBUG)
1060: #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   (PetscFree(m7) || PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1061: #else
1062: #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   ((m7)=0,(m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
1063: #endif

1065: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,int,const char[],const char[],const char[],void**);
1066: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[],const char[]);
1067: PETSC_EXTERN PetscErrorCode PetscMallocSet(PetscErrorCode (*)(size_t,int,const char[],const char[],const char[],void**),PetscErrorCode (*)(void*,int,const char[],const char[],const char[]));
1068: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);

1070: /*
1071:     PetscLogDouble variables are used to contain double precision numbers
1072:   that are not used in the numerical computations, but rather in logging,
1073:   timing etc.
1074: */
1075: typedef double PetscLogDouble;
1076: #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE

1078: /*
1079:    Routines for tracing memory corruption/bleeding with default PETSc  memory allocation
1080: */
1081: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1082: PETSC_EXTERN PetscErrorCode PetscMallocDumpLog(FILE *);
1083: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1084: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1085: PETSC_EXTERN PetscErrorCode PetscMallocDebug(PetscBool);
1086: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[],const char[]);
1087: PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLog(void);
1088: PETSC_EXTERN PetscErrorCode PetscMallocGetDumpLog(PetscBool*);

1090: /*E
1091:     PetscDataType - Used for handling different basic data types.

1093:    Level: beginner

1095:    Developer comment: It would be nice if we could always just use MPI Datatypes, why can we not?

1097: .seealso: PetscBinaryRead(), PetscBinaryWrite(), PetscDataTypeToMPIDataType(),
1098:           PetscDataTypeGetSize()

1100: E*/
1101: typedef enum {PETSC_INT = 0,PETSC_DOUBLE = 1,PETSC_COMPLEX = 2, PETSC_LONG = 3 ,PETSC_SHORT = 4,PETSC_FLOAT = 5,
1102:               PETSC_CHAR = 6,PETSC_BIT_LOGICAL = 7,PETSC_ENUM = 8,PETSC_BOOL=9, PETSC___FLOAT128 = 10} PetscDataType;
1103: PETSC_EXTERN const char *PetscDataTypes[];

1105: #if defined(PETSC_USE_COMPLEX)
1106: #define  PETSC_SCALAR  PETSC_COMPLEX
1107: #else
1108: #if defined(PETSC_USE_REAL_SINGLE)
1109: #define  PETSC_SCALAR  PETSC_FLOAT
1110: #elif defined(PETSC_USE_REAL___FLOAT128)
1111: #define  PETSC_SCALAR  PETSC___FLOAT128
1112: #else
1113: #define  PETSC_SCALAR  PETSC_DOUBLE
1114: #endif
1115: #endif
1116: #if defined(PETSC_USE_REAL_SINGLE)
1117: #define  PETSC_REAL  PETSC_FLOAT
1118: #elif defined(PETSC_USE_REAL___FLOAT128)
1119: #define  PETSC_REAL  PETSC___FLOAT128
1120: #else
1121: #define  PETSC_REAL  PETSC_DOUBLE
1122: #endif
1123: #define  PETSC_FORTRANADDR  PETSC_LONG

1125: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1126: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1127: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);

1129: /*
1130:     Basic memory and string operations. These are usually simple wrappers
1131:    around the basic Unix system calls, but a few of them have additional
1132:    functionality and/or error checking.
1133: */
1134: PETSC_EXTERN PetscErrorCode PetscBitMemcpy(void*,PetscInt,const void*,PetscInt,PetscInt,PetscDataType);
1135: PETSC_EXTERN PetscErrorCode PetscMemmove(void*,void *,size_t);
1136: PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1137: PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1138: PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],int*,char ***);
1139: PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1140: PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool  *);
1141: PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool  *);
1142: PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1143: PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1144: PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1145: PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1146: PETSC_EXTERN PetscErrorCode PetscStrncat(char[],const char[],size_t);
1147: PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1148: PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1149: PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1150: PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1151: PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1152: PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1153: PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1154: PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1155: PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1156: PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1157: PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1158: PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);

1160: /*S
1161:     PetscToken - 'Token' used for managing tokenizing strings

1163:   Level: intermediate

1165: .seealso: PetscTokenCreate(), PetscTokenFind(), PetscTokenDestroy()
1166: S*/
1167: typedef struct _p_PetscToken* PetscToken;

1169: PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1170: PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1171: PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);

1173: /*
1174:    These are  MPI operations for MPI_Allreduce() etc
1175: */
1176: PETSC_EXTERN MPI_Op PetscMaxSum_Op;
1177: #if (defined(PETSC_USE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128)
1178: PETSC_EXTERN MPI_Op MPIU_SUM;
1179: #else
1180: #define MPIU_SUM MPI_SUM
1181: #endif
1182: #if defined(PETSC_USE_REAL___FLOAT128)
1183: PETSC_EXTERN MPI_Op MPIU_MAX;
1184: PETSC_EXTERN MPI_Op MPIU_MIN;
1185: #else
1186: #define MPIU_MAX MPI_MAX
1187: #define MPIU_MIN MPI_MIN
1188: #endif
1189: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);

1191: PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1192: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);

1194: /*S
1195:      PetscObject - any PETSc object, PetscViewer, Mat, Vec, KSP etc

1197:    Level: beginner

1199:    Note: This is the base class from which all PETSc objects are derived from.

1201: .seealso:  PetscObjectDestroy(), PetscObjectView(), PetscObjectGetName(), PetscObjectSetName(), PetscObjectReference(), PetscObjectDereferenc()
1202: S*/
1203: typedef struct _p_PetscObject* PetscObject;

1205: /*S
1206:      PetscFList - Linked list of functions, possibly stored in dynamic libraries, accessed
1207:       by string name

1209:    Level: advanced

1211: .seealso:  PetscFListAdd(), PetscFListDestroy(), PetscOpFlist
1212: S*/
1213: typedef struct _n_PetscFList *PetscFList;

1215: /*S
1216:      PetscOpFList - Linked list of operations on objects, implemented by functions, possibly stored in dynamic libraries, 
1217:                     accessed by string op name together with string argument types.

1219:    Level: advanced

1221: .seealso:  PetscFList, PetscOpFListAdd(), PetscOpFListFind(), PetscOpFListDestroy()
1222: S*/
1223: typedef struct _n_PetscOpFList *PetscOpFList;

1225: /*E
1226:   PetscFileMode - Access mode for a file.

1228:   Level: beginner

1230:   FILE_MODE_READ - open a file at its beginning for reading

1232:   FILE_MODE_WRITE - open a file at its beginning for writing (will create if the file does not exist)

1234:   FILE_MODE_APPEND - open a file at end for writing

1236:   FILE_MODE_UPDATE - open a file for updating, meaning for reading and writing

1238:   FILE_MODE_APPEND_UPDATE - open a file for updating, meaning for reading and writing, at the end

1240: .seealso: PetscViewerFileSetMode()
1241: E*/
1242: typedef enum {FILE_MODE_READ, FILE_MODE_WRITE, FILE_MODE_APPEND, FILE_MODE_UPDATE, FILE_MODE_APPEND_UPDATE} PetscFileMode;

1244:  #include petscviewer.h
1245:  #include petscoptions.h

1247: #define PETSC_SMALLEST_CLASSID  1211211
1248: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1249: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1250: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);

1252: /*
1253:    Routines that get memory usage information from the OS
1254: */
1255: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1256: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1257: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1258: PETSC_EXTERN PetscErrorCode PetscMemoryShowUsage(PetscViewer,const char[]);

1260: PETSC_EXTERN PetscErrorCode PetscInfoAllow(PetscBool ,const char []);
1261: PETSC_EXTERN PetscErrorCode PetscGetTime(PetscLogDouble*);
1262: PETSC_EXTERN PetscErrorCode PetscGetCPUTime(PetscLogDouble*);
1263: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);

1265: /* 
1266:    Initialization of PETSc
1267: */
1268: PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1269: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1270: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1271: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1272: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1273: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1274: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1275: PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1276: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1277: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);

1279: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1280: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(const char[]);

1282: PETSC_EXTERN MPI_Comm PETSC_COMM_LOCAL_WORLD;
1283: PETSC_EXTERN PetscErrorCode PetscHMPIMerge(PetscMPIInt,PetscErrorCode (*)(void*),void*);
1284: PETSC_EXTERN PetscErrorCode PetscHMPISpawn(PetscMPIInt);
1285: PETSC_EXTERN PetscErrorCode PetscHMPIFinalize(void);
1286: PETSC_EXTERN PetscErrorCode PetscHMPIRun(MPI_Comm,PetscErrorCode (*)(MPI_Comm,void *),void*);
1287: PETSC_EXTERN PetscErrorCode PetscHMPIRunCtx(MPI_Comm,PetscErrorCode (*)(MPI_Comm,void*,void *),void*);
1288: PETSC_EXTERN PetscErrorCode PetscHMPIFree(MPI_Comm,void*);
1289: PETSC_EXTERN PetscErrorCode PetscHMPIMalloc(MPI_Comm,size_t,void**);

1291: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1292: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1293: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1294: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);

1296: /*
1297:      These are so that in extern C code we can caste function pointers to non-extern C
1298:    function pointers. Since the regular C++ code expects its function pointers to be C++
1299: */
1300: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1301: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1302: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);

1304: /*
1305:    PetscTryMethod - Queries an object for a method, if it exists then calls it.
1306:               These are intended to be used only inside PETSc functions.

1308:    Level: developer
1309:    
1310: .seealso: PetscUseMethod()
1311: */
1312: #define  PetscTryMethod(obj,A,B,C) \
1313:   0;{ PetscErrorCode (*f)B, __ierr; \
1314:     __PetscObjectQueryFunction((PetscObject)obj,A,(PetscVoidStarFunction)&f);CHKERRQ(__ierr); \
1315:     if (f) {__(*f)C;CHKERRQ(__ierr);}\
1316:   }

1318: /*
1319:    PetscUseMethod - Queries an object for a method, if it exists then calls it, otherwise generates an error.
1320:               These are intended to be used only inside PETSc functions.

1322:    Level: developer
1323:    
1324: .seealso: PetscTryMethod()
1325: */
1326: #define  PetscUseMethod(obj,A,B,C) \
1327:   0;{ PetscErrorCode (*f)B, __ierr; \
1328:     __PetscObjectQueryFunction((PetscObject)obj,A,(PetscVoidStarFunction)&f);CHKERRQ(__ierr); \
1329:     if (f) {__(*f)C;CHKERRQ(__ierr);}\
1330:     else SETERRQ1(((PetscObject)obj)->comm,PETSC_ERR_SUP,"Cannot locate function %s in object",A); \
1331:   }

1333: /*
1334:     Functions that can act on any PETSc object.
1335: */
1336: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1337: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1338: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1339: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1340: PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1341: PETSC_EXTERN PetscErrorCode PetscObjectSetPrecision(PetscObject,PetscPrecision);
1342: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1343: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1344: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1345: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer,const char[]);
1346: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1347: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1348: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1349: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1350: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1351: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1352: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt *);
1353: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1354: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1355: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1356: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject *);
1357: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction(PetscObject,const char[],const char[],void (*)(void));
1358: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1359: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1360: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt *);
1361: PETSC_EXTERN PetscErrorCode PetscObjectAddOptionsHandler(PetscObject,PetscErrorCode (*)(PetscObject,void*),PetscErrorCode (*)(PetscObject,void*),void*);
1362: PETSC_EXTERN PetscErrorCode PetscObjectProcessOptionsHandlers(PetscObject);
1363: PETSC_EXTERN PetscErrorCode PetscObjectDestroyOptionsHandlers(PetscObject);
1364: PETSC_EXTERN PetscErrorCode PetscObjectsGetGlobalNumbering(MPI_Comm,PetscInt,PetscObject*,PetscInt*,PetscInt*);

1366: /*MC
1367:    PetscObjectComposeFunctionDynamic - Associates a function with a given PETSc object. 
1368:                        
1369:     Synopsis:
1370:     PetscErrorCode PetscObjectComposeFunctionDynamic(PetscObject obj,const char name[],const char fname[],void *ptr)

1372:    Logically Collective on PetscObject

1374:    Input Parameters:
1375: +  obj - the PETSc object; this must be cast with a (PetscObject), for example, 
1376:          PetscObjectCompose((PetscObject)mat,...);
1377: .  name - name associated with the child function
1378: .  fname - name of the function
1379: -  ptr - function pointer (or PETSC_NULL if using dynamic libraries)

1381:    Level: advanced


1384:    Notes:
1385:    To remove a registered routine, pass in a PETSC_NULL rname and fnc().

1387:    PetscObjectComposeFunctionDynamic() can be used with any PETSc object (such as
1388:    Mat, Vec, KSP, SNES, etc.) or any user-provided object. 

1390:    The composed function must be wrapped in a EXTERN_C_BEGIN/END for this to
1391:    work in C++/complex with dynamic link libraries (./configure options --with-shared-libraries --with-dynamic-loading)
1392:    enabled.

1394:    Concepts: objects^composing functions
1395:    Concepts: composing functions
1396:    Concepts: functions^querying
1397:    Concepts: objects^querying
1398:    Concepts: querying objects

1400: .seealso: PetscObjectQueryFunction()
1401: M*/
1402: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1403: #define PetscObjectComposeFunctionDynamic(a,b,c,d) PetscObjectComposeFunction(a,b,c,0)
1404: #else
1405: #define PetscObjectComposeFunctionDynamic(a,b,c,d) PetscObjectComposeFunction(a,b,c,(PetscVoidFunction)(d))
1406: #endif

1408: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction(PetscObject,const char[],void (**)(void));
1409: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1410: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1411: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1412: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1413: PETSC_EXTERN PetscErrorCode PetscObjectAMSPublish(PetscObject);
1414: PETSC_EXTERN PetscErrorCode PetscObjectUnPublish(PetscObject);
1415: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1416: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1417: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1418: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1419: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1420: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1421: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1422: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);

1424: /*
1425:     Defines PETSc error handling.
1426: */
1427:  #include petscerror.h

1429: /*S
1430:      PetscOList - Linked list of PETSc objects, each accessable by string name

1432:    Level: developer

1434:    Notes: Used by PetscObjectCompose() and PetscObjectQuery() 

1436: .seealso:  PetscOListAdd(), PetscOListDestroy(), PetscOListFind(), PetscObjectCompose(), PetscObjectQuery(), PetscFList
1437: S*/
1438: typedef struct _n_PetscOList *PetscOList;

1440: PETSC_EXTERN PetscErrorCode PetscOListDestroy(PetscOList*);
1441: PETSC_EXTERN PetscErrorCode PetscOListFind(PetscOList,const char[],PetscObject*);
1442: PETSC_EXTERN PetscErrorCode PetscOListReverseFind(PetscOList,PetscObject,char**,PetscBool*);
1443: PETSC_EXTERN PetscErrorCode PetscOListAdd(PetscOList *,const char[],PetscObject);
1444: PETSC_EXTERN PetscErrorCode PetscOListRemoveReference(PetscOList *,const char[]);
1445: PETSC_EXTERN PetscErrorCode PetscOListDuplicate(PetscOList,PetscOList *);

1447: /*
1448:     Dynamic library lists. Lists of names of routines in objects or in dynamic 
1449:   link libraries that will be loaded as needed.
1450: */
1451: PETSC_EXTERN PetscErrorCode PetscFListAdd(PetscFList*,const char[],const char[],void (*)(void));
1452: PETSC_EXTERN PetscErrorCode PetscFListDestroy(PetscFList*);
1453: PETSC_EXTERN PetscErrorCode PetscFListFind(PetscFList,MPI_Comm,const char[],PetscBool,void (**)(void));
1454: PETSC_EXTERN PetscErrorCode PetscFListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFList,const char[]);
1455: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1456: #define    PetscFListAddDynamic(a,b,p,c) PetscFListAdd(a,b,p,0)
1457: #else
1458: #define    PetscFListAddDynamic(a,b,p,c) PetscFListAdd(a,b,p,(void (*)(void))c)
1459: #endif
1460: PETSC_EXTERN PetscErrorCode PetscFListDuplicate(PetscFList,PetscFList *);
1461: PETSC_EXTERN PetscErrorCode PetscFListView(PetscFList,PetscViewer);
1462: PETSC_EXTERN PetscErrorCode PetscFListConcat(const char [],const char [],char []);
1463: PETSC_EXTERN PetscErrorCode PetscFListGet(PetscFList,char ***,int*);

1465: /*
1466:     Multiple dispatch operation function lists. Lists of names of routines with corresponding
1467:     argument type names with function pointers or in dynamic link libraries that will be loaded 
1468:     as needed.  Search on the op name and argument type names.
1469: */
1470: PETSC_EXTERN PetscErrorCode PetscOpFListAdd(MPI_Comm, PetscOpFList*,const char[],PetscVoidFunction, const char[], PetscInt, char*[]);
1471: PETSC_EXTERN PetscErrorCode PetscOpFListDestroy(PetscOpFList*);
1472: PETSC_EXTERN PetscErrorCode PetscOpFListFind(MPI_Comm, PetscOpFList, PetscVoidFunction*, const char[], PetscInt, char*[]);
1473: PETSC_EXTERN PetscErrorCode PetscOpFListView(PetscOpFList,PetscViewer);

1475: /*S
1476:      PetscDLLibrary - Linked list of dynamics libraries to search for functions

1478:    Level: advanced

1480:    --with-shared-libraries --with-dynamic-loading must be used with ./configure to use dynamic libraries

1482: .seealso:  PetscDLLibraryOpen()
1483: S*/
1484: typedef struct _n_PetscDLLibrary *PetscDLLibrary;
1485: PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1486: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1487: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1488: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1489: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1490: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1491: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1492: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);

1494: /*
1495:   PetscShell support.  Needs to be better documented.  
1496:   Logically it is an extension of PetscDLLXXX, PetscObjectCompose, etc.
1497: */
1498:  #include petscshell.h

1500: /*
1501:      Useful utility routines
1502: */
1503: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1504: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1505: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1506: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1507: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1508: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);

1510: /*
1511:     PetscNot - negates a logical type value and returns result as a PetscBool 

1513:     Notes: This is useful in cases like 
1514: $     int        *a;
1515: $     PetscBool  flag = PetscNot(a) 
1516:      where !a does not return a PetscBool  because we cannot provide a cast from int to PetscBool  in C.
1517: */
1518: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)

1520: /*
1521:     Defines basic graphics available from PETSc. 
1522: */
1523:  #include petscdraw.h

1525: /*
1526:     Defines the base data structures for all PETSc objects
1527: */
1528: #include "petsc-private/petscimpl.h"


1531: /*MC
1532:     PetscErrorPrintf - Prints error messages.

1534:    Synopsis:
1535:      PetscErrorCode (*PetscErrorPrintf)(const char format[],...);

1537:     Not Collective

1539:     Input Parameters:
1540: .   format - the usual printf() format string 

1542:    Options Database Keys:
1543: +    -error_output_stdout - cause error messages to be printed to stdout instead of the
1544:          (default) stderr
1545: -    -error_output_none to turn off all printing of error messages (does not change the way the 
1546:           error is handled.)

1548:    Notes: Use
1549: $     PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the 
1550: $                        error is handled.) and
1551: $     PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on
1552: $        of you can use your own function

1554:           Use
1555:      PETSC_STDERR = FILE* obtained from a file open etc. to have stderr printed to the file. 
1556:      PETSC_STDOUT = FILE* obtained from a file open etc. to have stdout printed to the file. 

1558:           Use
1559:       PetscPushErrorHandler() to provide your own error handler that determines what kind of messages to print

1561:    Level: developer

1563:     Fortran Note:
1564:     This routine is not supported in Fortran.

1566:     Concepts: error messages^printing
1567:     Concepts: printing^error messages

1569: .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscHelpPrintf(), PetscPrintf(), PetscErrorHandlerPush(), PetscVFPrintf(), PetscHelpPrintf()
1570: M*/
1571: PETSC_EXTERN PetscErrorCode (*PetscErrorPrintf)(const char[],...);

1573: /*MC
1574:     PetscHelpPrintf - Prints help messages.

1576:    Synopsis:
1577:      PetscErrorCode (*PetscHelpPrintf)(const char format[],...);

1579:     Not Collective

1581:     Input Parameters:
1582: .   format - the usual printf() format string 

1584:    Level: developer

1586:     Fortran Note:
1587:     This routine is not supported in Fortran.

1589:     Concepts: help messages^printing
1590:     Concepts: printing^help messages

1592: .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1593: M*/
1594: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);

1596: /*
1597:      Defines PETSc profiling.
1598: */
1599:  #include petsclog.h

1601: /*
1602:           For locking, unlocking and destroying AMS memories associated with  PETSc objects. ams.h is included in petscviewer.h
1603: */
1604: #if defined(PETSC_HAVE_AMS)
1605: PETSC_EXTERN PetscBool PetscAMSPublishAll;
1606: #define PetscObjectTakeAccess(obj)  ((((PetscObject)(obj))->amem == -1) ? 0 : AMS_Memory_take_access(((PetscObject)(obj))->amem))
1607: #define PetscObjectGrantAccess(obj) ((((PetscObject)(obj))->amem == -1) ? 0 : AMS_Memory_grant_access(((PetscObject)(obj))->amem))
1608: #define PetscObjectDepublish(obj)   ((((PetscObject)(obj))->amem == -1) ? 0 : AMS_Memory_destroy(((PetscObject)(obj))->amem));((PetscObject)(obj))->amem = -1;
1609: #else
1610: #define PetscObjectTakeAccess(obj)   0
1611: #define PetscObjectGrantAccess(obj)  0
1612: #define PetscObjectDepublish(obj)      0
1613: #endif

1615: /*
1616:       Simple PETSc parallel IO for ASCII printing
1617: */
1618: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1619: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1620: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1621: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1622: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1623: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1624: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);

1626: /* These are used internally by PETSc ASCII IO routines*/
1627: #include <stdarg.h>
1628: PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list);
1629: PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE*,const char[],va_list);
1630: PETSC_EXTERN PetscErrorCode PetscVFPrintfDefault(FILE*,const char[],va_list);

1632: #if defined(PETSC_HAVE_MATLAB_ENGINE)
1633: PETSC_EXTERN PetscErrorCode PetscVFPrintf_Matlab(FILE*,const char[],va_list);
1634: #endif

1636: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1637: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1638: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);

1640: #if defined(PETSC_HAVE_POPEN)
1641: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1642: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*);
1643: #endif

1645: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1646: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1647: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm);
1648: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1649: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1650: PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1651: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);

1653: PETSC_EXTERN PetscErrorCode PetscPopUpSelect(MPI_Comm,const char*,const char*,int,const char**,int*);

1655: /*S
1656:      PetscContainer - Simple PETSc object that contains a pointer to any required data

1658:    Level: advanced

1660: .seealso:  PetscObject, PetscContainerCreate()
1661: S*/
1662: PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1663: typedef struct _p_PetscContainer*  PetscContainer;
1664: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void **);
1665: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void *);
1666: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1667: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer *);
1668: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void*));

1670: /*
1671:    For use in debuggers 
1672: */
1673: PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1674: PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1675: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1676: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1677: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);

1679: #if defined(PETSC_HAVE_MEMORY_H)
1680: #include <memory.h>
1681: #endif
1682: #if defined(PETSC_HAVE_STDLIB_H)
1683: #include <stdlib.h>
1684: #endif
1685: #if defined(PETSC_HAVE_STRINGS_H)
1686: #include <strings.h>
1687: #endif
1688: #if defined(PETSC_HAVE_STRING_H)
1689: #include <string.h>
1690: #endif

1692: #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1693: #include <xmmintrin.h>
1694: #endif
1695: #if defined(PETSC_HAVE_STDINT_H)
1696: #include <stdint.h>
1697: #endif

1701: /*@C
1702:    PetscMemcpy - Copies n bytes, beginning at location b, to the space
1703:    beginning at location a. The two memory regions CANNOT overlap, use
1704:    PetscMemmove() in that case.

1706:    Not Collective

1708:    Input Parameters:
1709: +  b - pointer to initial memory space
1710: -  n - length (in bytes) of space to copy

1712:    Output Parameter:
1713: .  a - pointer to copy space

1715:    Level: intermediate

1717:    Compile Option:
1718:     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used 
1719:                                   for memory copies on double precision values.
1720:     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used 
1721:                                   for memory copies on double precision values.
1722:     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used 
1723:                                   for memory copies on double precision values.

1725:    Note:
1726:    This routine is analogous to memcpy().

1728:    Developer Note: this is inlined for fastest performance

1730:   Concepts: memory^copying
1731:   Concepts: copying^memory
1732:   
1733: .seealso: PetscMemmove()

1735: @*/
1736: PETSC_STATIC_INLINE PetscErrorCode  PetscMemcpy(void *a,const void *b,size_t n)
1737: {
1738: #if defined(PETSC_USE_DEBUG)
1739:   unsigned long al = (unsigned long) a,bl = (unsigned long) b;
1740:   unsigned long nl = (unsigned long) n;
1742:   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1743:   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1744: #else
1746: #endif
1747:   if (a != b) {
1748: #if defined(PETSC_USE_DEBUG)
1749:     if ((al > bl && (al - bl) < nl) || (bl - al) < nl)  SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1750:               or make sure your copy regions and lengths are correct. \n\
1751:               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1752: #endif
1753: #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1754:    if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1755:       size_t len = n/sizeof(PetscScalar);
1756: #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1757:       PetscBLASInt one = 1,blen = PetscBLASIntCast(len);
1758:       BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one);
1759: #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1760:       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1761: #else
1762:       size_t      i;
1763:       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1764:       for (i=0; i<len; i++) y[i] = x[i];
1765: #endif
1766:     } else {
1767:       memcpy((char*)(a),(char*)(b),n);
1768:     }
1769: #else
1770:     memcpy((char*)(a),(char*)(b),n);
1771: #endif
1772:   }
1773:   return(0);
1774: }

1776: /*@C
1777:    PetscMemzero - Zeros the specified memory.

1779:    Not Collective

1781:    Input Parameters:
1782: +  a - pointer to beginning memory location
1783: -  n - length (in bytes) of memory to initialize

1785:    Level: intermediate

1787:    Compile Option:
1788:    PETSC_PREFER_BZERO - on certain machines (the IBM RS6000) the bzero() routine happens
1789:   to be faster than the memset() routine. This flag causes the bzero() routine to be used.

1791:    Developer Note: this is inlined for fastest performance

1793:    Concepts: memory^zeroing
1794:    Concepts: zeroing^memory

1796: .seealso: PetscMemcpy()
1797: @*/
1798: PETSC_STATIC_INLINE PetscErrorCode  PetscMemzero(void *a,size_t n)
1799: {
1800:   if (n > 0) {
1801: #if defined(PETSC_USE_DEBUG)
1802:     if (!a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer");
1803: #endif
1804: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1805:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1806:       size_t      i,len = n/sizeof(PetscScalar);
1807:       PetscScalar *x = (PetscScalar*)a;
1808:       for (i=0; i<len; i++) x[i] = 0.0;
1809:     } else {
1810: #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1811:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1812:       PetscInt len = n/sizeof(PetscScalar);
1813:       fortranzero_(&len,(PetscScalar*)a);
1814:     } else {
1815: #endif
1816: #if defined(PETSC_PREFER_BZERO)
1817:       bzero((char *)a,n);
1818: #else
1819:       memset((char*)a,0,n);
1820: #endif
1821: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1822:     }
1823: #endif
1824:   }
1825:   return 0;
1826: }

1828: /*MC
1829:    PetscPrefetchBlock - Prefetches a block of memory

1831:    Synopsis:
1832:     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)

1834:    Not Collective

1836:    Input Parameters:
1837: +  a - pointer to first element to fetch (any type but usually PetscInt or PetscScalar)
1838: .  n - number of elements to fetch
1839: .  rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors)
1840: -  t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note

1842:    Level: developer

1844:    Notes:
1845:    The last two arguments (rw and t) must be compile-time constants.

1847:    Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality.  Not all architectures offer
1848:    equivalent locality hints, but the following macros are always defined to their closest analogue.
1849: +  PETSC_PREFETCH_HINT_NTA - Non-temporal.  Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched).
1850: .  PETSC_PREFETCH_HINT_T0 - Fetch to all levels of cache and evict to the closest level.  Use this when the memory will be reused regularly despite necessary eviction from L1.
1851: .  PETSC_PREFETCH_HINT_T1 - Fetch to level 2 and higher (not L1).
1852: -  PETSC_PREFETCH_HINT_T2 - Fetch to high-level cache only.  (On many systems, T0 and T1 are equivalent.)

1854:    This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid
1855:    address).

1857:    Concepts: memory
1858: M*/
1859: #define PetscPrefetchBlock(a,n,rw,t) do {                               \
1860:     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
1861:     for ( ; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
1862:   } while (0)

1864: /*
1865:     Allows accessing MATLAB Engine
1866: */
1867:  #include petscmatlab.h

1869: /*
1870:       Determine if some of the kernel computation routines use
1871:    Fortran (rather than C) for the numerical calculations. On some machines
1872:    and compilers (like complex numbers) the Fortran version of the routines
1873:    is faster than the C/C++ versions. The flag --with-fortran-kernels
1874:    should be used with ./configure to turn these on.
1875: */
1876: #if defined(PETSC_USE_FORTRAN_KERNELS)

1878: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1879: #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1880: #endif

1882: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
1883: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
1884: #endif

1886: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1887: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
1888: #endif

1890: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1891: #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
1892: #endif

1894: #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
1895: #define PETSC_USE_FORTRAN_KERNEL_NORM
1896: #endif

1898: #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
1899: #define PETSC_USE_FORTRAN_KERNEL_MAXPY
1900: #endif

1902: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1903: #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
1904: #endif

1906: #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1907: #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
1908: #endif

1910: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
1911: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
1912: #endif

1914: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1915: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
1916: #endif

1918: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
1919: #define PETSC_USE_FORTRAN_KERNEL_MDOT
1920: #endif

1922: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
1923: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
1924: #endif

1926: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
1927: #define PETSC_USE_FORTRAN_KERNEL_AYPX
1928: #endif

1930: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
1931: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
1932: #endif

1934: #endif

1936: /*
1937:     Macros for indicating code that should be compiled with a C interface,
1938:    rather than a C++ interface. Any routines that are dynamically loaded
1939:    (such as the PCCreate_XXX() routines) must be wrapped so that the name
1940:    mangler does not change the functions symbol name. This just hides the 
1941:    ugly extern "C" {} wrappers.
1942: */
1943: #if defined(__cplusplus)
1945: #define EXTERN_C_BEGIN extern "C" {
1946: #define EXTERN_C_END }
1947: #else
1949: #define EXTERN_C_BEGIN 
1950: #define EXTERN_C_END 
1951: #endif

1953: /* --------------------------------------------------------------------*/

1955: /*MC
1956:     MPI_Comm - the basic object used by MPI to determine which processes are involved in a 
1957:         communication

1959:    Level: beginner

1961:    Note: This manual page is a place-holder because MPICH does not have a manual page for MPI_Comm

1963: .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
1964: M*/

1966: /*MC
1967:     PetscScalar - PETSc type that represents either a double precision real number, a double precision
1968:        complex number, a single precision real number, a long double or an int - if the code is configured 
1969:        with --with-scalar-type=real,complex --with-precision=single,double,longdouble,int,matsingle


1972:    Level: beginner

1974: .seealso: PetscReal, PassiveReal, PassiveScalar, MPIU_SCALAR, PetscInt
1975: M*/

1977: /*MC
1978:     PetscReal - PETSc type that represents a real number version of PetscScalar

1980:    Level: beginner

1982: .seealso: PetscScalar, PassiveReal, PassiveScalar
1983: M*/

1985: /*MC
1986:     PassiveScalar - PETSc type that represents a PetscScalar
1987:    Level: beginner

1989:     This is the same as a PetscScalar except in code that is automatically differentiated it is
1990:    treated as a constant (not an indendent or dependent variable)

1992: .seealso: PetscReal, PassiveReal, PetscScalar
1993: M*/

1995: /*MC
1996:     PassiveReal - PETSc type that represents a PetscReal

1998:    Level: beginner

2000:     This is the same as a PetscReal except in code that is automatically differentiated it is
2001:    treated as a constant (not an indendent or dependent variable)

2003: .seealso: PetscScalar, PetscReal, PassiveScalar
2004: M*/

2006: /*MC
2007:     MPIU_SCALAR - MPI datatype corresponding to PetscScalar

2009:    Level: beginner

2011:     Note: In MPI calls that require an MPI datatype that matches a PetscScalar or array of PetscScalars
2012:           pass this value

2014: .seealso: PetscReal, PassiveReal, PassiveScalar, PetscScalar, MPIU_INT
2015: M*/

2017: #if defined(PETSC_HAVE_MPIIO)
2018: #if !defined(PETSC_WORDS_BIGENDIAN)
2019: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2020: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2021: #else
2022: #define MPIU_File_write_all(a,b,c,d,e) MPI_File_write_all(a,b,c,d,e) 
2023: #define MPIU_File_read_all(a,b,c,d,e) MPI_File_read_all(a,b,c,d,e) 
2024: #endif
2025: #endif

2027: /* the following petsc_static_inline require petscerror.h */

2029: /* Limit MPI to 32-bits */
2030: #define PETSC_MPI_INT_MAX  2147483647
2031: #define PETSC_MPI_INT_MIN -2147483647
2032: /* Limit BLAS to 32-bits */
2033: #define PETSC_BLAS_INT_MAX  2147483647
2034: #define PETSC_BLAS_INT_MIN -2147483647
2035: /* On 32 bit systems HDF5 is limited by size of integer, because hsize_t is defined as size_t */
2036: #define PETSC_HDF5_INT_MAX  2147483647
2037: #define PETSC_HDF5_INT_MIN -2147483647

2039: #if defined(PETSC_USE_64BIT_INDICES)
2040: #define PetscMPIIntCheck(a)  if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Message too long for MPI")
2041: #define PetscBLASIntCheck(a)  if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK")
2042: #define PetscMPIIntCast(a) (PetscMPIInt)(a);PetscMPIIntCheck(a)
2043: #define PetscBLASIntCast(a) (PetscBLASInt)(a);PetscBLASIntCheck(a)

2045: #if (PETSC_SIZEOF_SIZE_T == 4)
2046: #define PetscHDF5IntCheck(a)  if ((a) > PETSC_HDF5_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for HDF5")
2047: #define PetscHDF5IntCast(a) (hsize_t)(a);PetscHDF5IntCheck(a)
2048: #else
2049: #define PetscHDF5IntCheck(a)
2050: #define PetscHDF5IntCast(a) a
2051: #endif

2053: #else
2054: #define PetscMPIIntCheck(a) 
2055: #define PetscBLASIntCheck(a) 
2056: #define PetscHDF5IntCheck(a)
2057: #define PetscMPIIntCast(a) a
2058: #define PetscBLASIntCast(a) a
2059: #define PetscHDF5IntCast(a) a
2060: #endif  

2062: /*
2063:      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2064: */
2065: #if defined(hz)
2066: #undef hz
2067: #endif

2069: /*  For arrays that contain filenames or paths */


2072: #if defined(PETSC_HAVE_LIMITS_H)
2073: #include <limits.h>
2074: #endif
2075: #if defined(PETSC_HAVE_SYS_PARAM_H)
2076: #include <sys/param.h>
2077: #endif
2078: #if defined(PETSC_HAVE_SYS_TYPES_H)
2079: #include <sys/types.h>
2080: #endif
2081: #if defined(MAXPATHLEN)
2082: #  define PETSC_MAX_PATH_LEN     MAXPATHLEN
2083: #elif defined(MAX_PATH)
2084: #  define PETSC_MAX_PATH_LEN     MAX_PATH
2085: #elif defined(_MAX_PATH)
2086: #  define PETSC_MAX_PATH_LEN     _MAX_PATH
2087: #else
2088: #  define PETSC_MAX_PATH_LEN     4096
2089: #endif

2091: /* Special support for C++ */
2092:  #include petscsys.hh

2094: /*MC

2096:     UsingFortran - Fortran can be used with PETSc in four distinct approaches

2098: $    1) classic Fortran 77 style
2099: $#include "finclude/petscXXX.h" to work with material from the XXX component of PETSc
2100: $       XXX variablename
2101: $      You cannot use this approach if you wish to use the Fortran 90 specific PETSc routines
2102: $      which end in F90; such as VecGetArrayF90()
2103: $
2104: $    2) classic Fortran 90 style
2105: $#include "finclude/petscXXX.h" 
2106: $#include "finclude/petscXXX.h90" to work with material from the XXX component of PETSc
2107: $       XXX variablename
2108: $
2109: $    3) Using Fortran modules
2110: $#include "finclude/petscXXXdef.h" 
2111: $         use petscXXXX
2112: $       XXX variablename
2113: $
2114: $    4) Use Fortran modules and Fortran data types for PETSc types
2115: $#include "finclude/petscXXXdef.h" 
2116: $         use petscXXXX
2117: $       type(XXX) variablename
2118: $      To use this approach you must ./configure PETSc with the additional
2119: $      option --with-fortran-datatypes You cannot use the type(XXX) declaration approach without using Fortran modules

2121:     Finally if you absolutely do not want to use any #include you can use either 

2123: $    3a) skip the #include BUT you cannot use any PETSc data type names like Vec, Mat, PetscInt, PetscErrorCode etc
2124: $        and you must declare the variables as integer, for example 
2125: $        integer variablename
2126: $
2127: $    4a) skip the #include, you use the object types like type(Vec) type(Mat) but cannot use the data type
2128: $        names like PetscErrorCode, PetscInt etc. again for those you must use integer

2130:    We recommend either 2 or 3. Approaches 2 and 3 provide type checking for most PETSc function calls; 4 has type checking 
2131: for only a few PETSc functions.

2133:    Fortran type checking with interfaces is strick, this means you cannot pass a scalar value when an array value
2134: is expected (even though it is legal Fortran). For example when setting a single value in a matrix with MatSetValues()
2135: you cannot have something like
2136: $      PetscInt row,col
2137: $      PetscScalar val
2138: $        ...
2139: $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)
2140: You must instead have 
2141: $      PetscInt row(1),col(1)
2142: $      PetscScalar val(1)
2143: $        ...
2144: $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)


2147:     See the example src/vec/vec/examples/tutorials/ex20f90.F90 for an example that can use all four approaches

2149:     Developer Notes: The finclude/petscXXXdef.h contain all the #defines (would be typedefs in C code) these
2150:      automatically include their predecessors; for example finclude/petscvecdef.h includes finclude/petscisdef.h

2152:      The finclude/petscXXXX.h contain all the parameter statements for that package. These automatically include
2153:      their finclude/petscXXXdef.h file but DO NOT automatically include their predecessors;  for example 
2154:      finclude/petscvec.h does NOT automatically include finclude/petscis.h

2156:      The finclude/ftn-custom/petscXXXdef.h90 are not intended to be used directly in code, they define the
2157:      Fortran data type type(XXX) (for example type(Vec)) when PETSc is ./configure with the --with-fortran-datatypes option.

2159:      The finclude/ftn-custom/petscXXX.h90 (not included directly by code) contain interface definitions for
2160:      the PETSc Fortran stubs that have different bindings then their C version (for example VecGetArrayF90).

2162:      The finclude/ftn-auto/petscXXX.h90 (not included directly by code) contain interface definitions generated
2163:      automatically by "make allfortranstubs".

2165:      The finclude/petscXXX.h90 includes the custom finclude/ftn-custom/petscXXX.h90 and if ./configure 
2166:      was run with --with-fortran-interfaces it also includes the finclude/ftn-auto/petscXXX.h90 These DO NOT automatically
2167:      include their predecessors

2169:     Level: beginner

2171: M*/

2173: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2174: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2175: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2176: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2177: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2178: PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);

2180: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2181: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2182: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2183: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2184: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2185: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt*,PetscInt*,PetscInt*);
2186: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2187: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2188: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2189: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2190: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2191: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2192: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2193: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt*,const PetscInt*,PetscInt,const PetscInt*,const PetscInt*,PetscInt*,PetscInt**,PetscInt**);

2195: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2196: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);

2198: /*J
2199:     PetscRandomType - String with the name of a PETSc randomizer
2200:        with an optional dynamic library name, for example
2201:        http://www.mcs.anl.gov/petsc/lib.a:myrandcreate()

2203:    Level: beginner

2205:    Notes: to use the SPRNG you must have ./configure PETSc
2206:    with the option --download-sprng

2208: .seealso: PetscRandomSetType(), PetscRandom
2209: J*/
2210: #define PetscRandomType char*
2211: #define PETSCRAND       "rand"
2212: #define PETSCRAND48     "rand48"
2213: #define PETSCSPRNG      "sprng"          

2215: /* Logging support */
2216: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;

2218: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(const char[]);

2220: /*S
2221:      PetscRandom - Abstract PETSc object that manages generating random numbers

2223:    Level: intermediate

2225:   Concepts: random numbers

2227: .seealso:  PetscRandomCreate(), PetscRandomGetValue(), PetscRandomType
2228: S*/
2229: typedef struct _p_PetscRandom*   PetscRandom;

2231: /* Dynamic creation and loading functions */
2232: PETSC_EXTERN PetscFList PetscRandomList;
2233: PETSC_EXTERN PetscBool PetscRandomRegisterAllCalled;

2235: PETSC_EXTERN PetscErrorCode PetscRandomRegisterAll(const char []);
2236: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],const char[],const char[],PetscErrorCode (*)(PetscRandom));
2237: PETSC_EXTERN PetscErrorCode PetscRandomRegisterDestroy(void);
2238: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, const PetscRandomType);
2239: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2240: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, const PetscRandomType*);
2241: PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom,char*);
2242: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);

2244: /*MC
2245:   PetscRandomRegisterDynamic - Adds a new PetscRandom component implementation

2247:   Synopsis:
2248:   PetscErrorCode PetscRandomRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(PetscRandom))

2250:   Not Collective

2252:   Input Parameters:
2253: + name        - The name of a new user-defined creation routine
2254: . path        - The path (either absolute or relative) of the library containing this routine
2255: . func_name   - The name of routine to create method context
2256: - create_func - The creation routine itself

2258:   Notes:
2259:   PetscRandomRegisterDynamic() may be called multiple times to add several user-defined randome number generators

2261:   If dynamic libraries are used, then the fourth input argument (routine_create) is ignored.

2263:   Sample usage:
2264: .vb
2265:     PetscRandomRegisterDynamic("my_rand","/home/username/my_lib/lib/libO/solaris/libmy.a", "MyPetscRandomtorCreate", MyPetscRandomtorCreate);
2266: .ve

2268:   Then, your random type can be chosen with the procedural interface via
2269: .vb
2270:     PetscRandomCreate(MPI_Comm, PetscRandom *);
2271:     PetscRandomSetType(PetscRandom,"my_random_name");
2272: .ve
2273:    or at runtime via the option
2274: .vb
2275:     -random_type my_random_name
2276: .ve

2278:   Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.

2280:          For an example of the code needed to interface your own random number generator see
2281:          src/sys/random/impls/rand/rand.c
2282:         
2283:   Level: advanced

2285: .keywords: PetscRandom, register
2286: .seealso: PetscRandomRegisterAll(), PetscRandomRegisterDestroy(), PetscRandomRegister()
2287: M*/
2288: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
2289: #define PetscRandomRegisterDynamic(a,b,c,d) PetscRandomRegister(a,b,c,0)
2290: #else
2291: #define PetscRandomRegisterDynamic(a,b,c,d) PetscRandomRegister(a,b,c,d)
2292: #endif

2294: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2295: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2296: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2297: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2298: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2299: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2300: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2301: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2302: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);

2304: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2305: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2306: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2307: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2308: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2309: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2310: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);

2312: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscDataType);
2313: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscDataType);
2314: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,void*,PetscInt,PetscDataType,PetscBool );
2315: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,void*,PetscInt,PetscDataType,PetscBool );
2316: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2317: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2318: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2319: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2320: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2321: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2322: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2323: PETSC_EXTERN PetscErrorCode PetscOpenSocket(char*,int,int*);
2324: PETSC_EXTERN PetscErrorCode PetscWebServe(MPI_Comm,int);

2326: /*
2327:    In binary files variables are stored using the following lengths,
2328:   regardless of how they are stored in memory on any one particular
2329:   machine. Use these rather then sizeof() in computing sizes for 
2330:   PetscBinarySeek().
2331: */
2332: #define PETSC_BINARY_INT_SIZE   (32/8)
2333: #define PETSC_BINARY_FLOAT_SIZE  (32/8)
2334: #define PETSC_BINARY_CHAR_SIZE  (8/8)
2335: #define PETSC_BINARY_SHORT_SIZE  (16/8)
2336: #define PETSC_BINARY_DOUBLE_SIZE  (64/8)
2337: #define PETSC_BINARY_SCALAR_SIZE  sizeof(PetscScalar)

2339: /*E
2340:   PetscBinarySeekType - argument to PetscBinarySeek()

2342:   Level: advanced

2344: .seealso: PetscBinarySeek(), PetscBinarySynchronizedSeek()
2345: E*/
2346: typedef enum {PETSC_BINARY_SEEK_SET = 0,PETSC_BINARY_SEEK_CUR = 1,PETSC_BINARY_SEEK_END = 2} PetscBinarySeekType;
2347: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2348: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2349: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);

2351: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2352: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool );
2353: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2354: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2355: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2356: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);

2358: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2359: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2360: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2361: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2362: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);

2364: PETSC_EXTERN PetscErrorCode PetscSSEIsEnabled(MPI_Comm,PetscBool  *,PetscBool  *);

2366: /*E
2367:   InsertMode - Whether entries are inserted or added into vectors or matrices

2369:   Level: beginner

2371: .seealso: VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2372:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(),
2373:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd()
2374: E*/
2375:  typedef enum {NOT_SET_VALUES, INSERT_VALUES, ADD_VALUES, MAX_VALUES, INSERT_ALL_VALUES, ADD_ALL_VALUES} InsertMode;

2377: /*MC
2378:     INSERT_VALUES - Put a value into a vector or matrix, overwrites any previous value

2380:     Level: beginner

2382: .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2383:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), ADD_VALUES,
2384:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES

2386: M*/

2388: /*MC
2389:     ADD_VALUES - Adds a value into a vector or matrix, if there previously was no value, just puts the
2390:                 value into that location

2392:     Level: beginner

2394: .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2395:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), INSERT_VALUES,
2396:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES

2398: M*/

2400: /*MC
2401:     MAX_VALUES - Puts the maximum of the scattered/gathered value and the current value into each location

2403:     Level: beginner

2405: .seealso: InsertMode, VecScatterBegin(), VecScatterEnd(), ADD_VALUES, INSERT_VALUES

2407: M*/

2409: /*S
2410:    PetscSubcomm - Context of MPI subcommunicators, used by PCREDUNDANT

2412:    Level: advanced

2414:    Concepts: communicator, create
2415: S*/
2416: typedef struct _n_PetscSubcomm* PetscSubcomm;

2418: struct _n_PetscSubcomm {
2419:   MPI_Comm   parent;      /* parent communicator */
2420:   MPI_Comm   dupparent;   /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2421:   MPI_Comm   comm;        /* this communicator */
2422:   PetscInt   n;           /* num of subcommunicators under the parent communicator */
2423:   PetscInt   color;       /* color of processors belong to this communicator */
2424: };

2426: typedef enum {PETSC_SUBCOMM_GENERAL=0,PETSC_SUBCOMM_CONTIGUOUS=1,PETSC_SUBCOMM_INTERLACED=2} PetscSubcommType;
2427: PETSC_EXTERN const char *PetscSubcommTypes[];

2429: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2430: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2431: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2432: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2433: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt,PetscMPIInt);

2435: #include <petscctable.h>


2438: /* Reset __FUNCT__ in case the user does not define it themselves */

2442: #endif