Actual source code: petscsys.h

  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
 12: */
 13: #include "petscconf.h"
 14: #include "petscfix.h"

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

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

 38: #if defined(__cplusplus)
 39: #  define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_CXX
 40: #else
 41: #  define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_C
 42: #endif

 44: #if defined(PETSC_USE_EXTERN_CXX) && defined(__cplusplus)
 45: #define PETSC_EXTERN_CXX_BEGIN extern "C" {
 46: #define PETSC_EXTERN_CXX_END  }
 47: #else
 48: #define PETSC_EXTERN_CXX_BEGIN
 49: #define PETSC_EXTERN_CXX_END
 50: #endif
 51: /* ========================================================================== */
 52: /* 
 53:    Current PETSc version number and release date. Also listed in
 54:     Web page
 55:     src/docs/tex/manual/intro.tex,
 56:     src/docs/tex/manual/manual.tex.
 57:     src/docs/website/index.html.
 58: */
 59:  #include petscversion.h
 60: #define PETSC_AUTHOR_INFO        "       The PETSc Team\n    petsc-maint@mcs.anl.gov\n http://www.mcs.anl.gov/petsc/\n"
 61: #if (PETSC_VERSION_RELEASE == 1)
 62: #define PetscGetVersion(version,len) PetscSNPrintf(version,len,"Petsc Release Version %d.%d.%d, Patch %d, %s ", \
 63:                                          PETSC_VERSION_MAJOR,PETSC_VERSION_MINOR, PETSC_VERSION_SUBMINOR, \
 64:                                          PETSC_VERSION_PATCH,PETSC_VERSION_PATCH_DATE)
 65: #else
 66: #define PetscGetVersion(version,len) PetscSNPrintf(version,len,"Petsc Development HG revision: %s  HG Date: %s", \
 67:                                         PETSC_VERSION_HG, PETSC_VERSION_DATE_HG)
 68: #endif

 70: /*MC
 71:     PetscGetVersion - Gets the PETSc version information in a string.

 73:     Input Parameter:
 74: .   len - length of the string

 76:     Output Parameter:
 77: .   version - version string

 79:     Level: developer

 81:     Usage:
 82:     char version[256];
 83:     PetscGetVersion(version,256);CHKERRQ(ierr)

 85:     Fortran Note:
 86:     This routine is not supported in Fortran.

 88: .seealso: PetscGetProgramName()

 90: M*/

 92: /* ========================================================================== */

 94: /*
 95:     Defines the interface to MPI allowing the use of all MPI functions.

 97:     PETSc does not use the C++ binding of MPI at ALL. The following flag
 98:     makes sure the C++ bindings are not included. The C++ bindings REQUIRE
 99:     putting mpi.h before ANY C++ include files, we cannot control this
100:     with all PETSc users. Users who want to use the MPI C++ bindings can include 
101:     mpicxx.h directly in their code
102: */
103: #define MPICH_SKIP_MPICXX 1
104: #define OMPI_SKIP_MPICXX 1
105: #include "mpi.h"

107: /*
108:     Yuck, we need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler 
109:     see the top of mpicxx.h in the MPICH2 distribution.

111:     The MPI STANDARD HAS TO BE CHANGED to prevent this nonsense.
112: */
113: #include <stdio.h>

115: /* MSMPI on 32bit windows requires this yukky hack - that breaks MPI standard compliance */
116: #if !defined(MPIAPI)
117: #define MPIAPI
118: #endif


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

124:     Level: beginner

126: .seealso: CHKERRQ, SETERRQ
127: M*/
128: typedef int PetscErrorCode;

130: /*MC

132:     PetscClassId - A unique id used to identify each PETSc class.
133:          (internal integer in the data structure used for error
134:          checking). These are all defined by an offset from the lowest
135:          one, PETSC_SMALLEST_CLASSID. 

137:     Level: advanced

139: .seealso: PetscClassIdRegister(), PetscLogEventRegister(), PetscHeaderCreate()
140: M*/
141: typedef int PetscClassId;


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

147:     Level: intermediate

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

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

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

159:     Developer Notes: The 64bit versions of MATLAB ship with BLAS and LAPACK that use 64 bit integers for sizes etc,
160:      if you run ./configure with the option
161:      --with-blas-lapack-lib=[/Applications/MATLAB_R2010b.app/bin/maci64/libmwblas.dylib,/Applications/MATLAB_R2010b.app/bin/maci64/libmwlapack.dylib]
162:      for example, you can change the int below to long int. Since MATLAB uses the MKL (Intel Math Libraries) it is likely one can
163:      purchase a 64 bit integer version of the MKL and use that with a  PetscBLASInt of long int.

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

168: .seealso: PetscMPIInt, PetscInt

170: M*/
171: typedef int PetscBLASInt;

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

176:     Level: intermediate

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

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

184:     PetscMPIInt b = PetscMPIIntCast(a) checks if the given PetscInt a will fit in a PetscMPIInt, if not it 
185:       generates a PETSC_ERR_ARG_OUTOFRANGE

187: .seealso: PetscBLASInt, PetscInt

189: M*/
190: typedef int PetscMPIInt;

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

195:     Level: intermediate

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

200:     PetscMPIInt b = PetscMPIIntCast(a) checks if the given PetscInt a will fit in a PetscMPIInt, if not it 
201:       generates a PETSC_ERR_ARG_OUTOFRANGE

203: .seealso: PetscOptionsGetEnum(), PetscOptionsEnum(), PetscBagRegisterEnum()
204: M*/
205: typedef enum { ENUM_DUMMY } PetscEnum;

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

212:    Level: intermediate

214: .seealso: PetscScalar, PetscBLASInt, PetscMPIInt
215: M*/
216: #if defined(PETSC_HAVE___INT64)
217: typedef __int64 Petsc64bitInt;
218: #else
219: typedef long long Petsc64bitInt;
220: #endif
221: #if defined(PETSC_USE_64BIT_INDICES)
222: typedef Petsc64bitInt PetscInt;
223: #define MPIU_INT MPI_LONG_LONG_INT
224: #else
225: typedef int PetscInt;
226: #define MPIU_INT MPI_INT
227: #endif


230: /*EC

232:     PetscPrecision - indicates what precision the object is using

234:     Level: advanced

236: .seealso: PetscObjectSetPrecision()
237: E*/
238: typedef enum { PETSC_PRECISION_SINGLE=4,PETSC_PRECISION_DOUBLE=8 } PetscPrecision;
239: extern const char *PetscPrecisions[];


242: /* 
243:     For the rare cases when one needs to send a size_t object with MPI
244: */
245: #if (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_INT)
246: #define MPIU_SIZE_T MPI_INT
247: #elif  (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_LONG)
248: #define MPIU_SIZE_T MPI_LONG
249: #elif  (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_LONG_LONG)
250: #define MPIU_SIZE_T MPI_LONG_LONG_INT
251: #else
252: #error "Unknown size for size_t! Send us a bugreport at petsc-maint@mcs.anl.gov"
253: #endif


256: /*
257:       You can use PETSC_STDOUT as a replacement of stdout. You can also change
258:     the value of PETSC_STDOUT to redirect all standard output elsewhere
259: */

261: extern FILE* PETSC_STDOUT;

263: /*
264:       You can use PETSC_STDERR as a replacement of stderr. You can also change
265:     the value of PETSC_STDERR to redirect all standard error elsewhere
266: */
267: extern FILE* PETSC_STDERR;

269: /*
270:       PETSC_ZOPEFD is used to send data to the PETSc webpage.  It can be used
271:     in conjunction with PETSC_STDOUT, or by itself.
272: */
273: extern FILE* PETSC_ZOPEFD;

275: #if !defined(PETSC_USE_EXTERN_CXX) && defined(__cplusplus)
276: /*MC
277:       PetscPolymorphicSubroutine - allows defining a C++ polymorphic version of 
278:             a PETSc function that remove certain optional arguments for a simplier user interface

280:    Synopsis:
281:    PetscPolymorphicSubroutine(Functionname,(arguments of C++ function),(arguments of C function))
282:  
283:      Not collective

285:    Level: developer

287:     Example:
288:       PetscPolymorphicSubroutine(VecNorm,(Vec x,PetscReal *r),(x,NORM_2,r)) generates the new routine
289:            PetscErrorCode VecNorm(Vec x,PetscReal *r) = VecNorm(x,NORM_2,r)

291: .seealso: PetscPolymorphicFunction()

293: M*/
294: #define PetscPolymorphicSubroutine(A,B,C) PETSC_STATIC_INLINE PetscErrorCode A B {return A C;}

296: /*MC
297:       PetscPolymorphicScalar - allows defining a C++ polymorphic version of 
298:             a PETSc function that replaces a PetscScalar * argument with a PetscScalar argument

300:    Synopsis:
301:    PetscPolymorphicScalar(Functionname,(arguments of C++ function),(arguments of C function))
302:  
303:    Not collective

305:    Level: developer

307:     Example:
308:       PetscPolymorphicScalar(VecAXPY,(PetscScalar _val,Vec x,Vec y),(&_Val,x,y)) generates the new routine
309:            PetscErrorCode VecAXPY(PetscScalar _val,Vec x,Vec y) = {PetscScalar _Val = _val; return VecAXPY(&_Val,x,y);}

311: .seealso: PetscPolymorphicFunction(),PetscPolymorphicSubroutine()

313: M*/
314: #define PetscPolymorphicScalar(A,B,C) PETSC_STATIC_INLINE PetscErrorCode A B {PetscScalar _Val = _val; return A C;}

316: /*MC
317:       PetscPolymorphicFunction - allows defining a C++ polymorphic version of 
318:             a PETSc function that remove certain optional arguments for a simplier user interface
319:             and returns the computed value (istead of an error code)

321:    Synopsis:
322:    PetscPolymorphicFunction(Functionname,(arguments of C++ function),(arguments of C function),return type,return variable name)
323:  
324:      Not collective

326:    Level: developer

328:     Example:
329:       PetscPolymorphicFunction(VecNorm,(Vec x,NormType t),(x,t,&r),PetscReal,r) generates the new routine
330:          PetscReal VecNorm(Vec x,NormType t) = {PetscReal r; VecNorm(x,t,&r); return r;}

332: .seealso: PetscPolymorphicSubroutine()

334: M*/
335: #define PetscPolymorphicFunction(A,B,C,D,E) PETSC_STATIC_INLINE D A B {D E; A C;return E;}

