Actual source code: petscsys.h

petsc-dev 2014-04-20
Report Typos and Errors
  1: /*
  2:    This is the main PETSc include file (for C and C++).  It is included by all
  3:    other PETSc include files, so it almost never has to be specifically included.
  4: */
  7: /* ========================================================================== */
  8: /*
  9:    petscconf.h is contained in ${PETSC_ARCH}/include/petscconf.h it is
 10:    found automatically by the compiler due to the -I${PETSC_DIR}/${PETSC_ARCH}/include
 11:    in the conf/variables definition of PETSC_INCLUDE. For --prefix installs the ${PETSC_ARCH}/
 12:    does not exist and petscconf.h is in the same directory as the other PETSc include files.
 13: */
 14: #include <petscconf.h>
 15: #include <petscfix.h>

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

 33: /* ========================================================================== */
 34: /*
 35:    This facilitates using the C version of PETSc from C++ and the C++ version from C.
 36: */
 37: #if defined(__cplusplus)
 38: #  define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_CXX
 39: #else
 40: #  define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_C
 41: #endif

 43: #if defined(__cplusplus)
 44: #  define PETSC_RESTRICT PETSC_CXX_RESTRICT
 45: #else
 46: #  define PETSC_RESTRICT PETSC_C_RESTRICT
 47: #endif

 49: #if defined(__cplusplus)
 50: #  define PETSC_STATIC_INLINE PETSC_CXX_STATIC_INLINE
 51: #else
 52: #  define PETSC_STATIC_INLINE PETSC_C_STATIC_INLINE
 53: #endif

 55: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES) /* For Win32 shared libraries */
 56: #  define  __declspec(dllexport)
 57: #  define PETSC_DLLIMPORT __declspec(dllimport)
 58: #  define PETSC_VISIBILITY_INTERNAL
 59: #elif defined(PETSC_USE_VISIBILITY)
 60: #  define  __attribute__((visibility ("default")))
 61: #  define PETSC_DLLIMPORT __attribute__((visibility ("default")))
 62: #  define PETSC_VISIBILITY_INTERNAL __attribute__((visibility ("hidden")))
 63: #else
 64: #  define 
 65: #  define PETSC_DLLIMPORT
 66: #  define PETSC_VISIBILITY_INTERNAL
 67: #endif

 69: #if defined(petsc_EXPORTS)      /* CMake defines this when building the shared library */
 70: #  define PETSC_VISIBILITY_PUBLIC 
 71: #else  /* Win32 users need this to import symbols from petsc.dll */
 72: #  define PETSC_VISIBILITY_PUBLIC PETSC_DLLIMPORT
 73: #endif

 75: #if defined(__cplusplus)
 76: #define PETSC_EXTERN extern "C" PETSC_VISIBILITY_PUBLIC
 77: #define PETSC_EXTERN_TYPEDEF extern "C"
 78: #define PETSC_INTERN extern "C" PETSC_VISIBILITY_INTERNAL
 79: #else
 80: #define PETSC_EXTERN extern PETSC_VISIBILITY_PUBLIC
 81: #define PETSC_EXTERN_TYPEDEF
 82: #define PETSC_INTERN extern PETSC_VISIBILITY_INTERNAL
 83: #endif

 85: #include <petscversion.h>
 86: #define PETSC_AUTHOR_INFO  "       The PETSc Team\n    petsc-maint@mcs.anl.gov\n http://www.mcs.anl.gov/petsc/\n"

 88: /* ========================================================================== */

 90: /*
 91:     Defines the interface to MPI allowing the use of all MPI functions.

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

107: /*
108:     Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler
109:     see the top of mpicxx.h in the MPICH2 distribution.
110: */
111: #include <stdio.h>

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

118: /* Support for Clang (>=3.2) matching type tag arguments with void* buffer types */
119: #if defined(__has_attribute)
120: #  if __has_attribute(argument_with_type_tag) && __has_attribute(pointer_with_type_tag) && __has_attribute(type_tag_for_datatype)
121: #    define PetscAttrMPIPointerWithType(bufno,typeno) __attribute__((pointer_with_type_tag(MPI,bufno,typeno)))
122: #    define PetscAttrMPITypeTag(type)                 __attribute__((type_tag_for_datatype(MPI,type)))
123: #    define PetscAttrMPITypeTagLayoutCompatible(type) __attribute__((type_tag_for_datatype(MPI,type,layout_compatible)))
124: #  endif
125: #endif
126: #if !defined(PetscAttrMPIPointerWithType)
127: #  define PetscAttrMPIPointerWithType(bufno,typeno)
128: #  define PetscAttrMPITypeTag(type)
129: #  define PetscAttrMPITypeTagLayoutCompatible(type)
130: #endif

132: /*MC
133:     PetscErrorCode - datatype used for return error code from almost all PETSc functions

135:     Level: beginner

137: .seealso: CHKERRQ, SETERRQ
138: M*/
139: typedef int PetscErrorCode;

141: /*MC

143:     PetscClassId - A unique id used to identify each PETSc class.

145:     Notes: Use PetscClassIdRegister() to obtain a new value for a new class being created. Usually
146:          XXXInitializePackage() calls it for each class it defines.

148:     Developer Notes: Internal integer stored in the _p_PetscObject data structure.
149:          These are all computed by an offset from the lowest one, PETSC_SMALLEST_CLASSID.

151:     Level: developer

153: .seealso: PetscClassIdRegister(), PetscLogEventRegister(), PetscHeaderCreate()
154: M*/
155: typedef int PetscClassId;


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

161:     Level: intermediate

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

166:     PetscMPIIntCast(a,&b) checks if the given PetscInt a will fit in a PetscMPIInt, if not it
167:       generates a PETSC_ERR_ARG_OUTOFRANGE error.

169: .seealso: PetscBLASInt, PetscInt

171: M*/
172: typedef int PetscMPIInt;

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

177:     Level: intermediate

179: .seealso: PetscOptionsGetEnum(), PetscOptionsEnum(), PetscBagRegisterEnum()
180: M*/
181: typedef enum { ENUM_DUMMY } PetscEnum;
182: PETSC_EXTERN MPI_Datatype MPIU_ENUM PetscAttrMPITypeTag(PetscEnum);

184: #if defined(PETSC_HAVE_STDINT_H)
185: #include <stdint.h>
186: #endif

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

193:    Level: intermediate

195: .seealso: PetscScalar, PetscBLASInt, PetscMPIInt
196: M*/
197: #if defined(PETSC_HAVE_STDINT_H)
198: typedef int64_t Petsc64bitInt;
199: #elif (PETSC_SIZEOF_LONG_LONG == 8)
200: typedef long long Petsc64bitInt;
201: #elif defined(PETSC_HAVE___INT64)
202: typedef __int64 Petsc64bitInt;
203: #else
204: typedef unknown64bit Petsc64bitInt
205: #endif
206: #if defined(PETSC_USE_64BIT_INDICES)
207: typedef Petsc64bitInt PetscInt;
208: #  if defined(PETSC_HAVE_MPI_INT64_T) /* MPI_INT64_T is not guaranteed to be a macro */
209: #    define MPIU_INT MPI_INT64_T
210: #  else
211: #    define MPIU_INT MPI_LONG_LONG_INT
212: #  endif
213: #else
214: typedef int PetscInt;
215: #define MPIU_INT MPI_INT
216: #endif
217: #if defined(PETSC_HAVE_MPI_INT64_T)
218: #  define MPIU_INT64 MPI_INT64_T
219: #else
220: #  define MPIU_INT64 MPI_LONG_LONG_INT
221: #endif


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

227:     Level: intermediate

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

233:     PetscErrorCode PetscBLASIntCast(a,&b) checks if the given PetscInt a will fit in a PetscBLASInt, if not it
234:       generates a PETSC_ERR_ARG_OUTOFRANGE error

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

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

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

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

248: .seealso: PetscMPIInt, PetscInt

250: M*/
251: #if defined(PETSC_HAVE_64BIT_BLAS_INDICES)
252: typedef Petsc64bitInt PetscBLASInt;
253: #else
254: typedef int PetscBLASInt;
255: #endif

257: /*EC

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

261:     Level: advanced

263: .seealso: PetscObjectSetPrecision()
264: E*/
265: typedef enum { PETSC_PRECISION_SINGLE=4,PETSC_PRECISION_DOUBLE=8 } PetscPrecision;
266: PETSC_EXTERN const char *PetscPrecisions[];

268: /*
269:     For the rare cases when one needs to send a size_t object with MPI
270: */
271: #if (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_INT)
272: #define MPIU_SIZE_T MPI_UNSIGNED
273: #elif  (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_LONG)
274: #define MPIU_SIZE_T MPI_UNSIGNED_LONG
275: #elif  (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_LONG_LONG)
276: #define MPIU_SIZE_T MPI_UNSIGNED_LONG_LONG
277: #else
278: #error "Unknown size for size_t! Send us a bugreport at petsc-maint@mcs.anl.gov"
279: #endif


282: /*
283:       You can use PETSC_STDOUT as a replacement of stdout. You can also change
284:     the value of PETSC_STDOUT to redirect all standard output elsewhere
285: */
286: PETSC_EXTERN FILE* PETSC_STDOUT;

288: /*
289:       You can use PETSC_STDERR as a replacement of stderr. You can also change
290:     the value of PETSC_STDERR to redirect all standard error elsewhere
291: */
292: PETSC_EXTERN FILE* PETSC_STDERR;