337: #else
338: #define PetscPolymorphicSubroutine(A,B,C)
339: #define PetscPolymorphicScalar(A,B,C)
340: #define PetscPolymorphicFunction(A,B,C,D,E)
341: #endif

343: /*MC
344:     PetscUnlikely - hints the compiler that the given condition is usually FALSE

346:     Synopsis:
347:     PetscBool  PetscUnlikely(PetscBool  cond)

349:     Not Collective

351:     Input Parameters:
352: .   cond - condition or expression

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

357:     Level: advanced

359: .seealso: PetscLikely(), CHKERRQ
360: M*/

362: /*MC
363:     PetscLikely - hints the compiler that the given condition is usually TRUE

365:     Synopsis:
366:     PetscBool  PetscUnlikely(PetscBool  cond)

368:     Not Collective

370:     Input Parameters:
371: .   cond - condition or expression

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

376:     Level: advanced

378: .seealso: PetscUnlikely()
379: M*/
380: #if defined(PETSC_HAVE_BUILTIN_EXPECT)
381: #  define PetscUnlikely(cond)   __builtin_expect(!!(cond),0)
382: #  define PetscLikely(cond)     __builtin_expect(!!(cond),1)
383: #else
384: #  define PetscUnlikely(cond)   (cond)
385: #  define PetscLikely(cond)     (cond)
386: #endif

388: /*
389:     Defines some elementary mathematics functions and constants.
390: */
391:  #include petscmath.h

393: /*
394:     Declare extern C stuff after including external header files
395: */

397: PETSC_EXTERN_CXX_BEGIN

399: /*
400:        Basic PETSc constants
401: */

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

406:    Level: beginner

408:    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 
409:       boolean values. It is not easy to have a simple macro that that will work properly in all circumstances with all three mechanisms.

411: E*/
412: typedef enum { PETSC_FALSE,PETSC_TRUE } PetscBool;
413: extern const char *PetscBools[];

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

418:    Level: beginner

420: $   PETSC_COPY_VALUES - the array values are copied into new space, the user is free to reuse or delete the passed in array
421: $   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 
422: $                       delete the array. The array MUST have been obtained with PetscMalloc(). Hence this mode cannot be used in Fortran.
423: $   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 
424:                         the array but the user must delete the array after the object is destroyed.

426: E*/
427: typedef enum { PETSC_COPY_VALUES, PETSC_OWN_POINTER, PETSC_USE_POINTER} PetscCopyMode;
428: extern const char *PetscCopyModes[];

430: /*MC
431:     PETSC_FALSE - False value of PetscBool 

433:     Level: beginner

435:     Note: Zero integer

437: .seealso: PetscBool , PETSC_TRUE
438: M*/

440: /*MC
441:     PETSC_TRUE - True value of PetscBool 

443:     Level: beginner

445:     Note: Nonzero integer

447: .seealso: PetscBool , PETSC_FALSE
448: M*/

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

453:    Level: beginner

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

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

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

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

466: M*/
467: #define PETSC_NULL           0

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

472:    Level: beginner

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

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

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

482: M*/
483: #define PETSC_IGNORE         PETSC_NULL

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

489:    Level: beginner

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

493: M*/
494: #define PETSC_DECIDE  -1

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

500:    Level: beginner


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

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

508: M*/
509: #define PETSC_DETERMINE PETSC_DECIDE

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

515:    Level: beginner

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

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

521: M*/
522: #define PETSC_DEFAULT  -2

524: /*MC
525:     PETSC_COMM_WORLD - the equivalent of the MPI_COMM_WORLD communicator which represents
526:            all the processs that PETSc knows about. 

528:    Level: beginner

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

535: .seealso: PETSC_COMM_SELF

537: M*/
538: extern MPI_Comm PETSC_COMM_WORLD;

540: /*MC
541:     PETSC_COMM_SELF - This is always MPI_COMM_SELF

543:    Level: beginner

545: .seealso: PETSC_COMM_WORLD

547: M*/
548: #define PETSC_COMM_SELF MPI_COMM_SELF

550: extern  PetscBool  PetscInitializeCalled;
551: extern  PetscBool  PetscFinalizeCalled;

553: extern PetscErrorCode  PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm));
554: extern PetscErrorCode  PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*);
555: extern PetscErrorCode  PetscCommDestroy(MPI_Comm*);

557: /*MC
558:    PetscMalloc - Allocates memory

560:    Synopsis:
561:    PetscErrorCode PetscMalloc(size_t m,void **result)

563:    Not Collective

565:    Input Parameter:
566: .  m - number of bytes to allocate

568:    Output Parameter:
569: .  result - memory allocated

571:    Level: beginner

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

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

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

580:   Concepts: memory allocation

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

585: /*MC
586:    PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment

588:    Synopsis:
589:    void *PetscAddrAlign(void *addr)

591:    Not Collective

593:    Input Parameters:
594: .  addr - address to align (any pointer type)

596:    Level: developer

598: .seealso: PetscMallocAlign()

600:   Concepts: memory allocation
601: M*/
602: #define PetscAddrAlign(a) (void*)((((PETSC_UINTPTR_T)(a))+(PETSC_MEMALIGN-1)) & ~(PETSC_MEMALIGN-1))

604: /*MC
605:    PetscMalloc2 - Allocates 2 chunks of  memory both aligned to PETSC_MEMALIGN

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

610:    Not Collective

612:    Input Parameter:
613: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
614: .  t1 - type of first memory elements 
615: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
616: -  t2 - type of second memory elements

618:    Output Parameter:
619: +  r1 - memory allocated in first chunk
620: -  r2 - memory allocated in second chunk

622:    Level: developer

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

626:   Concepts: memory allocation

628: M*/
629: #if defined(PETSC_USE_DEBUG)
630: #define PetscMalloc2(m1,t1,r1,m2,t2,r2) (PetscMalloc((m1)*sizeof(t1),r1) || PetscMalloc((m2)*sizeof(t2),r2))
631: #else
632: #define PetscMalloc2(m1,t1,r1,m2,t2,r2) ((*(r2) = 0,PetscMalloc((m1)*sizeof(t1)+(m2)*sizeof(t2)+(PETSC_MEMALIGN-1),r1)) \
633:                                          || (*(r2) = (t2*)PetscAddrAlign(*(r1)+m1),0))
634: #endif

636: /*MC
637:    PetscMalloc3 - Allocates 3 chunks of  memory  all aligned to PETSC_MEMALIGN

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

642:    Not Collective

644:    Input Parameter:
645: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
646: .  t1 - type of first memory elements 
647: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
648: .  t2 - type of second memory elements
649: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
650: -  t3 - type of third memory elements

652:    Output Parameter:
653: +  r1 - memory allocated in first chunk
654: .  r2 - memory allocated in second chunk
655: -  r3 - memory allocated in third chunk

657:    Level: developer

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

661:   Concepts: memory allocation

663: M*/
664: #if defined(PETSC_USE_DEBUG)
665: #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))
666: #else
667: #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)) \
668:                                                   || (*(r2) = (t2*)PetscAddrAlign(*(r1)+m1),*(r3) = (t3*)PetscAddrAlign(*(r2)+m2),0))
669: #endif

671: /*MC
672:    PetscMalloc4 - Allocates 4 chunks of  memory  all aligned to PETSC_MEMALIGN

674:    Synopsis:
675:    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)

677:    Not Collective

679:    Input Parameter:
680: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
681: .  t1 - type of first memory elements 
682: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
683: .  t2 - type of second memory elements
684: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
685: .  t3 - type of third memory elements
686: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
687: -  t4 - type of fourth memory elements

689:    Output Parameter:
690: +  r1 - memory allocated in first chunk
691: .  r2 - memory allocated in second chunk
692: .  r3 - memory allocated in third chunk
693: -  r4 - memory allocated in fourth chunk

695:    Level: developer

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

699:   Concepts: memory allocation

701: M*/
702: #if defined(PETSC_USE_DEBUG)
703: #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))
704: #else
705: #define PetscMalloc4(m1,t1,r1,m2,t2,r2,m3,t3,r3,m4,t4,r4)               \
706:   ((*(r2) = 0, *(r3) = 0, *(r4) = 0,PetscMalloc((m1)*sizeof(t1)+(m2)*sizeof(t2)+(m3)*sizeof(t3)+(m4)*sizeof(t4)+3*(PETSC_MEMALIGN-1),r1)) \
707:    || (*(r2) = (t2*)PetscAddrAlign(*(r1)+m1),*(r3) = (t3*)PetscAddrAlign(*(r2)+m2),*(r4) = (t4*)PetscAddrAlign(*(r3)+m3),0))
708: #endif

710: /*MC
711:    PetscMalloc5 - Allocates 5 chunks of  memory all aligned to PETSC_MEMALIGN

713:    Synopsis:
714:    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)

716:    Not Collective

718:    Input Parameter:
719: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
720: .  t1 - type of first memory elements 
721: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
722: .  t2 - type of second memory elements
723: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
724: .  t3 - type of third memory elements
725: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
726: .  t4 - type of fourth memory elements
727: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
728: -  t5 - type of fifth memory elements

730:    Output Parameter:
731: +  r1 - memory allocated in first chunk
732: .  r2 - memory allocated in second chunk
733: .  r3 - memory allocated in third chunk
734: .  r4 - memory allocated in fourth chunk
735: -  r5 - memory allocated in fifth chunk

737:    Level: developer

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

741:   Concepts: memory allocation

743: M*/
744: #if defined(PETSC_USE_DEBUG)
745: #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))
746: #else
747: #define PetscMalloc5(m1,t1,r1,m2,t2,r2,m3,t3,r3,m4,t4,r4,m5,t5,r5)      \
748:   ((*(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)) \
749:    || (*(r2) = (t2*)PetscAddrAlign(*(r1)+m1),*(r3) = (t3*)PetscAddrAlign(*(r2)+m2),*(r4) = (t4*)PetscAddrAlign(*(r3)+m3),*(r5) = (t5*)PetscAddrAlign(*(r4)+m4),0))
750: #endif