294: /*MC
295:     PetscUnlikely - hints the compiler that the given condition is usually FALSE

297:     Synopsis:
298:     #include <petscsys.h>
299:     PetscBool  PetscUnlikely(PetscBool  cond)

301:     Not Collective

303:     Input Parameters:
304: .   cond - condition or expression

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

309:     Level: advanced

311: .seealso: PetscLikely(), CHKERRQ
312: M*/

314: /*MC
315:     PetscLikely - hints the compiler that the given condition is usually TRUE

317:     Synopsis:
318:     #include <petscsys.h>
319:     PetscBool  PetscUnlikely(PetscBool  cond)

321:     Not Collective

323:     Input Parameters:
324: .   cond - condition or expression

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

329:     Level: advanced

331: .seealso: PetscUnlikely()
332: M*/
333: #if defined(PETSC_HAVE_BUILTIN_EXPECT)
334: #  define PetscUnlikely(cond)   __builtin_expect(!!(cond),0)
335: #  define PetscLikely(cond)     __builtin_expect(!!(cond),1)
336: #else
337: #  define PetscUnlikely(cond)   (cond)
338: #  define PetscLikely(cond)     (cond)
339: #endif

341: /*
342:     Declare extern C stuff after including external header files
343: */


346: /*
347:        Basic PETSc constants
348: */

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

353:    Level: beginner

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

358: E*/
359: typedef enum { PETSC_FALSE,PETSC_TRUE } PetscBool;
360: PETSC_EXTERN const char *const PetscBools[];
361: PETSC_EXTERN MPI_Datatype MPIU_BOOL PetscAttrMPITypeTag(PetscBool);

363: /*
364:     Defines some elementary mathematics functions and constants.
365: */
366: #include <petscmath.h>

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

371:    Level: beginner

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

379: E*/
380: typedef enum { PETSC_COPY_VALUES, PETSC_OWN_POINTER, PETSC_USE_POINTER} PetscCopyMode;
381: PETSC_EXTERN const char *const PetscCopyModes[];

383: /*MC
384:     PETSC_FALSE - False value of PetscBool

386:     Level: beginner

388:     Note: Zero integer

390: .seealso: PetscBool , PETSC_TRUE
391: M*/

393: /*MC
394:     PETSC_TRUE - True value of PetscBool

396:     Level: beginner

398:     Note: Nonzero integer

400: .seealso: PetscBool , PETSC_FALSE
401: M*/

403: /*MC
404:     PETSC_NULL - standard way of passing in a null or array or pointer. This is deprecated in C/C++ simply use NULL

406:    Level: beginner

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

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

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

416: M*/
417: #define PETSC_NULL           NULL

419: /*MC
420:     PETSC_IGNORE - same as NULL, means PETSc will ignore this argument

422:    Level: beginner

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

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

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

432: M*/
433: #define PETSC_IGNORE         NULL

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

439:    Level: beginner

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

443: M*/
444: #define PETSC_DECIDE  -1

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

450:    Level: beginner


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

456: .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_IGNORE, VecSetSizes()

458: M*/
459: #define PETSC_DETERMINE PETSC_DECIDE

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

465:    Level: beginner

467:    Fortran Notes: You need to use PETSC_DEFAULT_INTEGER or PETSC_DEFAULT_REAL.

469: .seealso: PETSC_DECIDE, PETSC_IGNORE, PETSC_DETERMINE

471: M*/
472: #define PETSC_DEFAULT  -2

474: /*MC
475:     PETSC_COMM_WORLD - the equivalent of the MPI_COMM_WORLD communicator which represents
476:            all the processs that PETSc knows about.

478:    Level: beginner

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

485: .seealso: PETSC_COMM_SELF

487: M*/
488: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;

490: /*MC
491:     PETSC_COMM_SELF - This is always MPI_COMM_SELF

493:    Level: beginner

495: .seealso: PETSC_COMM_WORLD

497: M*/
498: #define PETSC_COMM_SELF MPI_COMM_SELF

500: PETSC_EXTERN PetscBool PetscBeganMPI;
501: PETSC_EXTERN PetscBool PetscInitializeCalled;
502: PETSC_EXTERN PetscBool PetscFinalizeCalled;
503: PETSC_EXTERN PetscBool PetscCUSPSynchronize;
504: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;

506: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm));
507: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*);
508: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm*);

510: /*MC
511:    PetscMalloc - Allocates memory

513:    Synopsis:
514:     #include <petscsys.h>
515:    PetscErrorCode PetscMalloc(size_t m,void **result)

517:    Not Collective

519:    Input Parameter:
520: .  m - number of bytes to allocate

522:    Output Parameter:
523: .  result - memory allocated

525:    Level: beginner

527:    Notes:
528:    Memory is always allocated at least double aligned

530:    It is safe to allocate size 0 and pass the resulting pointer (which may or may not be NULL) to PetscFree().

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

534:   Concepts: memory allocation

536: M*/
537: #define PetscMalloc(a,b)  ((*PetscTrMalloc)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__,(void**)(b)))

539: /*MC
540:    PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment

542:    Synopsis:
543:     #include <petscsys.h>
544:    void *PetscAddrAlign(void *addr)

546:    Not Collective

548:    Input Parameters:
549: .  addr - address to align (any pointer type)

551:    Level: developer

553: .seealso: PetscMallocAlign()

555:   Concepts: memory allocation
556: M*/
557: #define PetscAddrAlign(a) (void*)((((PETSC_UINTPTR_T)(a))+(PETSC_MEMALIGN-1)) & ~(PETSC_MEMALIGN-1))

559: /*MC
560:    PetscMalloc1 - Allocates an array of memory aligned to PETSC_MEMALIGN

562:    Synopsis:
563:     #include <petscsys.h>
564:    PetscErrorCode PetscMalloc1(size_t m1,type **r1)

566:    Not Collective

568:    Input Parameter:
569: .  m1 - number of elements to allocate in 1st chunk  (may be zero)

571:    Output Parameter:
572: .  r1 - memory allocated in first chunk

574:    Level: developer

576: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc1(), PetscMalloc2()

578:   Concepts: memory allocation

580: M*/
581: #define PetscMalloc1(m1,r1) ((m1) ? PetscMalloc((m1)*sizeof(**(r1)),r1) : (*(r1) = 0,0))

583: /*MC
584:    PetscCalloc1 - Allocates a cleared (zeroed) array of memory aligned to PETSC_MEMALIGN

586:    Synopsis:
587:     #include <petscsys.h>
588:    PetscErrorCode PetscCalloc1(size_t m1,type **r1)

590:    Not Collective

592:    Input Parameter:
593: .  m1 - number of elements to allocate in 1st chunk  (may be zero)

595:    Output Parameter:
596: .  r1 - memory allocated in first chunk

598:    Level: developer

600: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc1(), PetscCalloc2()

602:   Concepts: memory allocation

604: M*/
605: #define PetscCalloc1(m1,r1) (PetscMalloc1((m1),r1) || PetscMemzero(*(r1),(m1)*sizeof(**(r1))))

607: /*MC
608:    PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN

610:    Synopsis:
611:     #include <petscsys.h>
612:    PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)

614:    Not Collective

616:    Input Parameter:
617: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
618: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

620:    Output Parameter:
621: +  r1 - memory allocated in first chunk
622: -  r2 - memory allocated in second chunk

624:    Level: developer

626: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc1(), PetscCalloc2()

628:   Concepts: memory allocation

630: M*/
631: #if !defined(PETSC_USE_MALLOC_COALESCED)
632: #define PetscMalloc2(m1,r1,m2,r2) (PetscMalloc((m1)*sizeof(**(r1)),r1) || PetscMalloc((m2)*sizeof(**(r2)),r2))
633: #else
634: #define PetscMalloc2(m1,r1,m2,r2) ((((m1)+(m2)) ? (*(r2) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(PETSC_MEMALIGN-1),r1)) : 0) \
635:                                    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),0) \
636:                                    || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0))
637: #endif

639: /*MC
640:    PetscCalloc2 - Allocates 2 cleared (zeroed) arrays of memory both aligned to PETSC_MEMALIGN

642:    Synopsis:
643:     #include <petscsys.h>
644:    PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)

646:    Not Collective

648:    Input Parameter:
649: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
650: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

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

656:    Level: developer

658: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc1(), PetscMalloc2()

660:   Concepts: memory allocation
661: M*/
662: #define PetscCalloc2(m1,r1,m2,r2) (PetscMalloc2((m1),(r1),(m2),(r2)) || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))))

664: /*MC
665:    PetscMalloc3 - Allocates 3 arrays of memory, all aligned to PETSC_MEMALIGN

667:    Synopsis:
668:     #include <petscsys.h>
669:    PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

671:    Not Collective

673:    Input Parameter:
674: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
675: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
676: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

678:    Output Parameter:
679: +  r1 - memory allocated in first chunk
680: .  r2 - memory allocated in second chunk
681: -  r3 - memory allocated in third chunk

683:    Level: developer

685: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc3(), PetscFree3()

687:   Concepts: memory allocation

689: M*/
690: #if !defined(PETSC_USE_MALLOC_COALESCED)
691: #define PetscMalloc3(m1,r1,m2,r2,m3,r3) (PetscMalloc((m1)*sizeof(**(r1)),r1) || PetscMalloc((m2)*sizeof(**(r2)),r2) || PetscMalloc((m3)*sizeof(**(r3)),r3))
692: #else
693: #define PetscMalloc3(m1,r1,m2,r2,m3,r3) ((((m1)+(m2)+(m3)) ? (*(r2) = 0,*(r3) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(m3)*sizeof(**(r3))+2*(PETSC_MEMALIGN-1),r1)) : 0) \
694:                                          || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),0) \
695:                                          || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0))
696: #endif