753: /*MC
754:    PetscMalloc6 - Allocates 6 chunks of  memory all aligned to PETSC_MEMALIGN

756:    Synopsis:
757:    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)

759:    Not Collective

761:    Input Parameter:
762: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
763: .  t1 - type of first memory elements 
764: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
765: .  t2 - type of second memory elements
766: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
767: .  t3 - type of third memory elements
768: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
769: .  t4 - type of fourth memory elements
770: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
771: .  t5 - type of fifth memory elements
772: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
773: -  t6 - type of sixth memory elements

775:    Output Parameter:
776: +  r1 - memory allocated in first chunk
777: .  r2 - memory allocated in second chunk
778: .  r3 - memory allocated in third chunk
779: .  r4 - memory allocated in fourth chunk
780: .  r5 - memory allocated in fifth chunk
781: -  r6 - memory allocated in sixth chunk

783:    Level: developer

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

787:   Concepts: memory allocation

789: M*/
790: #if defined(PETSC_USE_DEBUG)
791: #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))
792: #else
793: #define PetscMalloc6(m1,t1,r1,m2,t2,r2,m3,t3,r3,m4,t4,r4,m5,t5,r5,m6,t6,r6) \
794:   ((*(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)) \
795:    || (*(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))
796: #endif

798: /*MC
799:    PetscMalloc7 - Allocates 7 chunks of  memory all aligned to PETSC_MEMALIGN

801:    Synopsis:
802:    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)

804:    Not Collective

806:    Input Parameter:
807: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
808: .  t1 - type of first memory elements 
809: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
810: .  t2 - type of second memory elements
811: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
812: .  t3 - type of third memory elements
813: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
814: .  t4 - type of fourth memory elements
815: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
816: .  t5 - type of fifth memory elements
817: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
818: .  t6 - type of sixth memory elements
819: .  m7 - number of elements to allocate in 7th chunk  (may be zero)
820: -  t7 - type of sixth memory elements

822:    Output Parameter:
823: +  r1 - memory allocated in first chunk
824: .  r2 - memory allocated in second chunk
825: .  r3 - memory allocated in third chunk
826: .  r4 - memory allocated in fourth chunk
827: .  r5 - memory allocated in fifth chunk
828: .  r6 - memory allocated in sixth chunk
829: -  r7 - memory allocated in seventh chunk

831:    Level: developer

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

835:   Concepts: memory allocation

837: M*/
838: #if defined(PETSC_USE_DEBUG)
839: #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))
840: #else
841: #define PetscMalloc7(m1,t1,r1,m2,t2,r2,m3,t3,r3,m4,t4,r4,m5,t5,r5,m6,t6,r6,m7,t7,r7) \
842:   ((*(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)) \
843:    || (*(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))
844: #endif

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

849:    Synopsis:
850:    PetscErrorCode PetscNew(struct type,((type *))result)

852:    Not Collective

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

857:    Output Parameter:
858: .  result - memory allocated

860:    Level: beginner

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

864:   Concepts: memory allocation

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

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

873:    Synopsis:
874:    PetscErrorCode PetscNewLog(PetscObject obj,struct type,((type *))result)

876:    Not Collective

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

882:    Output Parameter:
883: .  result - memory allocated

885:    Level: developer

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

889:   Concepts: memory allocation

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

894: /*MC
895:    PetscFree - Frees memory

897:    Synopsis:
898:    PetscErrorCode PetscFree(void *memory)

900:    Not Collective

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

905:    Level: beginner

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

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

911:   Concepts: memory allocation

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

916: /*MC
917:    PetscFreeVoid - Frees memory

919:    Synopsis:
920:    void PetscFreeVoid(void *memory)

922:    Not Collective

924:    Input Parameter:
925: .   memory - memory to free

927:    Level: beginner

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

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

933:   Concepts: memory allocation

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


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

942:    Synopsis:
943:    PetscErrorCode PetscFree2(void *memory1,void *memory2)

945:    Not Collective

947:    Input Parameter:
948: +   memory1 - memory to free
949: -   memory2 - 2nd memory to free

951:    Level: developer

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

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

957:   Concepts: memory allocation

959: M*/
960: #if defined(PETSC_USE_DEBUG)
961: #define PetscFree2(m1,m2)   (PetscFree(m2) || PetscFree(m1))
962: #else
963: #define PetscFree2(m1,m2)   ((m2)=0, PetscFree(m1))
964: #endif

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

969:    Synopsis:
970:    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)

972:    Not Collective

974:    Input Parameter:
975: +   memory1 - memory to free
976: .   memory2 - 2nd memory to free
977: -   memory3 - 3rd memory to free

979:    Level: developer

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

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

985:   Concepts: memory allocation

987: M*/
988: #if defined(PETSC_USE_DEBUG)
989: #define PetscFree3(m1,m2,m3)   (PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
990: #else
991: #define PetscFree3(m1,m2,m3)   ((m3)=0,(m2)=0,PetscFree(m1))
992: #endif

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

997:    Synopsis:
998:    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)

1000:    Not Collective

1002:    Input Parameter:
1003: +   m1 - memory to free
1004: .   m2 - 2nd memory to free
1005: .   m3 - 3rd memory to free
1006: -   m4 - 4th memory to free

1008:    Level: developer

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

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

1014:   Concepts: memory allocation

1016: M*/
1017: #if defined(PETSC_USE_DEBUG)
1018: #define PetscFree4(m1,m2,m3,m4)   (PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1019: #else
1020: #define PetscFree4(m1,m2,m3,m4)   ((m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
1021: #endif

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

1026:    Synopsis:
1027:    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)

1029:    Not Collective

1031:    Input Parameter:
1032: +   m1 - memory to free
1033: .   m2 - 2nd memory to free
1034: .   m3 - 3rd memory to free
1035: .   m4 - 4th memory to free
1036: -   m5 - 5th memory to free

1038:    Level: developer

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

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

1044:   Concepts: memory allocation

1046: M*/
1047: #if defined(PETSC_USE_DEBUG)
1048: #define PetscFree5(m1,m2,m3,m4,m5)   (PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1049: #else
1050: #define PetscFree5(m1,m2,m3,m4,m5)   ((m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
1051: #endif


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

1057:    Synopsis:
1058:    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)

1060:    Not Collective

1062:    Input Parameter:
1063: +   m1 - memory to free
1064: .   m2 - 2nd memory to free
1065: .   m3 - 3rd memory to free
1066: .   m4 - 4th memory to free
1067: .   m5 - 5th memory to free
1068: -   m6 - 6th memory to free


1071:    Level: developer

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

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

1077:   Concepts: memory allocation

1079: M*/
1080: #if defined(PETSC_USE_DEBUG)
1081: #define PetscFree6(m1,m2,m3,m4,m5,m6)   (PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1082: #else
1083: #define PetscFree6(m1,m2,m3,m4,m5,m6)   ((m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
1084: #endif

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

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

1092:    Not Collective

1094:    Input Parameter:
1095: +   m1 - memory to free
1096: .   m2 - 2nd memory to free
1097: .   m3 - 3rd memory to free
1098: .   m4 - 4th memory to free
1099: .   m5 - 5th memory to free
1100: .   m6 - 6th memory to free
1101: -   m7 - 7th memory to free


1104:    Level: developer

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

1108: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1109:           PetscMalloc7()

1111:   Concepts: memory allocation

1113: M*/
1114: #if defined(PETSC_USE_DEBUG)
1115: #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   (PetscFree(m7) || PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1116: #else
1117: #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   ((m7)=0,(m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1))
1118: #endif

1120: extern  PetscErrorCode (*PetscTrMalloc)(size_t,int,const char[],const char[],const char[],void**);
1121: extern  PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[],const char[]);
1122: extern PetscErrorCode   PetscMallocSet(PetscErrorCode (*)(size_t,int,const char[],const char[],const char[],void**),PetscErrorCode (*)(void*,int,const char[],const char[],const char[]));
1123: extern PetscErrorCode   PetscMallocClear(void);

1125: /*
1126:    Routines for tracing memory corruption/bleeding with default PETSc  memory allocation
1127: */
1128: extern PetscErrorCode    PetscMallocDump(FILE *);
1129: extern PetscErrorCode    PetscMallocDumpLog(FILE *);
1130: extern PetscErrorCode    PetscMallocGetCurrentUsage(PetscLogDouble *);
1131: extern PetscErrorCode    PetscMallocGetMaximumUsage(PetscLogDouble *);
1132: extern PetscErrorCode    PetscMallocDebug(PetscBool);
1133: extern PetscErrorCode    PetscMallocValidate(int,const char[],const char[],const char[]);
1134: extern PetscErrorCode    PetscMallocSetDumpLog(void);


1137: /*E
1138:     PetscDataType - Used for handling different basic data types.

1140:    Level: beginner

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

1144: .seealso: PetscBinaryRead(), PetscBinaryWrite(), PetscDataTypeToMPIDataType(),
1145:           PetscDataTypeGetSize()

1147: E*/
1148: typedef enum {PETSC_INT = 0,PETSC_DOUBLE = 1,PETSC_COMPLEX = 2, PETSC_LONG = 3 ,PETSC_SHORT = 4,PETSC_FLOAT = 5,
1149:               PETSC_CHAR = 6,PETSC_BIT_LOGICAL = 7,PETSC_ENUM = 8,PETSC_BOOL=9, PETSC_LONG_DOUBLE = 10} PetscDataType;
1150: extern const char *PetscDataTypes[];

1152: #if defined(PETSC_USE_COMPLEX)
1153: #define  PETSC_SCALAR  PETSC_COMPLEX
1154: #else
1155: #if defined(PETSC_USE_REAL_SINGLE)
1156: #define  PETSC_SCALAR  PETSC_FLOAT
1157: #elif defined(PETSC_USE_REAL_LONG_DOUBLE)
1158: #define  PETSC_SCALAR  PETSC_LONG_DOUBLE
1159: #else
1160: #define  PETSC_SCALAR  PETSC_DOUBLE
1161: #endif
1162: #endif
1163: #if defined(PETSC_USE_REAL_SINGLE)
1164: #define  PETSC_REAL  PETSC_FLOAT
1165: #elif defined(PETSC_USE_REAL_LONG_DOUBLE)
1166: #define  PETSC_REAL  PETSC_LONG_DOUBLE
1167: #else
1168: #define  PETSC_REAL  PETSC_DOUBLE
1169: #endif
1170: #define  PETSC_FORTRANADDR  PETSC_LONG

1172: extern PetscErrorCode  PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1173: extern PetscErrorCode  PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1174: extern PetscErrorCode  PetscDataTypeGetSize(PetscDataType,size_t*);

1176: /*
1177:     Basic memory and string operations. These are usually simple wrappers
1178:    around the basic Unix system calls, but a few of them have additional
1179:    functionality and/or error checking.
1180: */
1181: extern PetscErrorCode    PetscBitMemcpy(void*,PetscInt,const void*,PetscInt,PetscInt,PetscDataType);
1182: extern PetscErrorCode    PetscMemmove(void*,void *,size_t);
1183: extern PetscErrorCode    PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1184: extern PetscErrorCode    PetscStrlen(const char[],size_t*);
1185: extern PetscErrorCode    PetscStrToArray(const char[],int*,char ***);
1186: extern PetscErrorCode    PetscStrToArrayDestroy(int,char **);
1187: extern PetscErrorCode    PetscStrcmp(const char[],const char[],PetscBool  *);
1188: extern PetscErrorCode    PetscStrgrt(const char[],const char[],PetscBool  *);
1189: extern PetscErrorCode    PetscStrcasecmp(const char[],const char[],PetscBool *);
1190: extern PetscErrorCode    PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1191: extern PetscErrorCode    PetscStrcpy(char[],const char[]);
1192: extern PetscErrorCode    PetscStrcat(char[],const char[]);
1193: extern PetscErrorCode    PetscStrncat(char[],const char[],size_t);
1194: extern PetscErrorCode    PetscStrncpy(char[],const char[],size_t);
1195: extern PetscErrorCode    PetscStrchr(const char[],char,char *[]);
1196: extern PetscErrorCode    PetscStrtolower(char[]);
1197: extern PetscErrorCode    PetscStrrchr(const char[],char,char *[]);
1198: extern PetscErrorCode    PetscStrstr(const char[],const char[],char *[]);
1199: extern PetscErrorCode    PetscStrrstr(const char[],const char[],char *[]);
1200: extern PetscErrorCode    PetscStrendswith(const char[],const char[],PetscBool*);
1201: extern PetscErrorCode    PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1202: extern PetscErrorCode    PetscStrallocpy(const char[],char *[]);
1203: extern PetscErrorCode    PetscStrreplace(MPI_Comm,const char[],char[],size_t);

1205: /*S
1206:     PetscToken - 'Token' used for managing tokenizing strings

1208:   Level: intermediate

1210: .seealso: PetscTokenCreate(), PetscTokenFind(), PetscTokenDestroy()
1211: S*/
1212: typedef struct _p_PetscToken* PetscToken;

1214: extern PetscErrorCode    PetscTokenCreate(const char[],const char,PetscToken*);
1215: extern PetscErrorCode    PetscTokenFind(PetscToken,char *[]);
1216: extern PetscErrorCode    PetscTokenDestroy(PetscToken);

1218: /*
1219:    These are  MPI operations for MPI_Allreduce() etc
1220: */
1221: extern  MPI_Op PetscMaxSum_Op;
1222: #if (defined(PETSC_USE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128)
1223: extern  MPI_Op MPIU_SUM;
1224: #else
1225: #define MPIU_SUM MPI_SUM
1226: #endif
1227: #if defined(PETSC_USE_REAL___FLOAT128)
1228: extern  MPI_Op MPIU_MAX;
1229: extern  MPI_Op MPIU_MIN;
1230: #else
1231: #define MPIU_MAX MPI_MAX
1232: #define MPIU_MIN MPI_MIN
1233: #endif
1234: extern PetscErrorCode  PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);

1236: extern PetscErrorCode MPILong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1237: extern PetscErrorCode MPILong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);

1239: /*S
1240:      PetscObject - any PETSc object, PetscViewer, Mat, Vec, KSP etc

1242:    Level: beginner

1244:    Note: This is the base class from which all objects appear.

1246: .seealso:  PetscObjectDestroy(), PetscObjectView(), PetscObjectGetName(), PetscObjectSetName(), PetscObjectReference(), PetscObjectDereferenc()
1247: S*/
1248: typedef struct _p_PetscObject* PetscObject;

1250: /*S
1251:      PetscFList - Linked list of functions, possibly stored in dynamic libraries, accessed
1252:       by string name

1254:    Level: advanced

1256: .seealso:  PetscFListAdd(), PetscFListDestroy()
1257: S*/
1258: typedef struct _n_PetscFList *PetscFList;

1260: /*E
1261:   PetscFileMode - Access mode for a file.

1263:   Level: beginner

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

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

1269:   FILE_MODE_APPEND - open a file at end for writing

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

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

1275: .seealso: PetscViewerFileSetMode()
1276: E*/
1277: typedef enum {FILE_MODE_READ, FILE_MODE_WRITE, FILE_MODE_APPEND, FILE_MODE_UPDATE, FILE_MODE_APPEND_UPDATE} PetscFileMode;

1279:  #include petscviewer.h
1280:  #include petscoptions.h

1282: #define PETSC_SMALLEST_CLASSID  1211211
1283: extern PetscClassId   PETSC_LARGEST_CLASSID;
1284: extern PetscClassId   PETSC_OBJECT_CLASSID;
1285: extern PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);

1287: /*
1288:    Routines that get memory usage information from the OS
1289: */
1290: extern PetscErrorCode  PetscMemoryGetCurrentUsage(PetscLogDouble *);
1291: extern PetscErrorCode  PetscMemoryGetMaximumUsage(PetscLogDouble *);
1292: extern PetscErrorCode  PetscMemorySetGetMaximumUsage(void);
1293: extern PetscErrorCode  PetscMemoryShowUsage(PetscViewer,const char[]);

1295: extern PetscErrorCode  PetscInfoAllow(PetscBool ,const char []);
1296: extern PetscErrorCode  PetscGetTime(PetscLogDouble*);
1297: extern PetscErrorCode  PetscGetCPUTime(PetscLogDouble*);
1298: extern PetscErrorCode  PetscSleep(PetscReal);

1300: /* 
1301:    Initialization of PETSc
1302: */
1303: extern PetscErrorCode  PetscInitialize(int*,char***,const char[],const char[]);
1304: PetscPolymorphicSubroutine(PetscInitialize,(int *argc,char ***args),(argc,args,PETSC_NULL,PETSC_NULL))
1305: extern PetscErrorCode  PetscInitializeNoArguments(void);
1306: extern PetscErrorCode  PetscInitialized(PetscBool  *);
1307: extern PetscErrorCode  PetscFinalized(PetscBool  *);
1308: extern PetscErrorCode  PetscFinalize(void);
1309: extern PetscErrorCode  PetscInitializeFortran(void);
1310: extern PetscErrorCode  PetscGetArgs(int*,char ***);
1311: extern PetscErrorCode  PetscGetArguments(char ***);
1312: extern PetscErrorCode  PetscFreeArguments(char **);

1314: extern PetscErrorCode  PetscEnd(void);
1315: extern PetscErrorCode  PetscSysInitializePackage(const char[]);

1317: extern MPI_Comm PETSC_COMM_LOCAL_WORLD;
1318: extern PetscErrorCode  PetscHMPIMerge(PetscMPIInt,PetscErrorCode (*)(void*),void*);
1319: extern PetscErrorCode  PetscHMPISpawn(PetscMPIInt);
1320: extern PetscErrorCode  PetscHMPIFinalize(void);
1321: extern PetscErrorCode  PetscHMPIRun(MPI_Comm,PetscErrorCode (*)(MPI_Comm,void *),void*);
1322: extern PetscErrorCode  PetscHMPIRunCtx(MPI_Comm,PetscErrorCode (*)(MPI_Comm,void*,void *),void*);
1323: extern PetscErrorCode  PetscHMPIFree(MPI_Comm,void*);
1324: extern PetscErrorCode  PetscHMPIMalloc(MPI_Comm,size_t,void**);

1326: extern PetscErrorCode  PetscPythonInitialize(const char[],const char[]);
1327: extern PetscErrorCode  PetscPythonFinalize(void);
1328: extern PetscErrorCode  PetscPythonPrintError(void);
1329: extern PetscErrorCode  PetscPythonMonitorSet(PetscObject,const char[]);

1331: /*
1332:      These are so that in extern C code we can caste function pointers to non-extern C
1333:    function pointers. Since the regular C++ code expects its function pointers to be 
1334:    C++.
1335: */
1336: typedef void (**PetscVoidStarFunction)(void);
1337: typedef void (*PetscVoidFunction)(void);
1338: typedef PetscErrorCode (*PetscErrorCodeFunction)(void);

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

1344:    Level: developer
1345:    
1346: .seealso: PetscUseMethod()
1347: */
1348: #define  PetscTryMethod(obj,A,B,C) \
1349:   0;{ PetscErrorCode (*f)B, __ierr; \
1350:     __PetscObjectQueryFunction((PetscObject)obj,A,(PetscVoidStarFunction)&f);CHKERRQ(__ierr); \
1351:     if (f) {__(*f)C;CHKERRQ(__ierr);}\
1352:   }

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

1358:    Level: developer
1359:    
1360: .seealso: PetscTryMethod()
1361: */
1362: #define  PetscUseMethod(obj,A,B,C) \
1363:   0;{ PetscErrorCode (*f)B, __ierr; \
1364:     __PetscObjectQueryFunction((PetscObject)obj,A,(PetscVoidStarFunction)&f);CHKERRQ(__ierr); \
1365:     if (f) {__(*f)C;CHKERRQ(__ierr);}\
1366:     else SETERRQ1(((PetscObject)obj)->comm,PETSC_ERR_SUP,"Cannot locate function %s in object",A); \
1367:   }