698: /*MC
699:    PetscCalloc3 - Allocates 3 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN

701:    Synopsis:
702:     #include <petscsys.h>
703:    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

705:    Not Collective

707:    Input Parameter:
708: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
709: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
710: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

712:    Output Parameter:
713: +  r1 - memory allocated in first chunk
714: .  r2 - memory allocated in second chunk
715: -  r3 - memory allocated in third chunk

717:    Level: developer

719: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc2(), PetscMalloc3(), PetscFree3()

721:   Concepts: memory allocation
722: M*/
723: #define PetscCalloc3(m1,r1,m2,r2,m3,r3)                                 \
724:   (PetscMalloc3((m1),(r1),(m2),(r2),(m3),(r3))                          \
725:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))))

727: /*MC
728:    PetscMalloc4 - Allocates 4 arrays of memory, all aligned to PETSC_MEMALIGN

730:    Synopsis:
731:     #include <petscsys.h>
732:    PetscErrorCode PetscMalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)

734:    Not Collective

736:    Input Parameter:
737: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
738: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
739: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
740: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

742:    Output Parameter:
743: +  r1 - memory allocated in first chunk
744: .  r2 - memory allocated in second chunk
745: .  r3 - memory allocated in third chunk
746: -  r4 - memory allocated in fourth chunk

748:    Level: developer

750: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc4(), PetscFree4()

752:   Concepts: memory allocation

754: M*/
755: #if !defined(PETSC_USE_MALLOC_COALESCED)
756: #define PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4) (PetscMalloc((m1)*sizeof(**(r1)),r1) || PetscMalloc((m2)*sizeof(**(r2)),r2) || PetscMalloc((m3)*sizeof(**(r3)),r3) || PetscMalloc((m4)*sizeof(**(r4)),r4))
757: #else
758: #define PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4)                           \
759:   ((((m1)+(m2)+(m3)+(m4)) ? (*(r2) = 0, *(r3) = 0, *(r4) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(m3)*sizeof(**(r3))+(m4)*sizeof(**(r4))+3*(PETSC_MEMALIGN-1),r1)) : 0) \
760:    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),0) \
761:    || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0) || (!(m4) ? (*(r4) = 0,0) : 0))
762: #endif

764: /*MC
765:    PetscCalloc4 - Allocates 4 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN

767:    Synopsis:
768:     #include <petscsys.h>
769:    PetscErrorCode PetscCalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)

771:    Not Collective

773:    Input Parameter:
774: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
775: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
776: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
777: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

779:    Output Parameter:
780: +  r1 - memory allocated in first chunk
781: .  r2 - memory allocated in second chunk
782: .  r3 - memory allocated in third chunk
783: -  r4 - memory allocated in fourth chunk

785:    Level: developer

787: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc4(), PetscFree4()

789:   Concepts: memory allocation

791: M*/
792: #define PetscCalloc4(m1,r1,m2,r2,m3,r3,m4,r4)                           \
793:   (PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4)                                \
794:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
795:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))))

797: /*MC
798:    PetscMalloc5 - Allocates 5 arrays of memory, all aligned to PETSC_MEMALIGN

800:    Synopsis:
801:     #include <petscsys.h>
802:    PetscErrorCode PetscMalloc5(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5)

804:    Not Collective

806:    Input Parameter:
807: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
808: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
809: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
810: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
811: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

813:    Output Parameter:
814: +  r1 - memory allocated in first chunk
815: .  r2 - memory allocated in second chunk
816: .  r3 - memory allocated in third chunk
817: .  r4 - memory allocated in fourth chunk
818: -  r5 - memory allocated in fifth chunk

820:    Level: developer

822: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc5(), PetscFree5()

824:   Concepts: memory allocation

826: M*/
827: #if !defined(PETSC_USE_MALLOC_COALESCED)
828: #define PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5) (PetscMalloc((m1)*sizeof(**(r1)),r1) || PetscMalloc((m2)*sizeof(**(r2)),r2) || PetscMalloc((m3)*sizeof(**(r3)),r3) || PetscMalloc((m4)*sizeof(**(r4)),r4) || PetscMalloc((m5)*sizeof(**(r5)),r5))
829: #else
830: #define PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5)      \
831:   ((((m1)+(m2)+(m3)+(m4)+(m5)) ? (*(r2) = 0, *(r3) = 0, *(r4) = 0,*(r5) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(m3)*sizeof(**(r3))+(m4)*sizeof(**(r4))+(m5)*sizeof(**(r5))+4*(PETSC_MEMALIGN-1),r1)) : 0) \
832:    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),*(void**)(r5) = PetscAddrAlign(*(r4)+(m4)),0) \
833:    || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0) || (!(m4) ? (*(r4) = 0,0) : 0) || (!(m5) ? (*(r5) = 0,0) : 0))
834: #endif

836: /*MC
837:    PetscCalloc5 - Allocates 5 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN

839:    Synopsis:
840:     #include <petscsys.h>
841:    PetscErrorCode PetscCalloc5(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5)

843:    Not Collective

845:    Input Parameter:
846: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
847: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
848: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
849: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
850: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

852:    Output Parameter:
853: +  r1 - memory allocated in first chunk
854: .  r2 - memory allocated in second chunk
855: .  r3 - memory allocated in third chunk
856: .  r4 - memory allocated in fourth chunk
857: -  r5 - memory allocated in fifth chunk

859:    Level: developer

861: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc5(), PetscFree5()

863:   Concepts: memory allocation

865: M*/
866: #define PetscCalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5)                     \
867:   (PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5)                          \
868:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
869:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))))

871: /*MC
872:    PetscMalloc6 - Allocates 6 arrays of memory, all aligned to PETSC_MEMALIGN

874:    Synopsis:
875:     #include <petscsys.h>
876:    PetscErrorCode PetscMalloc6(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5,size_t m6,type **r6)

878:    Not Collective

880:    Input Parameter:
881: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
882: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
883: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
884: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
885: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
886: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

888:    Output Parameter:
889: +  r1 - memory allocated in first chunk
890: .  r2 - memory allocated in second chunk
891: .  r3 - memory allocated in third chunk
892: .  r4 - memory allocated in fourth chunk
893: .  r5 - memory allocated in fifth chunk
894: -  r6 - memory allocated in sixth chunk

896:    Level: developer

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

900:   Concepts: memory allocation

902: M*/
903: #if !defined(PETSC_USE_MALLOC_COALESCED)
904: #define PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) (PetscMalloc((m1)*sizeof(**(r1)),r1) || PetscMalloc((m2)*sizeof(**(r2)),r2) || PetscMalloc((m3)*sizeof(**(r3)),r3) || PetscMalloc((m4)*sizeof(**(r4)),r4) || PetscMalloc((m5)*sizeof(**(r5)),r5) || PetscMalloc((m6)*sizeof(**(r6)),r6))
905: #else
906: #define PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) \
907:   ((((m1)+(m2)+(m3)+(m4)+(m5)+(m6)) ? (*(r2) = 0, *(r3) = 0, *(r4) = 0,*(r5) = 0,*(r6) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(m3)*sizeof(**(r3))+(m4)*sizeof(**(r4))+(m5)*sizeof(**(r5))+(m6)*sizeof(**(r6))+5*(PETSC_MEMALIGN-1),r1)) : 0) \
908:    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),*(void**)(r5) = PetscAddrAlign(*(r4)+(m4)),*(void**)(r6) = PetscAddrAlign(*(r5)+(m5)),0) \
909:    || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0) || (!(m4) ? (*(r4) = 0,0) : 0) || (!(m5) ? (*(r5) = 0,0) : 0) || (!(m6) ? (*(r6) = 0,0) : 0))
910: #endif

912: /*MC
913:    PetscCalloc6 - Allocates 6 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN

915:    Synopsis:
916:     #include <petscsys.h>
917:    PetscErrorCode PetscCalloc6(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5,size_t m6,type **r6)

919:    Not Collective

921:    Input Parameter:
922: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
923: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
924: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
925: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
926: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
927: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

929:    Output Parameter:
930: +  r1 - memory allocated in first chunk
931: .  r2 - memory allocated in second chunk
932: .  r3 - memory allocated in third chunk
933: .  r4 - memory allocated in fourth chunk
934: .  r5 - memory allocated in fifth chunk
935: -  r6 - memory allocated in sixth chunk

937:    Level: developer

939: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscMalloc6(), PetscFree6()

941:   Concepts: memory allocation
942: M*/
943: #define PetscCalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6)               \
944:   (PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6)                    \
945:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
946:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))) || PetscMemzero(*(r6),(m6)*sizeof(**(r6))))