1369: /*
1370:     Functions that can act on any PETSc object.
1371: */
1372: extern PetscErrorCode  PetscObjectDestroy(PetscObject*);
1373: extern PetscErrorCode  PetscObjectGetComm(PetscObject,MPI_Comm *);
1374: extern PetscErrorCode  PetscObjectGetClassId(PetscObject,PetscClassId *);
1375: extern PetscErrorCode  PetscObjectGetClassName(PetscObject,const char *[]);
1376: extern PetscErrorCode  PetscObjectSetType(PetscObject,const char []);
1377: extern PetscErrorCode  PetscObjectSetPrecision(PetscObject,PetscPrecision);
1378: extern PetscErrorCode  PetscObjectGetType(PetscObject,const char *[]);
1379: extern PetscErrorCode  PetscObjectSetName(PetscObject,const char[]);
1380: extern PetscErrorCode  PetscObjectGetName(PetscObject,const char*[]);
1381: extern PetscErrorCode  PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer,const char[]);
1382: extern PetscErrorCode  PetscObjectSetTabLevel(PetscObject,PetscInt);
1383: extern PetscErrorCode  PetscObjectGetTabLevel(PetscObject,PetscInt*);
1384: extern PetscErrorCode  PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1385: extern PetscErrorCode  PetscObjectReference(PetscObject);
1386: extern PetscErrorCode  PetscObjectGetReference(PetscObject,PetscInt*);
1387: extern PetscErrorCode  PetscObjectDereference(PetscObject);
1388: extern PetscErrorCode  PetscObjectGetNewTag(PetscObject,PetscMPIInt *);
1389: extern PetscErrorCode  PetscObjectView(PetscObject,PetscViewer);
1390: extern PetscErrorCode  PetscObjectCompose(PetscObject,const char[],PetscObject);
1391: extern PetscErrorCode  PetscObjectRemoveReference(PetscObject,const char[]);
1392: extern PetscErrorCode  PetscObjectQuery(PetscObject,const char[],PetscObject *);
1393: extern PetscErrorCode  PetscObjectComposeFunction(PetscObject,const char[],const char[],void (*)(void));
1394: extern PetscErrorCode  PetscObjectSetFromOptions(PetscObject);
1395: extern PetscErrorCode  PetscObjectSetUp(PetscObject);
1396: extern PetscErrorCode  PetscCommGetNewTag(MPI_Comm,PetscMPIInt *);
1397: extern PetscErrorCode  PetscObjectAddOptionsHandler(PetscObject,PetscErrorCode (*)(PetscObject,void*),PetscErrorCode (*)(PetscObject,void*),void*);
1398: extern PetscErrorCode  PetscObjectProcessOptionsHandlers(PetscObject);
1399: extern PetscErrorCode  PetscObjectDestroyOptionsHandlers(PetscObject);

1401: /*MC
1402:    PetscObjectComposeFunctionDynamic - Associates a function with a given PETSc object. 
1403:                        
1404:     Synopsis:
1405:     PetscErrorCode PetscObjectComposeFunctionDynamic(PetscObject obj,const char name[],const char fname[],void *ptr)

1407:    Logically Collective on PetscObject

1409:    Input Parameters:
1410: +  obj - the PETSc object; this must be cast with a (PetscObject), for example, 
1411:          PetscObjectCompose((PetscObject)mat,...);
1412: .  name - name associated with the child function
1413: .  fname - name of the function
1414: -  ptr - function pointer (or PETSC_NULL if using dynamic libraries)

1416:    Level: advanced


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

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

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

1429:    Concepts: objects^composing functions
1430:    Concepts: composing functions
1431:    Concepts: functions^querying
1432:    Concepts: objects^querying
1433:    Concepts: querying objects

1435: .seealso: PetscObjectQueryFunction()
1436: M*/
1437: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1438: #define PetscObjectComposeFunctionDynamic(a,b,c,d) PetscObjectComposeFunction(a,b,c,0)
1439: #else
1440: #define PetscObjectComposeFunctionDynamic(a,b,c,d) PetscObjectComposeFunction(a,b,c,(PetscVoidFunction)(d))
1441: #endif

1443: extern PetscErrorCode  PetscObjectQueryFunction(PetscObject,const char[],void (**)(void));
1444: extern PetscErrorCode  PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1445: extern PetscErrorCode  PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1446: extern PetscErrorCode  PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1447: extern PetscErrorCode  PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1448: extern PetscErrorCode  PetscObjectAMSPublish(PetscObject);
1449: extern PetscErrorCode  PetscObjectUnPublish(PetscObject);
1450: extern PetscErrorCode  PetscObjectChangeTypeName(PetscObject,const char[]);
1451: extern PetscErrorCode  PetscObjectRegisterDestroy(PetscObject);
1452: extern PetscErrorCode  PetscObjectRegisterDestroyAll(void);
1453: extern PetscErrorCode  PetscObjectName(PetscObject);
1454: extern PetscErrorCode  PetscTypeCompare(PetscObject,const char[],PetscBool *);
1455: extern PetscErrorCode  PetscTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1456: extern PetscErrorCode  PetscRegisterFinalize(PetscErrorCode (*)(void));
1457: extern PetscErrorCode  PetscRegisterFinalizeAll(void);

1459: /*
1460:     Defines PETSc error handling.
1461: */
1462:  #include petscerror.h

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

1467:    Level: developer

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

1471: .seealso:  PetscOListAdd(), PetscOListDestroy(), PetscOListFind(), PetscObjectCompose(), PetscObjectQuery() 
1472: S*/
1473: typedef struct _n_PetscOList *PetscOList;

1475: extern PetscErrorCode  PetscOListDestroy(PetscOList*);
1476: extern PetscErrorCode  PetscOListFind(PetscOList,const char[],PetscObject*);
1477: extern PetscErrorCode  PetscOListReverseFind(PetscOList,PetscObject,char**,PetscBool*);
1478: extern PetscErrorCode  PetscOListAdd(PetscOList *,const char[],PetscObject);
1479: extern PetscErrorCode  PetscOListRemoveReference(PetscOList *,const char[]);
1480: extern PetscErrorCode  PetscOListDuplicate(PetscOList,PetscOList *);

1482: /*
1483:     Dynamic library lists. Lists of names of routines in objects or in dynamic 
1484:   link libraries that will be loaded as needed.
1485: */
1486: extern PetscErrorCode  PetscFListAdd(PetscFList*,const char[],const char[],void (*)(void));
1487: extern PetscErrorCode  PetscFListDestroy(PetscFList*);
1488: extern PetscErrorCode  PetscFListFind(PetscFList,MPI_Comm,const char[],PetscBool,void (**)(void));
1489: extern PetscErrorCode  PetscFListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFList,const char[]);
1490: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1491: #define    PetscFListAddDynamic(a,b,p,c) PetscFListAdd(a,b,p,0)
1492: #else
1493: #define    PetscFListAddDynamic(a,b,p,c) PetscFListAdd(a,b,p,(void (*)(void))c)
1494: #endif
1495: extern PetscErrorCode  PetscFListDuplicate(PetscFList,PetscFList *);
1496: extern PetscErrorCode  PetscFListView(PetscFList,PetscViewer);
1497: extern PetscErrorCode  PetscFListConcat(const char [],const char [],char []);
1498: extern PetscErrorCode  PetscFListGet(PetscFList,char ***,int*);

1500: /*S
1501:      PetscDLLibrary - Linked list of dynamics libraries to search for functions

1503:    Level: advanced

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

1507: .seealso:  PetscDLLibraryOpen()
1508: S*/
1509: typedef struct _n_PetscDLLibrary *PetscDLLibrary;
1510: extern  PetscDLLibrary DLLibrariesLoaded;
1511: extern PetscErrorCode  PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1512: extern PetscErrorCode  PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1513: extern PetscErrorCode  PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1514: extern PetscErrorCode  PetscDLLibraryPrintPath(PetscDLLibrary);
1515: extern PetscErrorCode  PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1516: extern PetscErrorCode  PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1517: extern PetscErrorCode  PetscDLLibraryClose(PetscDLLibrary);
1518: extern PetscErrorCode  PetscDLLibraryCCAAppend(MPI_Comm,PetscDLLibrary *,const char[]);

1520: /*
1521:   PetscFwk support.  Needs to be documented.  
1522:   Logically it is an extension of PetscDLLXXX, PetscObjectCompose, etc.
1523: */
1524:  #include petscfwk.h

1526: /*
1527:      Useful utility routines
1528: */
1529: extern PetscErrorCode  PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1530: extern PetscErrorCode  PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1531: extern PetscErrorCode  PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1532: PetscPolymorphicSubroutine(PetscSequentialPhaseBegin,(MPI_Comm comm),(comm,1))
1533: PetscPolymorphicSubroutine(PetscSequentialPhaseBegin,(void),(PETSC_COMM_WORLD,1))
1534: extern PetscErrorCode  PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1535: PetscPolymorphicSubroutine(PetscSequentialPhaseEnd,(MPI_Comm comm),(comm,1))
1536: PetscPolymorphicSubroutine(PetscSequentialPhaseEnd,(void),(PETSC_COMM_WORLD,1))
1537: extern PetscErrorCode  PetscBarrier(PetscObject);
1538: extern PetscErrorCode  PetscMPIDump(FILE*);

1540: /*
1541:     PetscNot - negates a logical type value and returns result as a PetscBool 

1543:     Notes: This is useful in cases like 
1544: $     int        *a;
1545: $     PetscBool  flag = PetscNot(a) 
1546:      where !a does not return a PetscBool  because we cannot provide a cast from int to PetscBool  in C.
1547: */
1548:  #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)

1550: /*
1551:     Defines basic graphics available from PETSc. 
1552: */
1553:  #include petscdraw.h

1555: /*
1556:     Defines the base data structures for all PETSc objects
1557: */
1558:  #include private/petscimpl.h


1561: /*MC
1562:     PetscErrorPrintf - Prints error messages.

1564:    Synopsis:
1565:      PetscErrorCode (*PetscErrorPrintf)(const char format[],...);

1567:     Not Collective

1569:     Input Parameters:
1570: .   format - the usual printf() format string 

1572:    Options Database Keys:
1573: +    -error_output_stdout - cause error messages to be printed to stdout instead of the
1574:          (default) stderr
1575: -    -error_output_none to turn off all printing of error messages (does not change the way the 
1576:           error is handled.)

1578:    Notes: Use
1579: $     PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the 
1580: $                        error is handled.) and
1581: $     PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on
1582: $        of you can use your own function

1584:           Use
1585:      PETSC_STDERR = FILE* obtained from a file open etc. to have stderr printed to the file. 
1586:      PETSC_STDOUT = FILE* obtained from a file open etc. to have stdout printed to the file. 

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

1591:    Level: developer

1593:     Fortran Note:
1594:     This routine is not supported in Fortran.

1596:     Concepts: error messages^printing
1597:     Concepts: printing^error messages

1599: .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscHelpPrintf(), PetscPrintf(), PetscErrorHandlerPush(), PetscVFPrintf(), PetscHelpPrintf()
1600: M*/
1601: extern  PetscErrorCode (*PetscErrorPrintf)(const char[],...);

1603: /*MC
1604:     PetscHelpPrintf - Prints help messages.

1606:    Synopsis:
1607:      PetscErrorCode (*PetscHelpPrintf)(const char format[],...);

1609:     Not Collective

1611:     Input Parameters:
1612: .   format - the usual printf() format string 

1614:    Level: developer

1616:     Fortran Note:
1617:     This routine is not supported in Fortran.

1619:     Concepts: help messages^printing
1620:     Concepts: printing^help messages

1622: .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1623: M*/
1624: extern  PetscErrorCode  (*PetscHelpPrintf)(MPI_Comm,const char[],...);

1626: /*
1627:      Defines PETSc profiling.
1628: */
1629:  #include petsclog.h

1631: /*
1632:           For locking, unlocking and destroying AMS memories associated with  PETSc objects. ams.h is included in petscviewer.h
1633: */
1634: #if defined(PETSC_HAVE_AMS)
1635: extern PetscBool  PetscAMSPublishAll;
1636: #define PetscObjectTakeAccess(obj)  ((((PetscObject)(obj))->amem == -1) ? 0 : AMS_Memory_take_access(((PetscObject)(obj))->amem))
1637: #define PetscObjectGrantAccess(obj) ((((PetscObject)(obj))->amem == -1) ? 0 : AMS_Memory_grant_access(((PetscObject)(obj))->amem))
1638: #define PetscObjectDepublish(obj)   ((((PetscObject)(obj))->amem == -1) ? 0 : AMS_Memory_destroy(((PetscObject)(obj))->amem));((PetscObject)(obj))->amem = -1;
1639: #else
1640: #define PetscObjectTakeAccess(obj)   0
1641: #define PetscObjectGrantAccess(obj)  0
1642: #define PetscObjectDepublish(obj)      0
1643: #endif

1645: /*
1646:       Simple PETSc parallel IO for ASCII printing
1647: */
1648: extern PetscErrorCode   PetscFixFilename(const char[],char[]);
1649: extern PetscErrorCode   PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1650: extern PetscErrorCode   PetscFClose(MPI_Comm,FILE*);
1651: extern PetscErrorCode   PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1652: extern PetscErrorCode   PetscPrintf(MPI_Comm,const char[],...);
1653: extern PetscErrorCode   PetscSNPrintf(char*,size_t,const char [],...);
1654: extern PetscErrorCode   PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);



1658: /* These are used internally by PETSc ASCII IO routines*/
1659: #include <stdarg.h>
1660: extern PetscErrorCode   PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list);
1661: extern PetscErrorCode   (*PetscVFPrintf)(FILE*,const char[],va_list);
1662: extern PetscErrorCode   PetscVFPrintfDefault(FILE*,const char[],va_list);

1664: #if defined(PETSC_HAVE_MATLAB_ENGINE)
1665: extern PetscErrorCode  PetscVFPrintf_Matlab(FILE*,const char[],va_list);
1666: #endif

1668: extern PetscErrorCode  PetscErrorPrintfDefault(const char [],...);
1669: extern PetscErrorCode  PetscErrorPrintfNone(const char [],...);
1670: extern PetscErrorCode  PetscHelpPrintfDefault(MPI_Comm,const char [],...);

1672: #if defined(PETSC_HAVE_POPEN)
1673: extern PetscErrorCode   PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1674: extern PetscErrorCode   PetscPClose(MPI_Comm,FILE*);
1675: #endif

1677: extern PetscErrorCode   PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1678: extern PetscErrorCode   PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1679: extern PetscErrorCode   PetscSynchronizedFlush(MPI_Comm);
1680: extern PetscErrorCode   PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1681: extern PetscErrorCode   PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1682: extern PetscErrorCode   PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1683: extern PetscErrorCode   PetscGetPetscDir(const char*[]);

1685: extern PetscErrorCode   PetscPopUpSelect(MPI_Comm,const char*,const char*,int,const char**,int*);

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

1690:    Level: advanced

1692: .seealso:  PetscObject, PetscContainerCreate()
1693: S*/
1694: extern PetscClassId  PETSC_CONTAINER_CLASSID;
1695: typedef struct _p_PetscContainer*  PetscContainer;
1696: extern PetscErrorCode  PetscContainerGetPointer(PetscContainer,void **);
1697: extern PetscErrorCode  PetscContainerSetPointer(PetscContainer,void *);
1698: extern PetscErrorCode  PetscContainerDestroy(PetscContainer*);
1699: extern PetscErrorCode  PetscContainerCreate(MPI_Comm,PetscContainer *);
1700: extern PetscErrorCode  PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void*));

1702: /*
1703:    For use in debuggers 
1704: */
1705: extern  PetscMPIInt PetscGlobalRank;
1706: extern  PetscMPIInt PetscGlobalSize;
1707: extern PetscErrorCode  PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1708: extern PetscErrorCode  PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1709: extern PetscErrorCode  PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);

1711: #if defined(PETSC_HAVE_MEMORY_H)
1712: #include <memory.h>
1713: #endif
1714: #if defined(PETSC_HAVE_STDLIB_H)
1715: #include <stdlib.h>
1716: #endif
1717: #if defined(PETSC_HAVE_STRINGS_H)
1718: #include <strings.h>
1719: #endif
1720: #if defined(PETSC_HAVE_STRING_H)
1721: #include <string.h>
1722: #endif


1725: #if defined(PETSC_HAVE_XMMINTRIN_H)
1726: #include <xmmintrin.h>
1727: #endif
1728: #if defined(PETSC_HAVE_STDINT_H)
1729: #include <stdint.h>
1730: #endif

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

1739:    Not Collective

1741:    Input Parameters:
1742: +  b - pointer to initial memory space
1743: -  n - length (in bytes) of space to copy

1745:    Output Parameter:
1746: .  a - pointer to copy space

1748:    Level: intermediate

1750:    Compile Option:
1751:     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used 
1752:                                   for memory copies on double precision values.
1753:     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used 
1754:                                   for memory copies on double precision values.
1755:     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used 
1756:                                   for memory copies on double precision values.

1758:    Note:
1759:    This routine is analogous to memcpy().

1761:    Developer Note: this is inlined for fastest performance

1763:   Concepts: memory^copying
1764:   Concepts: copying^memory
1765:   
1766: .seealso: PetscMemmove()

1768: @*/
1769: PETSC_STATIC_INLINE PetscErrorCode  PetscMemcpy(void *a,const void *b,size_t n)
1770: {
1771: #if defined(PETSC_USE_DEBUG)
1772:   unsigned long al = (unsigned long) a,bl = (unsigned long) b;
1773:   unsigned long nl = (unsigned long) n;
1775:   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1776:   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1777: #else
1779: #endif
1780:   if (a != b) {
1781: #if defined(PETSC_USE_DEBUG)
1782:     if ((al > bl && (al - bl) < nl) || (bl - al) < nl) {
1783:       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1784:               or make sure your copy regions and lengths are correct. \n\
1785:               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1786:     }
1787: #endif
1788: #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1789:    if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1790:       size_t len = n/sizeof(PetscScalar);
1791: #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1792:       PetscBLASInt one = 1,blen = PetscBLASIntCast(len);
1793:       BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one);
1794: #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1795:       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1796: #else
1797:       size_t      i;
1798:       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1799:       for (i=0; i<len; i++) y[i] = x[i];
1800: #endif
1801:     } else {
1802:       memcpy((char*)(a),(char*)(b),n);
1803:     }
1804: #else
1805:     memcpy((char*)(a),(char*)(b),n);
1806: #endif
1807:   }
1808:   return(0);
1809: }

1811: /*@C
1812:    PetscMemzero - Zeros the specified memory.

1814:    Not Collective

1816:    Input Parameters:
1817: +  a - pointer to beginning memory location
1818: -  n - length (in bytes) of memory to initialize

1820:    Level: intermediate

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

1826:    Developer Note: this is inlined for fastest performance

1828:    Concepts: memory^zeroing
1829:    Concepts: zeroing^memory

1831: .seealso: PetscMemcpy()
1832: @*/
1833: PETSC_STATIC_INLINE PetscErrorCode  PetscMemzero(void *a,size_t n)
1834: {
1835:   if (n > 0) {
1836: #if defined(PETSC_USE_DEBUG)
1837:     if (!a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer");
1838: #endif
1839: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1840:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1841:       size_t      i,len = n/sizeof(PetscScalar);
1842:       PetscScalar *x = (PetscScalar*)a;
1843:       for (i=0; i<len; i++) x[i] = 0.0;
1844:     } else {
1845: #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1846:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1847:       PetscInt len = n/sizeof(PetscScalar);
1848:       fortranzero_(&len,(PetscScalar*)a);
1849:     } else {
1850: #endif
1851: #if defined(PETSC_PREFER_BZERO)
1852:       bzero((char *)a,n);
1853: #else
1854:       memset((char*)a,0,n);
1855: #endif
1856: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1857:     }
1858: #endif
1859:   }
1860:   return 0;
1861: }