948: /*MC
949:    PetscMalloc7 - Allocates 7 arrays of memory, all aligned to PETSC_MEMALIGN

951:    Synopsis:
952:     #include <petscsys.h>
953:    PetscErrorCode PetscMalloc7(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5,size_t m6,type **r6,size_t m7,type **r7)

955:    Not Collective

957:    Input Parameter:
958: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
959: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
960: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
961: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
962: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
963: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
964: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

966:    Output Parameter:
967: +  r1 - memory allocated in first chunk
968: .  r2 - memory allocated in second chunk
969: .  r3 - memory allocated in third chunk
970: .  r4 - memory allocated in fourth chunk
971: .  r5 - memory allocated in fifth chunk
972: .  r6 - memory allocated in sixth chunk
973: -  r7 - memory allocated in seventh chunk

975:    Level: developer

977: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc7(), PetscFree7()

979:   Concepts: memory allocation

981: M*/
982: #if !defined(PETSC_USE_MALLOC_COALESCED)
983: #define PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) (PetscMalloc((m1)*sizeof(**(r1)),r1) || PetscMalloc((m2)*sizeof(**(r2)),r2) || PetscMalloc((m3)*sizeof(**(r3)),r3) || PetscMalloc((m4)*sizeof(**(r4)),r4) || PetscMalloc((m5)*sizeof(**(r5)),r5) || PetscMalloc((m6)*sizeof(**(r6)),r6) || PetscMalloc((m7)*sizeof(**(r7)),r7))
984: #else
985: #define PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) \
986:   ((((m1)+(m2)+(m3)+(m4)+(m5)+(m6)+(m7)) ? (*(r2) = 0, *(r3) = 0, *(r4) = 0,*(r5) = 0,*(r6) = 0,*(r7) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(m3)*sizeof(**(r3))+(m4)*sizeof(**(r4))+(m5)*sizeof(**(r5))+(m6)*sizeof(**(r6))+(m7)*sizeof(**(r7))+6*(PETSC_MEMALIGN-1),r1)) : 0) \
987:    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),*(void**)(r5) = PetscAddrAlign(*(r4)+(m4)),*(void**)(r6) = PetscAddrAlign(*(r5)+(m5)),*(void**)(r7) = PetscAddrAlign(*(r6)+(m6)),0) \
988:    || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0) || (!(m4) ? (*(r4) = 0,0) : 0) || (!(m5) ? (*(r5) = 0,0) : 0) || (!(m6) ? (*(r6) = 0,0) : 0) || (!(m7) ? (*(r7) = 0,0) : 0))
989: #endif

991: /*MC
992:    PetscCalloc7 - Allocates 7 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN

994:    Synopsis:
995:     #include <petscsys.h>
996:    PetscErrorCode PetscCalloc7(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4,size_t m5,type **r5,size_t m6,type **r6,size_t m7,type **r7)

998:    Not Collective

1000:    Input Parameter:
1001: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
1002: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
1003: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
1004: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
1005: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
1006: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
1007: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

1009:    Output Parameter:
1010: +  r1 - memory allocated in first chunk
1011: .  r2 - memory allocated in second chunk
1012: .  r3 - memory allocated in third chunk
1013: .  r4 - memory allocated in fourth chunk
1014: .  r5 - memory allocated in fifth chunk
1015: .  r6 - memory allocated in sixth chunk
1016: -  r7 - memory allocated in seventh chunk

1018:    Level: developer

1020: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscMalloc7(), PetscFree7()

1022:   Concepts: memory allocation
1023: M*/
1024: #define PetscCalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7)         \
1025:   (PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7)              \
1026:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
1027:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))) || PetscMemzero(*(r6),(m6)*sizeof(**(r6))) \
1028:    || PetscMemzero(*(r7),(m7)*sizeof(**(r7))))

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

1033:    Synopsis:
1034:     #include <petscsys.h>
1035:    PetscErrorCode PetscNew(type **result)

1037:    Not Collective

1039:    Output Parameter:
1040: .  result - memory allocated, sized to match pointer type

1042:    Level: beginner

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

1046:   Concepts: memory allocation

1048: M*/
1049: #define PetscNew(b)      PetscCalloc1(1,(b))

1051: /*MC
1052:    PetscNewLog - Allocates memory of a type matching pointer, zeros the memory! Aligned to PETSC_MEMALIGN. Associates the memory allocated
1053:          with the given object using PetscLogObjectMemory().

1055:    Synopsis:
1056:     #include <petscsys.h>
1057:    PetscErrorCode PetscNewLog(PetscObject obj,type **result)

1059:    Not Collective

1061:    Input Parameter:
1062: .  obj - object memory is logged to

1064:    Output Parameter:
1065: .  result - memory allocated, sized to match pointer type

1067:    Level: developer

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

1071:   Concepts: memory allocation

1073: M*/
1074: #define PetscNewLog(o,b) (PetscNew((b)) || ((o) ? PetscLogObjectMemory((PetscObject)o,sizeof(**(b))) : 0))

1076: /*MC
1077:    PetscFree - Frees memory

1079:    Synopsis:
1080:     #include <petscsys.h>
1081:    PetscErrorCode PetscFree(void *memory)

1083:    Not Collective

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

1088:    Level: beginner

1090:    Notes:
1091:    Memory must have been obtained with PetscNew() or PetscMalloc().
1092:    It is safe to call PetscFree() on a NULL pointer.

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

1096:   Concepts: memory allocation

1098: M*/
1099: #define PetscFree(a)   ((*PetscTrFree)((void*)(a),__LINE__,PETSC_FUNCTION_NAME,__FILE__) || ((a) = 0,0))

1101: /*MC
1102:    PetscFreeVoid - Frees memory

1104:    Synopsis:
1105:     #include <petscsys.h>
1106:    void PetscFreeVoid(void *memory)

1108:    Not Collective

1110:    Input Parameter:
1111: .   memory - memory to free

1113:    Level: beginner

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

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

1119:   Concepts: memory allocation

1121: M*/
1122: #define PetscFreeVoid(a) ((*PetscTrFree)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__),(a) = 0)


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

1128:    Synopsis:
1129:     #include <petscsys.h>
1130:    PetscErrorCode PetscFree2(void *memory1,void *memory2)

1132:    Not Collective

1134:    Input Parameter:
1135: +   memory1 - memory to free
1136: -   memory2 - 2nd memory to free

1138:    Level: developer

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

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

1144:   Concepts: memory allocation

1146: M*/
1147: #if defined(PETSC_USE_DEBUG)
1148: #define PetscFree2(m1,m2)   (PetscFree(m2) || PetscFree(m1))
1149: #else
1150: #define PetscFree2(m1,m2)   ((m1) ? ((m2)=0,PetscFree(m1)) : ((m1)=0,PetscFree(m2)))
1151: #endif

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

1156:    Synopsis:
1157:     #include <petscsys.h>
1158:    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)

1160:    Not Collective

1162:    Input Parameter:
1163: +   memory1 - memory to free
1164: .   memory2 - 2nd memory to free
1165: -   memory3 - 3rd memory to free

1167:    Level: developer

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

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

1173:   Concepts: memory allocation

1175: M*/
1176: #if defined(PETSC_USE_DEBUG)
1177: #define PetscFree3(m1,m2,m3)   (PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1178: #else
1179: #define PetscFree3(m1,m2,m3)   ((m1) ? ((m3)=0,(m2)=0,PetscFree(m1)) : ((m2) ? ((m3)=0,(m1)=0,PetscFree(m2)) : ((m2)=0,(m1)=0,PetscFree(m3))))
1180: #endif

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

1185:    Synopsis:
1186:     #include <petscsys.h>
1187:    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)

1189:    Not Collective

1191:    Input Parameter:
1192: +   m1 - memory to free
1193: .   m2 - 2nd memory to free
1194: .   m3 - 3rd memory to free
1195: -   m4 - 4th memory to free

1197:    Level: developer

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

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

1203:   Concepts: memory allocation

1205: M*/
1206: #if defined(PETSC_USE_DEBUG)
1207: #define PetscFree4(m1,m2,m3,m4)   (PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1208: #else
1209: #define PetscFree4(m1,m2,m3,m4)   ((m1) ? ((m4)=0,(m3)=0,(m2)=0,PetscFree(m1)) : ((m2) ? ((m4)=0,(m3)=0,(m1)=0,PetscFree(m2)) : ((m3) ? ((m4)=0,(m2)=0,(m1)=0,PetscFree(m3)) : ((m3)=0,(m2)=0,(m1)=0,PetscFree(m4)))))
1210: #endif

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

1215:    Synopsis:
1216:     #include <petscsys.h>
1217:    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)

1219:    Not Collective

1221:    Input Parameter:
1222: +   m1 - memory to free
1223: .   m2 - 2nd memory to free
1224: .   m3 - 3rd memory to free
1225: .   m4 - 4th memory to free
1226: -   m5 - 5th memory to free

1228:    Level: developer

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

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

1234:   Concepts: memory allocation

1236: M*/
1237: #if defined(PETSC_USE_DEBUG)
1238: #define PetscFree5(m1,m2,m3,m4,m5)   (PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1239: #else
1240: #define PetscFree5(m1,m2,m3,m4,m5)   ((m1) ? ((m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1)) : ((m2) ? ((m5)=0,(m4)=0,(m3)=0,(m1)=0,PetscFree(m2)) : ((m3) ? ((m5)=0,(m4)=0,(m2)=0,(m1)=0,PetscFree(m3)) : \
1241:                                      ((m4) ? ((m5)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m4)) : ((m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m5))))))
1242: #endif


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

1248:    Synopsis:
1249:     #include <petscsys.h>
1250:    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)

1252:    Not Collective

1254:    Input Parameter:
1255: +   m1 - memory to free
1256: .   m2 - 2nd memory to free
1257: .   m3 - 3rd memory to free
1258: .   m4 - 4th memory to free
1259: .   m5 - 5th memory to free
1260: -   m6 - 6th memory to free


1263:    Level: developer

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

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

1269:   Concepts: memory allocation

1271: M*/
1272: #if defined(PETSC_USE_DEBUG)
1273: #define PetscFree6(m1,m2,m3,m4,m5,m6)   (PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1274: #else
1275: #define PetscFree6(m1,m2,m3,m4,m5,m6)   ((m1) ? ((m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1)) : ((m2) ? ((m6)=0,(m5)=0,(m4)=0,(m3)=0,(m1)=0,PetscFree(m2)) : \
1276:                                         ((m3) ? ((m6)=0,(m5)=0,(m4)=0,(m2)=0,(m1)=0,PetscFree(m3)) : ((m4) ? ((m6)=0,(m5)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m4)) : \
1277:                                         ((m5) ? ((m6)=0,(m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m5)) : ((m5)=0,(m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m6)))))))
1278: #endif

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

1283:    Synopsis:
1284:     #include <petscsys.h>
1285:    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)

1287:    Not Collective

1289:    Input Parameter:
1290: +   m1 - memory to free
1291: .   m2 - 2nd memory to free
1292: .   m3 - 3rd memory to free
1293: .   m4 - 4th memory to free
1294: .   m5 - 5th memory to free
1295: .   m6 - 6th memory to free
1296: -   m7 - 7th memory to free


1299:    Level: developer

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

1303: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1304:           PetscMalloc7()

1306:   Concepts: memory allocation

1308: M*/
1309: #if defined(PETSC_USE_DEBUG)
1310: #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   (PetscFree(m7) || PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1311: #else
1312: #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   ((m1) ? ((m7)=0,(m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1)) : ((m2) ? ((m7)=0,(m6)=0,(m5)=0,(m4)=0,(m3)=0,(m1)=0,PetscFree(m2)) : \
1313:                                            ((m3) ? ((m7)=0,(m6)=0,(m5)=0,(m4)=0,(m2)=0,(m1)=0,PetscFree(m3)) : ((m4) ? ((m7)=0,(m6)=0,(m5)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m4)) : \
1314:                                            ((m5) ? ((m7)=0,(m6)=0,(m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m5)) : ((m6) ? ((m7)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m6)) : \
1315:                                                    ((m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m7))))))))
1316: #endif

1318: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,int,const char[],const char[],void**);
1319: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1320: PETSC_EXTERN PetscErrorCode PetscMallocSet(PetscErrorCode (*)(size_t,int,const char[],const char[],void**),PetscErrorCode (*)(void*,int,const char[],const char[]));
1321: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);

1323: /*
1324:     PetscLogDouble variables are used to contain double precision numbers
1325:   that are not used in the numerical computations, but rather in logging,
1326:   timing etc.
1327: */
1328: typedef double PetscLogDouble;
1329: #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE

1331: /*
1332:    Routines for tracing memory corruption/bleeding with default PETSc  memory allocation
1333: */
1334: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1335: PETSC_EXTERN PetscErrorCode PetscMallocDumpLog(FILE *);
1336: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1337: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1338: PETSC_EXTERN PetscErrorCode PetscMallocDebug(PetscBool);
1339: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*);
1340: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1341: PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLog(void);
1342: PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLogThreshold(PetscLogDouble);
1343: PETSC_EXTERN PetscErrorCode PetscMallocGetDumpLog(PetscBool*);

1345: /*E
1346:     PetscDataType - Used for handling different basic data types.

1348:    Level: beginner

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

1352: .seealso: PetscBinaryRead(), PetscBinaryWrite(), PetscDataTypeToMPIDataType(),
1353:           PetscDataTypeGetSize()

1355: E*/
1356: typedef enum {PETSC_INT = 0,PETSC_DOUBLE = 1,PETSC_COMPLEX = 2, PETSC_LONG = 3 ,PETSC_SHORT = 4,PETSC_FLOAT = 5,
1357:               PETSC_CHAR = 6,PETSC_BIT_LOGICAL = 7,PETSC_ENUM = 8,PETSC_BOOL=9, PETSC___FLOAT128 = 10,PETSC_OBJECT = 11, PETSC_FUNCTION = 12, PETSC_STRING = 12} PetscDataType;
1358: PETSC_EXTERN const char *const PetscDataTypes[];

1360: #if defined(PETSC_USE_COMPLEX)
1361: #define  PETSC_SCALAR  PETSC_COMPLEX
1362: #else
1363: #if defined(PETSC_USE_REAL_SINGLE)
1364: #define  PETSC_SCALAR  PETSC_FLOAT
1365: #elif defined(PETSC_USE_REAL___FLOAT128)
1366: #define  PETSC_SCALAR  PETSC___FLOAT128
1367: #else
1368: #define  PETSC_SCALAR  PETSC_DOUBLE
1369: #endif
1370: #endif
1371: #if defined(PETSC_USE_REAL_SINGLE)
1372: #define  PETSC_REAL  PETSC_FLOAT
1373: #elif defined(PETSC_USE_REAL___FLOAT128)
1374: #define  PETSC_REAL  PETSC___FLOAT128
1375: #else
1376: #define  PETSC_REAL  PETSC_DOUBLE
1377: #endif
1378: #define  PETSC_FORTRANADDR  PETSC_LONG

1380: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1381: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1382: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1383: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);

1385: /*
1386:     Basic memory and string operations. These are usually simple wrappers
1387:    around the basic Unix system calls, but a few of them have additional
1388:    functionality and/or error checking.
1389: */
1390: PETSC_EXTERN PetscErrorCode PetscBitMemcpy(void*,PetscInt,const void*,PetscInt,PetscInt,PetscDataType);
1391: PETSC_EXTERN PetscErrorCode PetscMemmove(void*,void *,size_t);
1392: PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1393: PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1394: PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***);
1395: PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1396: PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool  *);
1397: PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool  *);
1398: PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1399: PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1400: PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1401: PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1402: PETSC_EXTERN PetscErrorCode PetscStrncat(char[],const char[],size_t);
1403: PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1404: PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1405: PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1406: PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1407: PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1408: PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1409: PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1410: PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1411: PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*);
1412: PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1413: PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1414: PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1415: PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1416: PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);

1418: PETSC_EXTERN void PetscStrcmpNoError(const char[],const char[],PetscBool  *);

1420: /*S
1421:     PetscToken - 'Token' used for managing tokenizing strings

1423:   Level: intermediate

1425: .seealso: PetscTokenCreate(), PetscTokenFind(), PetscTokenDestroy()
1426: S*/
1427: typedef struct _p_PetscToken* PetscToken;

1429: PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1430: PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1431: PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);

1433: PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1434: PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);

1436: /*
1437:    These are  MPI operations for MPI_Allreduce() etc
1438: */
1439: PETSC_EXTERN MPI_Op PetscMaxSum_Op;
1440: #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128)
1441: PETSC_EXTERN MPI_Op MPIU_SUM;
1442: #else
1443: #define MPIU_SUM MPI_SUM
1444: #endif
1445: #if defined(PETSC_USE_REAL___FLOAT128)
1446: PETSC_EXTERN MPI_Op MPIU_MAX;
1447: PETSC_EXTERN MPI_Op MPIU_MIN;
1448: #else
1449: #define MPIU_MAX MPI_MAX
1450: #define MPIU_MIN MPI_MIN
1451: #endif
1452: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);

1454: PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1455: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);

1457: /*S
1458:      PetscObject - any PETSc object, PetscViewer, Mat, Vec, KSP etc

1460:    Level: beginner

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

1464: .seealso:  PetscObjectDestroy(), PetscObjectView(), PetscObjectGetName(), PetscObjectSetName(), PetscObjectReference(), PetscObjectDereferenc()
1465: S*/
1466: typedef struct _p_PetscObject* PetscObject;

1468: /*MC
1469:     PetscObjectId - unique integer Id for a PetscObject

1471:     Level: developer

1473:     Notes: Unlike pointer values, object ids are never reused.

1475: .seealso: PetscObjectState, PetscObjectGetId()
1476: M*/
1477: typedef Petsc64bitInt PetscObjectId;

1479: /*MC
1480:     PetscObjectState - integer state for a PetscObject

1482:     Level: developer

1484:     Notes:
1485:     Object state is always-increasing and (for objects that track state) can be used to determine if an object has
1486:     changed since the last time you interacted with it.  It is 64-bit so that it will not overflow for a very long time.

1488: .seealso: PetscObjectId, PetscObjectStateGet(), PetscObjectStateIncrease(), PetscObjectStateSet()
1489: M*/
1490: typedef Petsc64bitInt PetscObjectState;

1492: /*S
1493:      PetscFunctionList - Linked list of functions, possibly stored in dynamic libraries, accessed
1494:       by string name

1496:    Level: advanced

1498: .seealso:  PetscFunctionListAdd(), PetscFunctionListDestroy(), PetscOpFlist
1499: S*/
1500: typedef struct _n_PetscFunctionList *PetscFunctionList;

1502: /*E
1503:   PetscFileMode - Access mode for a file.

1505:   Level: beginner

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

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

1511:   FILE_MODE_APPEND - open a file at end for writing

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

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

1517: .seealso: PetscViewerFileSetMode()
1518: E*/
1519: typedef enum {FILE_MODE_READ, FILE_MODE_WRITE, FILE_MODE_APPEND, FILE_MODE_UPDATE, FILE_MODE_APPEND_UPDATE} PetscFileMode;
1520: extern const char *const PetscFileModes[];