1863: /*MC
1864:    PetscPrefetchBlock - Prefetches a block of memory

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

1869:    Not Collective

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

1877:    Level: developer

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

1882:    Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality.  Not all architectures offer
1883:    equivalent locality hints, but the following macros are always defined to their closest analogue.
1884: +  PETSC_PREFETCH_HINT_NTA - Non-temporal.  Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched).
1885: .  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.
1886: .  PETSC_PREFETCH_HINT_T1 - Fetch to level 2 and higher (not L1).
1887: -  PETSC_PREFETCH_HINT_T2 - Fetch to high-level cache only.  (On many systems, T0 and T1 are equivalent.)

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

1892:    Concepts: memory
1893: M*/
1894: #define PetscPrefetchBlock(a,n,rw,t) do {                               \
1895:     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
1896:     for ( ; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
1897:   } while (0)

1899: /*
1900:     Allows accessing MATLAB Engine
1901: */
1902:  #include petscmatlab.h

1904: /*
1905:       Determine if some of the kernel computation routines use
1906:    Fortran (rather than C) for the numerical calculations. On some machines
1907:    and compilers (like complex numbers) the Fortran version of the routines
1908:    is faster than the C/C++ versions. The flag --with-fortran-kernels
1909:    should be used with ./configure to turn these on.
1910: */
1911: #if defined(PETSC_USE_FORTRAN_KERNELS)

1913: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1914: #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1915: #endif

1917: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
1918: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
1919: #endif

1921: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1922: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
1923: #endif

1925: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1926: #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
1927: #endif

1929: #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
1930: #define PETSC_USE_FORTRAN_KERNEL_NORM
1931: #endif

1933: #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
1934: #define PETSC_USE_FORTRAN_KERNEL_MAXPY
1935: #endif

1937: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1938: #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
1939: #endif

1941: #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1942: #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
1943: #endif

1945: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
1946: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
1947: #endif

1949: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1950: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
1951: #endif

1953: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
1954: #define PETSC_USE_FORTRAN_KERNEL_MDOT
1955: #endif

1957: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
1958: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
1959: #endif

1961: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
1962: #define PETSC_USE_FORTRAN_KERNEL_AYPX
1963: #endif

1965: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
1966: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
1967: #endif

1969: #endif

1971: /*
1972:     Macros for indicating code that should be compiled with a C interface,
1973:    rather than a C++ interface. Any routines that are dynamically loaded
1974:    (such as the PCCreate_XXX() routines) must be wrapped so that the name
1975:    mangler does not change the functions symbol name. This just hides the 
1976:    ugly extern "C" {} wrappers.
1977: */
1978: #if defined(__cplusplus)
1979: #define EXTERN_C_BEGIN extern "C" {
1980: #define EXTERN_C_END }
1981: #else
1982: #define EXTERN_C_BEGIN 
1983: #define EXTERN_C_END 
1984: #endif

1986: /* --------------------------------------------------------------------*/

1988: /*MC
1989:     MPI_Comm - the basic object used by MPI to determine which processes are involved in a 
1990:         communication

1992:    Level: beginner

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

1996: .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
1997: M*/

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


2005:    Level: beginner

2007: .seealso: PetscReal, PassiveReal, PassiveScalar, MPIU_SCALAR, PetscInt
2008: M*/

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

2013:    Level: beginner

2015: .seealso: PetscScalar, PassiveReal, PassiveScalar
2016: M*/

2018: /*MC
2019:     PassiveScalar - PETSc type that represents a PetscScalar
2020:    Level: beginner

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

2025: .seealso: PetscReal, PassiveReal, PetscScalar
2026: M*/

2028: /*MC
2029:     PassiveReal - PETSc type that represents a PetscReal

2031:    Level: beginner

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

2036: .seealso: PetscScalar, PetscReal, PassiveScalar
2037: M*/

2039: /*MC
2040:     MPIU_SCALAR - MPI datatype corresponding to PetscScalar

2042:    Level: beginner

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

2047: .seealso: PetscReal, PassiveReal, PassiveScalar, PetscScalar, MPIU_INT
2048: M*/

2050: #if defined(PETSC_HAVE_MPIIO)
2051: #if !defined(PETSC_WORDS_BIGENDIAN)
2052: extern PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2053: extern PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2054: #else
2055: #define MPIU_File_write_all(a,b,c,d,e) MPI_File_write_all(a,b,c,d,e) 
2056: #define MPIU_File_read_all(a,b,c,d,e) MPI_File_read_all(a,b,c,d,e) 
2057: #endif
2058: #endif

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

2062: /* Limit MPI to 32-bits */
2063: #define PETSC_MPI_INT_MAX  2147483647
2064: #define PETSC_MPI_INT_MIN -2147483647
2065: /* Limit BLAS to 32-bits */
2066: #define PETSC_BLAS_INT_MAX  2147483647
2067: #define PETSC_BLAS_INT_MIN -2147483647
2068: /* On 32 bit systems HDF5 is limited by size of integer, because hsize_t is defined as size_t */
2069: #define PETSC_HDF5_INT_MAX  2147483647
2070: #define PETSC_HDF5_INT_MIN -2147483647

2072: #if defined(PETSC_USE_64BIT_INDICES)
2073: #define PetscMPIIntCheck(a)  if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Message too long for MPI")
2074: #define PetscBLASIntCheck(a)  if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK")
2075: #define PetscMPIIntCast(a) (PetscMPIInt)(a);PetscMPIIntCheck(a)
2076: #define PetscBLASIntCast(a) (PetscBLASInt)(a);PetscBLASIntCheck(a)

2078: #if (PETSC_SIZEOF_SIZE_T == 4)
2079: #define PetscHDF5IntCheck(a)  if ((a) > PETSC_HDF5_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for HDF5")
2080: #define PetscHDF5IntCast(a) (hsize_t)(a);PetscHDF5IntCheck(a)
2081: #else
2082: #define PetscHDF5IntCheck(a)
2083: #define PetscHDF5IntCast(a) a
2084: #endif

2086: #else
2087: #define PetscMPIIntCheck(a) 
2088: #define PetscBLASIntCheck(a) 
2089: #define PetscHDF5IntCheck(a)
2090: #define PetscMPIIntCast(a) a
2091: #define PetscBLASIntCast(a) a
2092: #define PetscHDF5IntCast(a) a
2093: #endif  


2096: /*
2097:      The IBM include files define hz, here we hide it so that it may be used
2098:    as a regular user variable.
2099: */
2100: #if defined(hz)
2101: #undef hz
2102: #endif

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


2107: #if defined(PETSC_HAVE_LIMITS_H)
2108: #include <limits.h>
2109: #endif
2110: #if defined(PETSC_HAVE_SYS_PARAM_H)
2111: #include <sys/param.h>
2112: #endif
2113: #if defined(PETSC_HAVE_SYS_TYPES_H)
2114: #include <sys/types.h>
2115: #endif
2116: #if defined(MAXPATHLEN)
2117: #  define PETSC_MAX_PATH_LEN     MAXPATHLEN
2118: #elif defined(MAX_PATH)
2119: #  define PETSC_MAX_PATH_LEN     MAX_PATH
2120: #elif defined(_MAX_PATH)
2121: #  define PETSC_MAX_PATH_LEN     _MAX_PATH
2122: #else
2123: #  define PETSC_MAX_PATH_LEN     4096
2124: #endif

2126: /* Special support for C++ */
2127:  #include petscsys.hh


2130: /*MC

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

2134: $    1) classic Fortran 77 style
2135: $#include "finclude/petscXXX.h" to work with material from the XXX component of PETSc
2136: $       XXX variablename
2137: $      You cannot use this approach if you wish to use the Fortran 90 specific PETSc routines
2138: $      which end in F90; such as VecGetArrayF90()
2139: $
2140: $    2) classic Fortran 90 style
2141: $#include "finclude/petscXXX.h" 
2142: $#include "finclude/petscXXX.h90" to work with material from the XXX component of PETSc
2143: $       XXX variablename
2144: $
2145: $    3) Using Fortran modules
2146: $#include "finclude/petscXXXdef.h" 
2147: $         use petscXXXX
2148: $       XXX variablename
2149: $
2150: $    4) Use Fortran modules and Fortran data types for PETSc types
2151: $#include "finclude/petscXXXdef.h" 
2152: $         use petscXXXX
2153: $       type(XXX) variablename
2154: $      To use this approach you must ./configure PETSc with the additional
2155: $      option --with-fortran-datatypes You cannot use the type(XXX) declaration approach without using Fortran modules

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

2159: $    3a) skip the #include BUT you cannot use any PETSc data type names like Vec, Mat, PetscInt, PetscErrorCode etc
2160: $        and you must declare the variables as integer, for example 
2161: $        integer variablename
2162: $
2163: $    4a) skip the #include, you use the object types like type(Vec) type(Mat) but cannot use the data type
2164: $        names like PetscErrorCode, PetscInt etc. again for those you must use integer

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

2169:    Fortran type checking with interfaces is strick, this means you cannot pass a scalar value when an array value
2170: is expected (even though it is legal Fortran). For example when setting a single value in a matrix with MatSetValues()
2171: you cannot have something like
2172: $      PetscInt row,col
2173: $      PetscScalar val
2174: $        ...
2175: $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)
2176: You must instead have 
2177: $      PetscInt row(1),col(1)
2178: $      PetscScalar val(1)
2179: $        ...
2180: $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)


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

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

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

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

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

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

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

2205:     Level: beginner

2207: M*/

2209: extern PetscErrorCode  PetscGetArchType(char[],size_t);
2210: extern PetscErrorCode  PetscGetHostName(char[],size_t);
2211: extern PetscErrorCode  PetscGetUserName(char[],size_t);
2212: extern PetscErrorCode  PetscGetProgramName(char[],size_t);
2213: extern PetscErrorCode  PetscSetProgramName(const char[]);
2214: extern PetscErrorCode  PetscGetDate(char[],size_t);

2216: extern PetscErrorCode  PetscSortInt(PetscInt,PetscInt[]);
2217: extern PetscErrorCode  PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2218: extern PetscErrorCode  PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2219: extern PetscErrorCode  PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2220: extern PetscErrorCode  PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2221: extern PetscErrorCode  PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2222: extern PetscErrorCode  PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2223: extern PetscErrorCode  PetscSortReal(PetscInt,PetscReal[]);
2224: extern PetscErrorCode  PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2225: extern PetscErrorCode  PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2226: extern PetscErrorCode  PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2227: extern PetscErrorCode  PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);