1522: /*
1523:     Defines PETSc error handling.
1524: */
1525: #include <petscerror.h>

1527: #define PETSC_SMALLEST_CLASSID  1211211
1528: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1529: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1530: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);

1532: /*
1533:    Routines that get memory usage information from the OS
1534: */
1535: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1536: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1537: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);

1539: PETSC_EXTERN PetscErrorCode PetscInfoAllow(PetscBool ,const char []);
1540: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);

1542: /*
1543:    Initialization of PETSc
1544: */
1545: PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1546: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1547: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1548: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1549: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1550: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1551: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1552: PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1553: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1554: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);

1556: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1557: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);

1559: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1560: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1561: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1562: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);

1564: /*
1565:      These are so that in extern C code we can caste function pointers to non-extern C
1566:    function pointers. Since the regular C++ code expects its function pointers to be C++
1567: */
1568: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1569: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1570: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);

1572: /*
1573:     Functions that can act on any PETSc object.
1574: */
1575: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1576: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1577: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1578: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1579: PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1580: PETSC_EXTERN PetscErrorCode PetscObjectSetPrecision(PetscObject,PetscPrecision);
1581: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1582: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1583: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1584: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1585: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1586: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1587: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1588: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1589: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1590: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt *);
1591: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1592: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1593: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject *);
1594: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void));
1595: #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d))
1596: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1597: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1598: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt *);
1599: PETSC_EXTERN PetscErrorCode PetscObjectAddOptionsHandler(PetscObject,PetscErrorCode (*)(PetscObject,void*),PetscErrorCode (*)(PetscObject,void*),void*);
1600: PETSC_EXTERN PetscErrorCode PetscObjectProcessOptionsHandlers(PetscObject);
1601: PETSC_EXTERN PetscErrorCode PetscObjectDestroyOptionsHandlers(PetscObject);
1602: PETSC_EXTERN PetscErrorCode PetscObjectsGetGlobalNumbering(MPI_Comm,PetscInt,PetscObject*,PetscInt*,PetscInt*);

1604: #include <petscviewertypes.h>
1605: #include <petscoptions.h>

1607: PETSC_EXTERN PetscErrorCode PetscMemoryShowUsage(PetscViewer,const char[]);
1608: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1609: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1610: #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1611: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1612: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1613: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1614: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1615: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1616: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1617: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1618: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1619: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,const char[],const char[]);
1620: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1621: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1622: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1623: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1624: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);

1626: #if defined(PETSC_HAVE_SAWS)
1627: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1628: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1629: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1630: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1631: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1632: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1633: PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1634: PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1635: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1636: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);

1638: #else
1639: #define PetscSAWsBlock()                        0
1640: #define PetscObjectSAWsViewOff(obj)             0
1641: #define PetscObjectSAWsSetBlock(obj,flg)        0
1642: #define PetscObjectSAWsBlock(obj)               0
1643: #define PetscObjectSAWsGrantAccess(obj)         0
1644: #define PetscObjectSAWsTakeAccess(obj)          0
1645: #define PetscStackViewSAWs()                    0
1646: #define PetscStackSAWsViewOff()                 0
1647: #define PetscStackSAWsTakeAccess()
1648: #define PetscStackSAWsGrantAccess()

1650: #endif

1652: typedef void* PetscDLHandle;
1653: typedef enum {PETSC_DL_DECIDE=0,PETSC_DL_NOW=1,PETSC_DL_LOCAL=2} PetscDLMode;
1654: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1655: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1656: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);


1659: #if defined(PETSC_USE_DEBUG)
1660: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1661: #endif
1662: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);

1664: /*S
1665:      PetscObjectList - Linked list of PETSc objects, each accessable by string name

1667:    Level: developer

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

1671: .seealso:  PetscObjectListAdd(), PetscObjectListDestroy(), PetscObjectListFind(), PetscObjectCompose(), PetscObjectQuery(), PetscFunctionList
1672: S*/
1673: typedef struct _n_PetscObjectList *PetscObjectList;

1675: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1676: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1677: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1678: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1679: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1680: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);

1682: /*
1683:     Dynamic library lists. Lists of names of routines in objects or in dynamic
1684:   link libraries that will be loaded as needed.
1685: */

1687: #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1688: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1689: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1690: #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1691: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1692: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[]);
1693: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1694: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1695: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);

1697: /*S
1698:      PetscDLLibrary - Linked list of dynamics libraries to search for functions

1700:    Level: advanced

1702: .seealso:  PetscDLLibraryOpen()
1703: S*/
1704: typedef struct _n_PetscDLLibrary *PetscDLLibrary;
1705: PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1706: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1707: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1708: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1709: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1710: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1711: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1712: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);

1714: /*
1715:      Useful utility routines
1716: */
1717: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1718: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1719: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1720: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1721: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1722: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);

1724: /*
1725:     PetscNot - negates a logical type value and returns result as a PetscBool

1727:     Notes: This is useful in cases like
1728: $     int        *a;
1729: $     PetscBool  flag = PetscNot(a)
1730:      where !a does not return a PetscBool  because we cannot provide a cast from int to PetscBool  in C.
1731: */
1732: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)

1734: #if defined(PETSC_HAVE_VALGRIND)
1735: #  include <valgrind/valgrind.h>
1736: #  define PETSC_RUNNING_ON_VALGRIND RUNNING_ON_VALGRIND
1737: #else
1738: #  define PETSC_RUNNING_ON_VALGRIND PETSC_FALSE
1739: #endif


1742: /*MC
1743:     PetscHelpPrintf - Prints help messages.

1745:    Synopsis:
1746:     #include <petscsys.h>
1747:      PetscErrorCode (*PetscHelpPrintf)(const char format[],...);

1749:     Not Collective

1751:     Input Parameters:
1752: .   format - the usual printf() format string

1754:    Level: developer

1756:     Fortran Note:
1757:     This routine is not supported in Fortran.

1759:     Concepts: help messages^printing
1760:     Concepts: printing^help messages

1762: .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1763: M*/
1764: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);

1766: /*
1767:      Defines PETSc profiling.
1768: */
1769: #include <petsclog.h>



1773: /*
1774:       Simple PETSc parallel IO for ASCII printing
1775: */
1776: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1777: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1778: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1779: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1780: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1781: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1782: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);

1784: /* These are used internally by PETSc ASCII IO routines*/
1785: #include <stdarg.h>
1786: PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list);
1787: PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE*,const char[],va_list);
1788: PETSC_EXTERN PetscErrorCode PetscVFPrintfDefault(FILE*,const char[],va_list);

1790: #if defined(PETSC_HAVE_MATLAB_ENGINE)
1791: PETSC_EXTERN PetscErrorCode PetscVFPrintf_Matlab(FILE*,const char[],va_list);
1792: #endif

1794: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1795: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1796: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);

1798: #if defined(PETSC_HAVE_POPEN)
1799: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1800: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*,int*);
1801: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1802: #endif

1804: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1805: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1806: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*);
1807: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1808: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1809: PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1810: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);

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

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

1817:    Level: advanced

1819: .seealso:  PetscObject, PetscContainerCreate()
1820: S*/
1821: PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1822: typedef struct _p_PetscContainer*  PetscContainer;
1823: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void **);
1824: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void *);
1825: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1826: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer *);
1827: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void*));

1829: /*
1830:    For use in debuggers
1831: */
1832: PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1833: PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1834: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1835: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1836: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);

1838: #include <stddef.h>
1839: #include <string.h>             /* for memcpy, memset */
1840: #if defined(PETSC_HAVE_STDLIB_H)
1841: #include <stdlib.h>
1842: #endif

1844: #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1845: #include <xmmintrin.h>
1846: #endif

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

1855:    Not Collective

1857:    Input Parameters:
1858: +  b - pointer to initial memory space
1859: -  n - length (in bytes) of space to copy

1861:    Output Parameter:
1862: .  a - pointer to copy space

1864:    Level: intermediate

1866:    Compile Option:
1867:     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1868:                                   for memory copies on double precision values.
1869:     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1870:                                   for memory copies on double precision values.
1871:     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1872:                                   for memory copies on double precision values.

1874:    Note:
1875:    This routine is analogous to memcpy().

1877:    Developer Note: this is inlined for fastest performance

1879:   Concepts: memory^copying
1880:   Concepts: copying^memory

1882: .seealso: PetscMemmove()

1884: @*/
1885: PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n)
1886: {
1887: #if defined(PETSC_USE_DEBUG)
1888:   unsigned long al = (unsigned long) a,bl = (unsigned long) b;
1889:   unsigned long nl = (unsigned long) n;
1891:   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1892:   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1893: #else
1895: #endif
1896:   if (a != b) {
1897: #if defined(PETSC_USE_DEBUG)
1898:     if ((al > bl && (al - bl) < nl) || (bl - al) < nl)  SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1899:               or make sure your copy regions and lengths are correct. \n\
1900:               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1901: #endif
1902: #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1903:    if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1904:       size_t len = n/sizeof(PetscScalar);
1905: #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1906:       PetscBLASInt   one = 1,blen;
1908:       PetscBLASIntCast(len,&blen);
1909:       PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one));
1910: #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1911:       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1912: #else
1913:       size_t      i;
1914:       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1915:       for (i=0; i<len; i++) y[i] = x[i];
1916: #endif
1917:     } else {
1918:       memcpy((char*)(a),(char*)(b),n);
1919:     }
1920: #else
1921:     memcpy((char*)(a),(char*)(b),n);
1922: #endif
1923:   }
1924:   return(0);
1925: }