2229: extern PetscErrorCode  PetscSetDisplay(void);
2230: extern PetscErrorCode  PetscGetDisplay(char[],size_t);

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

2237:    Level: beginner

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

2242: .seealso: PetscRandomSetType(), PetscRandom
2243: J*/
2244: #define PetscRandomType char*
2245: #define PETSCRAND       "rand"
2246: #define PETSCRAND48     "rand48"
2247: #define PETSCSPRNG      "sprng"          

2249: /* Logging support */
2250: extern  PetscClassId PETSC_RANDOM_CLASSID;

2252: extern PetscErrorCode  PetscRandomInitializePackage(const char[]);

2254: /*S
2255:      PetscRandom - Abstract PETSc object that manages generating random numbers

2257:    Level: intermediate

2259:   Concepts: random numbers

2261: .seealso:  PetscRandomCreate(), PetscRandomGetValue(), PetscRandomType
2262: S*/
2263: typedef struct _p_PetscRandom*   PetscRandom;

2265: /* Dynamic creation and loading functions */
2266: extern PetscFList PetscRandomList;
2267: extern PetscBool  PetscRandomRegisterAllCalled;

2269: extern PetscErrorCode  PetscRandomRegisterAll(const char []);
2270: extern PetscErrorCode  PetscRandomRegister(const char[],const char[],const char[],PetscErrorCode (*)(PetscRandom));
2271: extern PetscErrorCode  PetscRandomRegisterDestroy(void);
2272: extern PetscErrorCode  PetscRandomSetType(PetscRandom, const PetscRandomType);
2273: extern PetscErrorCode  PetscRandomSetFromOptions(PetscRandom);
2274: extern PetscErrorCode  PetscRandomGetType(PetscRandom, const PetscRandomType*);
2275: extern PetscErrorCode  PetscRandomViewFromOptions(PetscRandom,char*);
2276: extern PetscErrorCode  PetscRandomView(PetscRandom,PetscViewer);

2278: /*MC
2279:   PetscRandomRegisterDynamic - Adds a new PetscRandom component implementation

2281:   Synopsis:
2282:   PetscErrorCode PetscRandomRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(PetscRandom))

2284:   Not Collective

2286:   Input Parameters:
2287: + name        - The name of a new user-defined creation routine
2288: . path        - The path (either absolute or relative) of the library containing this routine
2289: . func_name   - The name of routine to create method context
2290: - create_func - The creation routine itself

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

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

2297:   Sample usage:
2298: .vb
2299:     PetscRandomRegisterDynamic("my_rand","/home/username/my_lib/lib/libO/solaris/libmy.a", "MyPetscRandomtorCreate", MyPetscRandomtorCreate);
2300: .ve

2302:   Then, your random type can be chosen with the procedural interface via
2303: .vb
2304:     PetscRandomCreate(MPI_Comm, PetscRandom *);
2305:     PetscRandomSetType(PetscRandom,"my_random_name");
2306: .ve
2307:    or at runtime via the option
2308: .vb
2309:     -random_type my_random_name
2310: .ve

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

2314:          For an example of the code needed to interface your own random number generator see
2315:          src/sys/random/impls/rand/rand.c
2316:         
2317:   Level: advanced

2319: .keywords: PetscRandom, register
2320: .seealso: PetscRandomRegisterAll(), PetscRandomRegisterDestroy(), PetscRandomRegister()
2321: M*/
2322: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
2323: #define PetscRandomRegisterDynamic(a,b,c,d) PetscRandomRegister(a,b,c,0)
2324: #else
2325: #define PetscRandomRegisterDynamic(a,b,c,d) PetscRandomRegister(a,b,c,d)
2326: #endif

2328: extern PetscErrorCode  PetscRandomCreate(MPI_Comm,PetscRandom*);
2329: extern PetscErrorCode  PetscRandomGetValue(PetscRandom,PetscScalar*);
2330: extern PetscErrorCode  PetscRandomGetValueReal(PetscRandom,PetscReal*);
2331: extern PetscErrorCode  PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2332: extern PetscErrorCode  PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2333: extern PetscErrorCode  PetscRandomSetSeed(PetscRandom,unsigned long);
2334: extern PetscErrorCode  PetscRandomGetSeed(PetscRandom,unsigned long *);
2335: extern PetscErrorCode  PetscRandomSeed(PetscRandom);
2336: extern PetscErrorCode  PetscRandomDestroy(PetscRandom*);

2338: extern PetscErrorCode  PetscGetFullPath(const char[],char[],size_t);
2339: extern PetscErrorCode  PetscGetRelativePath(const char[],char[],size_t);
2340: extern PetscErrorCode  PetscGetWorkingDirectory(char[],size_t);
2341: extern PetscErrorCode  PetscGetRealPath(const char[],char[]);
2342: extern PetscErrorCode  PetscGetHomeDirectory(char[],size_t);
2343: extern PetscErrorCode  PetscTestFile(const char[],char,PetscBool *);
2344: extern PetscErrorCode  PetscTestDirectory(const char[],char,PetscBool *);

2346: extern PetscErrorCode  PetscBinaryRead(int,void*,PetscInt,PetscDataType);
2347: extern PetscErrorCode  PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscDataType);
2348: extern PetscErrorCode  PetscBinarySynchronizedWrite(MPI_Comm,int,void*,PetscInt,PetscDataType,PetscBool );
2349: extern PetscErrorCode  PetscBinaryWrite(int,void*,PetscInt,PetscDataType,PetscBool );
2350: extern PetscErrorCode  PetscBinaryOpen(const char[],PetscFileMode,int *);
2351: extern PetscErrorCode  PetscBinaryClose(int);
2352: extern PetscErrorCode  PetscSharedTmp(MPI_Comm,PetscBool  *);
2353: extern PetscErrorCode  PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2354: extern PetscErrorCode  PetscGetTmp(MPI_Comm,char[],size_t);
2355: extern PetscErrorCode  PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2356: extern PetscErrorCode  PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2357: extern PetscErrorCode  PetscOpenSocket(char*,int,int*);
2358: extern PetscErrorCode  PetscWebServe(MPI_Comm,int);

2360: /*
2361:    In binary files variables are stored using the following lengths,
2362:   regardless of how they are stored in memory on any one particular
2363:   machine. Use these rather then sizeof() in computing sizes for 
2364:   PetscBinarySeek().
2365: */
2366: #define PETSC_BINARY_INT_SIZE   (32/8)
2367: #define PETSC_BINARY_FLOAT_SIZE  (32/8)
2368: #define PETSC_BINARY_CHAR_SIZE  (8/8)
2369: #define PETSC_BINARY_SHORT_SIZE  (16/8)
2370: #define PETSC_BINARY_DOUBLE_SIZE  (64/8)
2371: #define PETSC_BINARY_SCALAR_SIZE  sizeof(PetscScalar)

2373: /*E
2374:   PetscBinarySeekType - argument to PetscBinarySeek()

2376:   Level: advanced

2378: .seealso: PetscBinarySeek(), PetscBinarySynchronizedSeek()
2379: E*/
2380: typedef enum {PETSC_BINARY_SEEK_SET = 0,PETSC_BINARY_SEEK_CUR = 1,PETSC_BINARY_SEEK_END = 2} PetscBinarySeekType;
2381: extern PetscErrorCode  PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2382: extern PetscErrorCode  PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2383: extern PetscErrorCode  PetscByteSwap(void *,PetscDataType,PetscInt);

2385: extern PetscErrorCode  PetscSetDebugTerminal(const char[]);
2386: extern PetscErrorCode  PetscSetDebugger(const char[],PetscBool );
2387: extern PetscErrorCode  PetscSetDefaultDebugger(void);
2388: extern PetscErrorCode  PetscSetDebuggerFromString(char*);
2389: extern PetscErrorCode  PetscAttachDebugger(void);
2390: extern PetscErrorCode  PetscStopForDebugger(void);

2392: extern PetscErrorCode  PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2393: extern PetscErrorCode  PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2394: extern PetscErrorCode  PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2395: extern PetscErrorCode  PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2396: extern PetscErrorCode  PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);

2398: extern PetscErrorCode  PetscSSEIsEnabled(MPI_Comm,PetscBool  *,PetscBool  *);

2400: /*E
2401:   InsertMode - Whether entries are inserted or added into vectors or matrices

2403:   Level: beginner

2405: .seealso: VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2406:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(),
2407:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd()
2408: E*/
2409:  typedef enum {NOT_SET_VALUES, INSERT_VALUES, ADD_VALUES, MAX_VALUES, INSERT_ALL_VALUES, ADD_ALL_VALUES} InsertMode;

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

2414:     Level: beginner

2416: .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2417:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), ADD_VALUES,
2418:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES

2420: M*/

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

2426:     Level: beginner

2428: .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2429:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), INSERT_VALUES,
2430:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES

2432: M*/

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

2437:     Level: beginner

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

2441: M*/

2443: /*S
2444:    PetscSubcomm - Context of MPI subcommunicators, used by PCREDUNDANT

2446:    Level: advanced

2448:    Concepts: communicator, create
2449: S*/
2450: typedef struct _n_PetscSubcomm* PetscSubcomm;

2452: struct _n_PetscSubcomm {
2453:   MPI_Comm   parent;      /* parent communicator */
2454:   MPI_Comm   dupparent;   /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2455:   MPI_Comm   comm;        /* this communicator */
2456:   PetscInt   n;           /* num of subcommunicators under the parent communicator */
2457:   PetscInt   color;       /* color of processors belong to this communicator */
2458: };

2460: typedef enum {PETSC_SUBCOMM_GENERAL=0,PETSC_SUBCOMM_CONTIGUOUS=1,PETSC_SUBCOMM_INTERLACED=2} PetscSubcommType;
2461: extern const char *PetscSubcommTypes[];

2463: extern PetscErrorCode  PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2464: extern PetscErrorCode  PetscSubcommDestroy(PetscSubcomm*);
2465: extern PetscErrorCode  PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2466: extern PetscErrorCode  PetscSubcommSetType(PetscSubcomm,const PetscSubcommType);
2467: extern PetscErrorCode  PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt,PetscMPIInt);

2469: PETSC_EXTERN_CXX_END

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

2475: #endif