1927: /*@C
1928:    PetscMemzero - Zeros the specified memory.

1930:    Not Collective

1932:    Input Parameters:
1933: +  a - pointer to beginning memory location
1934: -  n - length (in bytes) of memory to initialize

1936:    Level: intermediate

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

1942:    Developer Note: this is inlined for fastest performance

1944:    Concepts: memory^zeroing
1945:    Concepts: zeroing^memory

1947: .seealso: PetscMemcpy()
1948: @*/
1949: PETSC_STATIC_INLINE PetscErrorCode  PetscMemzero(void *a,size_t n)
1950: {
1951:   if (n > 0) {
1952: #if defined(PETSC_USE_DEBUG)
1953:     if (!a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer");
1954: #endif
1955: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1956:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1957:       size_t      i,len = n/sizeof(PetscScalar);
1958:       PetscScalar *x = (PetscScalar*)a;
1959:       for (i=0; i<len; i++) x[i] = 0.0;
1960:     } else {
1961: #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1962:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1963:       PetscInt len = n/sizeof(PetscScalar);
1964:       fortranzero_(&len,(PetscScalar*)a);
1965:     } else {
1966: #endif
1967: #if defined(PETSC_PREFER_BZERO)
1968:       bzero((char *)a,n);
1969: #else
1970:       memset((char*)a,0,n);
1971: #endif
1972: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1973:     }
1974: #endif
1975:   }
1976:   return 0;
1977: }

1979: /*MC
1980:    PetscPrefetchBlock - Prefetches a block of memory

1982:    Synopsis:
1983:     #include <petscsys.h>
1984:     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)

1986:    Not Collective

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

1994:    Level: developer

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

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

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

2009:    Concepts: memory
2010: M*/
2011: #define PetscPrefetchBlock(a,n,rw,t) do {                               \
2012:     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
2013:     for ( ; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
2014:   } while (0)

2016: /*
2017:       Determine if some of the kernel computation routines use
2018:    Fortran (rather than C) for the numerical calculations. On some machines
2019:    and compilers (like complex numbers) the Fortran version of the routines
2020:    is faster than the C/C++ versions. The flag --with-fortran-kernels
2021:    should be used with ./configure to turn these on.
2022: */
2023: #if defined(PETSC_USE_FORTRAN_KERNELS)

2025: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
2026: #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
2027: #endif

2029: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
2030: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
2031: #endif

2033: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
2034: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
2035: #endif

2037: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
2038: #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
2039: #endif

2041: #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
2042: #define PETSC_USE_FORTRAN_KERNEL_NORM
2043: #endif

2045: #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
2046: #define PETSC_USE_FORTRAN_KERNEL_MAXPY
2047: #endif

2049: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
2050: #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
2051: #endif

2053: #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
2054: #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
2055: #endif

2057: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
2058: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
2059: #endif

2061: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
2062: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
2063: #endif

2065: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
2066: #define PETSC_USE_FORTRAN_KERNEL_MDOT
2067: #endif

2069: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
2070: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
2071: #endif

2073: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
2074: #define PETSC_USE_FORTRAN_KERNEL_AYPX
2075: #endif

2077: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
2078: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
2079: #endif

2081: #endif

2083: /*
2084:     Macros for indicating code that should be compiled with a C interface,
2085:    rather than a C++ interface. Any routines that are dynamically loaded
2086:    (such as the PCCreate_XXX() routines) must be wrapped so that the name
2087:    mangler does not change the functions symbol name. This just hides the
2088:    ugly extern "C" {} wrappers.
2089: */
2090: #if defined(__cplusplus)
2091: #define EXTERN_C_BEGIN extern "C" {
2092: #define EXTERN_C_END }
2093: #else
2094: #define EXTERN_C_BEGIN
2095: #define EXTERN_C_END
2096: #endif

2098: /* --------------------------------------------------------------------*/

2100: /*MC
2101:     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
2102:         communication

2104:    Level: beginner

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

2108: .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
2109: M*/

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


2117:    Level: beginner

2119: .seealso: PetscReal, PassiveReal, PassiveScalar, MPIU_SCALAR, PetscInt
2120: M*/

2122: /*MC
2123:     PetscComplex - PETSc type that represents a complex number with precision matching that of PetscReal.

2125:    Synopsis:
2126:    #define PETSC_DESIRE_COMPLEX
2127:    #include <petscsys.h>
2128:    PetscComplex number = 1. + 2.*PETSC_i;

2130:    Level: beginner

2132:    Note:
2133:    Complex numbers are automatically available if PETSc was configured --with-scalar-type=complex (in which case
2134:    PetscComplex will match PetscScalar), otherwise the macro PETSC_DESIRE_COMPLEX must be defined before including any
2135:    PETSc headers. If the compiler supports complex numbers, PetscComplex and associated variables and functions will be
2136:    defined and PETSC_HAVE_COMPLEX will be set.

2138: .seealso: PetscReal, PetscComplex, MPIU_COMPLEX, PetscInt, PETSC_i
2139: M*/

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

2144:    Level: beginner

2146: .seealso: PetscScalar, PassiveReal, PassiveScalar
2147: M*/

2149: /*MC
2150:     PassiveScalar - PETSc type that represents a PetscScalar
2151:    Level: beginner

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

2156: .seealso: PetscReal, PassiveReal, PetscScalar
2157: M*/

2159: /*MC
2160:     PassiveReal - PETSc type that represents a PetscReal

2162:    Level: beginner

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

2167: .seealso: PetscScalar, PetscReal, PassiveScalar
2168: M*/

2170: /*MC
2171:     MPIU_SCALAR - MPI datatype corresponding to PetscScalar

2173:    Level: beginner

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

2178: .seealso: PetscReal, PassiveReal, PassiveScalar, PetscScalar, MPIU_INT
2179: M*/

2181: #if defined(PETSC_HAVE_MPIIO)
2182: #if !defined(PETSC_WORDS_BIGENDIAN)
2183: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2184: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2185: #else
2186: #define MPIU_File_write_all(a,b,c,d,e) MPI_File_write_all(a,b,c,d,e)
2187: #define MPIU_File_read_all(a,b,c,d,e) MPI_File_read_all(a,b,c,d,e)
2188: #endif
2189: #endif

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

2193: /* Limit MPI to 32-bits */
2194: #define PETSC_MPI_INT_MAX  2147483647
2195: #define PETSC_MPI_INT_MIN -2147483647
2196: /* Limit BLAS to 32-bits */
2197: #define PETSC_BLAS_INT_MAX  2147483647
2198: #define PETSC_BLAS_INT_MIN -2147483647

2202: PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
2203: {
2205: #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
2206:   if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK");
2207: #endif
2208:   *b =  (PetscBLASInt)(a);
2209:   return(0);
2210: }

2214: PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
2215: {
2217: #if defined(PETSC_USE_64BIT_INDICES)
2218:   if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for MPI");
2219: #endif
2220:   *b =  (PetscMPIInt)(a);
2221:   return(0);
2222: }


2225: /*
2226:      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2227: */
2228: #if defined(hz)
2229: #undef hz
2230: #endif

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


2235: #if defined(PETSC_HAVE_LIMITS_H)
2236: #include <limits.h>
2237: #endif
2238: #if defined(PETSC_HAVE_SYS_PARAM_H)
2239: #include <sys/param.h>
2240: #endif
2241: #if defined(PETSC_HAVE_SYS_TYPES_H)
2242: #include <sys/types.h>
2243: #endif
2244: #if defined(MAXPATHLEN)
2245: #  define PETSC_MAX_PATH_LEN     MAXPATHLEN
2246: #elif defined(MAX_PATH)
2247: #  define PETSC_MAX_PATH_LEN     MAX_PATH
2248: #elif defined(_MAX_PATH)
2249: #  define PETSC_MAX_PATH_LEN     _MAX_PATH
2250: #else
2251: #  define PETSC_MAX_PATH_LEN     4096
2252: #endif

2254: /*MC

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

2258: $    1) classic Fortran 77 style
2259: $#include "finclude/petscXXX.h" to work with material from the XXX component of PETSc
2260: $       XXX variablename
2261: $      You cannot use this approach if you wish to use the Fortran 90 specific PETSc routines
2262: $      which end in F90; such as VecGetArrayF90()
2263: $
2264: $    2) classic Fortran 90 style
2265: $#include "finclude/petscXXX.h"
2266: $#include "finclude/petscXXX.h90" to work with material from the XXX component of PETSc
2267: $       XXX variablename
2268: $
2269: $    3) Using Fortran modules
2270: $#include "finclude/petscXXXdef.h"
2271: $         use petscXXXX
2272: $       XXX variablename
2273: $
2274: $    4) Use Fortran modules and Fortran data types for PETSc types
2275: $#include "finclude/petscXXXdef.h"
2276: $         use petscXXXX
2277: $       type(XXX) variablename
2278: $      To use this approach you must ./configure PETSc with the additional
2279: $      option --with-fortran-datatypes You cannot use the type(XXX) declaration approach without using Fortran modules

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

2283: $    3a) skip the #include BUT you cannot use any PETSc data type names like Vec, Mat, PetscInt, PetscErrorCode etc
2284: $        and you must declare the variables as integer, for example
2285: $        integer variablename
2286: $
2287: $    4a) skip the #include, you use the object types like type(Vec) type(Mat) but cannot use the data type
2288: $        names like PetscErrorCode, PetscInt etc. again for those you must use integer

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

2293:    Fortran type checking with interfaces is strick, this means you cannot pass a scalar value when an array value
2294: is expected (even though it is legal Fortran). For example when setting a single value in a matrix with MatSetValues()
2295: you cannot have something like
2296: $      PetscInt row,col
2297: $      PetscScalar val
2298: $        ...
2299: $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)
2300: You must instead have
2301: $      PetscInt row(1),col(1)
2302: $      PetscScalar val(1)
2303: $        ...
2304: $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)


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

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

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

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

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

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

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

2329:     Level: beginner

2331: M*/

2333: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2334: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2335: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2336: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2337: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2338: PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2339: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);

2341: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2342: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2343: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2344: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2345: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2346: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2347: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt*,PetscInt*,PetscInt*);
2348: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2349: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2350: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2351: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2352: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2353: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2354: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2355: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2356: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2357: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2358: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt*,const PetscInt*,PetscInt,const PetscInt*,const PetscInt*,PetscInt*,PetscInt**,PetscInt**);

2360: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2361: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);

2363: /*J
2364:     PetscRandomType - String with the name of a PETSc randomizer

2366:    Level: beginner

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

2371: .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2372: J*/
2373: typedef const char* PetscRandomType;
2374: #define PETSCRAND       "rand"
2375: #define PETSCRAND48     "rand48"
2376: #define PETSCSPRNG      "sprng"

2378: /* Logging support */
2379: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;

2381: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);

2383: /*S
2384:      PetscRandom - Abstract PETSc object that manages generating random numbers

2386:    Level: intermediate

2388:   Concepts: random numbers

2390: .seealso:  PetscRandomCreate(), PetscRandomGetValue(), PetscRandomType
2391: S*/
2392: typedef struct _p_PetscRandom*   PetscRandom;

2394: /* Dynamic creation and loading functions */
2395: PETSC_EXTERN PetscFunctionList PetscRandomList;
2396: PETSC_EXTERN PetscBool         PetscRandomRegisterAllCalled;

2398: PETSC_EXTERN PetscErrorCode PetscRandomRegisterAll(void);
2399: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2400: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2401: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2402: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType*);
2403: PETSC_STATIC_INLINE PetscErrorCode PetscRandomViewFromOptions(PetscRandom A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
2404: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);

2406: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2407: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2408: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2409: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2410: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2411: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2412: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2413: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2414: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);

2416: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2417: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2418: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2419: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2420: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2421: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2422: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);

2424: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscDataType);
2425: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscDataType);
2426: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,void*,PetscInt,PetscDataType,PetscBool );
2427: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,void*,PetscInt,PetscDataType,PetscBool );
2428: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2429: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2430: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2431: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2432: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2433: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2434: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2435: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);
2436: PETSC_EXTERN PetscErrorCode PetscWebServe(MPI_Comm,int);

2438: /*
2439:    In binary files variables are stored using the following lengths,
2440:   regardless of how they are stored in memory on any one particular
2441:   machine. Use these rather then sizeof() in computing sizes for
2442:   PetscBinarySeek().
2443: */
2444: #define PETSC_BINARY_INT_SIZE   (32/8)
2445: #define PETSC_BINARY_FLOAT_SIZE  (32/8)
2446: #define PETSC_BINARY_CHAR_SIZE  (8/8)
2447: #define PETSC_BINARY_SHORT_SIZE  (16/8)
2448: #define PETSC_BINARY_DOUBLE_SIZE  (64/8)
2449: #define PETSC_BINARY_SCALAR_SIZE  sizeof(PetscScalar)

2451: /*E
2452:   PetscBinarySeekType - argument to PetscBinarySeek()

2454:   Level: advanced

2456: .seealso: PetscBinarySeek(), PetscBinarySynchronizedSeek()
2457: E*/
2458: typedef enum {PETSC_BINARY_SEEK_SET = 0,PETSC_BINARY_SEEK_CUR = 1,PETSC_BINARY_SEEK_END = 2} PetscBinarySeekType;
2459: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2460: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2461: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);

2463: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2464: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool );
2465: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2466: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2467: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2468: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);

2470: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2471: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2472: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2473: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2474: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2475: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscInt,const PetscMPIInt*,const void*,PetscInt*,PetscMPIInt**,void*) PetscAttrMPIPointerWithType(6,3);

2477: /*E
2478:     PetscBuildTwoSidedType - algorithm for setting up two-sided communication

2480: $  PETSC_BUILDTWOSIDED_ALLREDUCE - classical algorithm using an MPI_Allreduce with
2481: $      a buffer of length equal to the communicator size. Not memory-scalable due to
2482: $      the large reduction size. Requires only MPI-1.
2483: $  PETSC_BUILDTWOSIDED_IBARRIER - nonblocking algorithm based on MPI_Issend and MPI_Ibarrier.
2484: $      Proved communication-optimal in Hoefler, Siebert, and Lumsdaine (2010). Requires MPI-3.

2486:    Level: developer

2488: .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedSetType(), PetscCommBuildTwoSidedGetType()
2489: E*/
2490: typedef enum {
2491:   PETSC_BUILDTWOSIDED_NOTSET = -1,
2492:   PETSC_BUILDTWOSIDED_ALLREDUCE = 0,
2493:   PETSC_BUILDTWOSIDED_IBARRIER = 1
2494:   /* Updates here must be accompanied by updates in finclude/petscsys.h and the string array in mpits.c */
2495: } PetscBuildTwoSidedType;
2496: PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2497: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2498: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);

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

2502: /*E
2503:   InsertMode - Whether entries are inserted or added into vectors or matrices

2505:   Level: beginner

2507: .seealso: VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2508:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(),
2509:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd()
2510: E*/
2511:  typedef enum {NOT_SET_VALUES, INSERT_VALUES, ADD_VALUES, MAX_VALUES, INSERT_ALL_VALUES, ADD_ALL_VALUES, INSERT_BC_VALUES, ADD_BC_VALUES} InsertMode;

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

2516:     Level: beginner

2518: .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2519:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), ADD_VALUES,
2520:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES

2522: M*/

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

2528:     Level: beginner

2530: .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2531:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), INSERT_VALUES,
2532:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES

2534: M*/

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

2539:     Level: beginner

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

2543: M*/

2545: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);

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

2550: /*S
2551:    PetscSubcomm - Context of MPI subcommunicators, used by PCREDUNDANT

2553:    Level: advanced

2555:    Concepts: communicator, create
2556: S*/
2557: typedef struct _n_PetscSubcomm* PetscSubcomm;

2559: struct _n_PetscSubcomm {
2560:   MPI_Comm    parent;           /* parent communicator */
2561:   MPI_Comm    dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2562:   MPI_Comm    comm;             /* this communicator */
2563:   PetscMPIInt n;                /* num of subcommunicators under the parent communicator */
2564:   PetscMPIInt color;            /* color of processors belong to this communicator */
2565:   PetscMPIInt *subsize;         /* size of subcommunicator[color] */
2566:   PetscSubcommType type;
2567: };

2569: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2570: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2571: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2572: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2573: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2574: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2575: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);

2577: /*S
2578:    PetscSegBuffer - a segmented extendable buffer

2580:    Level: developer

2582: .seealso: PetscSegBufferCreate(), PetscSegBufferGet(), PetscSegBufferExtract(), PetscSegBufferDestroy()
2583: S*/
2584: typedef struct _n_PetscSegBuffer *PetscSegBuffer;
2585: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2586: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2587: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2588: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2589: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2590: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2591: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2592: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);

2595:  * prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as
2596:  * possible. */
2597: PETSC_STATIC_INLINE PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg,PetscInt count,PetscInt *PETSC_RESTRICT *slot) {return PetscSegBufferGet(seg,count,(void**)slot);}

2599: extern PetscSegBuffer PetscCitationsList;
2602: /*@C
2603:       PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code.

2605:      Not Collective - only what is registered on rank 0 of PETSC_COMM_WORLD will be printed

2607:      Input Parameters:
2608: +      cite - the bibtex item, formated to displayed on multiple lines nicely
2609: -      set - a boolean variable initially set to PETSC_FALSE; this is used to insure only a single registration of the citation

2611:    Level: intermediate

2613:      Options Database:
2614: .     -citations [filenmae]   - print out the bibtex entries for the given computation
2615: @*/
2616: PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2617: {
2618:   size_t         len;
2619:   char           *vstring;

2623:   if (set && *set) return(0);
2624:   PetscStrlen(cit,&len);
2625:   PetscSegBufferGet(PetscCitationsList,(PetscInt)len,&vstring);
2626:   PetscMemcpy(vstring,cit,len);
2627:   if (set) *set = PETSC_TRUE;
2628:   return(0);
2629: }

2631: PETSC_EXTERN PetscErrorCode PetscURLShorten(const char[],char[],size_t);
2632: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm,char[],char[],size_t);
2633: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm,const char[],char[],size_t);
2634: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm,const char[],const char []);

2636: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm,char[],char[],size_t);
2637: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm,const char[],char[],char[],size_t);

2639: PETSC_EXTERN PetscErrorCode PetscTextBelt(MPI_Comm,const char[],const char[],PetscBool*);

2641: PETSC_EXTERN PetscErrorCode PetscPullJSONValue(const char[],const char[],char[],size_t,PetscBool*);
2642: PETSC_EXTERN PetscErrorCode PetscPushJSONValue(char[],const char[],const char[],size_t);

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

2648: #endif