Actual source code: petscsys.h

petsc-3.6.1 2015-08-06
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__DEFAULT_SOURCE) && !defined(_DEFAULT_SOURCE)
 29: #define _DEFAULT_SOURCE
 30: #endif
 31: #if defined(PETSC__GNU_SOURCE) && !defined(_GNU_SOURCE)
 32: #define _GNU_SOURCE
 33: #endif
 34: #endif

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

 46: #if defined(__cplusplus)
 47: #  define PETSC_RESTRICT PETSC_CXX_RESTRICT
 48: #else
 49: #  define PETSC_RESTRICT PETSC_C_RESTRICT
 50: #endif

 52: #if defined(__cplusplus)
 53: #  define PETSC_STATIC_INLINE PETSC_CXX_STATIC_INLINE
 54: #else
 55: #  define PETSC_STATIC_INLINE PETSC_C_STATIC_INLINE
 56: #endif

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

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

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

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

 91: /* ========================================================================== */

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

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

110: /*
111:    Perform various sanity checks that the correct mpi.h is being included at compile time.
112:    This usually happens because
113:       * either an unexpected mpi.h is in the default compiler path (i.e. in /usr/include) or
114:       * an extra include path -I/something (which contains the unexpected mpi.h) is being passed to the compiler
115: */
116: #if defined(PETSC_HAVE_MPIUNI)
117: #  if !defined(__MPIUNI_H)
118: #    error "PETSc was configured with --with-mpi=0 but now appears to be compiling using a different mpi.h"
119: #  endif
120: #elif defined(PETSC_HAVE_MPICH_NUMVERSION)
121: #  if !defined(MPICH_NUMVERSION)
122: #    error "PETSc was configured with MPICH but now appears to be compiling using a non-MPICH mpi.h"
123: #  elif MPICH_NUMVERSION != PETSC_HAVE_MPICH_NUMVERSION
124: #    error "PETSc was configured with one MPICH mpi.h version but now appears to be compiling using a different MPICH mpi.h version"
125: #  endif
126: #elif defined(PETSC_HAVE_OMPI_MAJOR_VERSION)
127: #  if !defined(OMPI_MAJOR_VERSION)
128: #    error "PETSc was configured with OpenMPI but now appears to be compiling using a non-OpenMPI mpi.h"
129: #  elif (OMPI_MAJOR_VERSION != PETSC_HAVE_OMPI_MAJOR_VERSION) || (OMPI_MINOR_VERSION != PETSC_HAVE_OMPI_MINOR_VERSION) || (OMPI_RELEASE_VERSION != PETSC_HAVE_OMPI_RELEASE_VERSION)
130: #    error "PETSc was configured with one OpenMPI mpi.h version but now appears to be compiling using a different OpenMPI mpi.h version"
131: #  endif
132: #endif

134: /*
135:     Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler
136:     see the top of mpicxx.h in the MPICH2 distribution.
137: */
138: #include <stdio.h>

140: /* MSMPI on 32bit windows requires this yukky hack - that breaks MPI standard compliance */
141: #if !defined(MPIAPI)
142: #define MPIAPI
143: #endif

145: /* Support for Clang (>=3.2) matching type tag arguments with void* buffer types */
146: #if defined(__has_attribute)
147: #  if __has_attribute(argument_with_type_tag) && __has_attribute(pointer_with_type_tag) && __has_attribute(type_tag_for_datatype)
148: #    define PetscAttrMPIPointerWithType(bufno,typeno) __attribute__((pointer_with_type_tag(MPI,bufno,typeno)))
149: #    define PetscAttrMPITypeTag(type)                 __attribute__((type_tag_for_datatype(MPI,type)))
150: #    define PetscAttrMPITypeTagLayoutCompatible(type) __attribute__((type_tag_for_datatype(MPI,type,layout_compatible)))
151: #  endif
152: #endif
153: #if !defined(PetscAttrMPIPointerWithType)
154: #  define PetscAttrMPIPointerWithType(bufno,typeno)
155: #  define PetscAttrMPITypeTag(type)
156: #  define PetscAttrMPITypeTagLayoutCompatible(type)
157: #endif

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

162:     Level: beginner

164: .seealso: CHKERRQ, SETERRQ
165: M*/
166: typedef int PetscErrorCode;

168: /*MC

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

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

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

178:     Level: developer

180: .seealso: PetscClassIdRegister(), PetscLogEventRegister(), PetscHeaderCreate()
181: M*/
182: typedef int PetscClassId;


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

188:     Level: intermediate

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

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

196: .seealso: PetscBLASInt, PetscInt

198: M*/
199: typedef int PetscMPIInt;

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

204:     Level: intermediate

206: .seealso: PetscOptionsGetEnum(), PetscOptionsEnum(), PetscBagRegisterEnum()
207: M*/
208: typedef enum { ENUM_DUMMY } PetscEnum;
209: PETSC_EXTERN MPI_Datatype MPIU_ENUM PetscAttrMPITypeTag(PetscEnum);

211: #if defined(PETSC_HAVE_STDINT_H)
212: #include <stdint.h>
213: #endif

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

220:    Level: intermediate

222: .seealso: PetscScalar, PetscBLASInt, PetscMPIInt
223: M*/
224: #if defined(PETSC_HAVE_STDINT_H)
225: typedef int64_t Petsc64bitInt;
226: #elif (PETSC_SIZEOF_LONG_LONG == 8)
227: typedef long long Petsc64bitInt;
228: #elif defined(PETSC_HAVE___INT64)
229: typedef __int64 Petsc64bitInt;
230: #else
231: typedef unknown64bit Petsc64bitInt
232: #endif
233: #if defined(PETSC_USE_64BIT_INDICES)
234: typedef Petsc64bitInt PetscInt;
235: #  if defined(PETSC_HAVE_MPI_INT64_T) /* MPI_INT64_T is not guaranteed to be a macro */
236: #    define MPIU_INT MPI_INT64_T
237: #  else
238: #    define MPIU_INT MPI_LONG_LONG_INT
239: #  endif
240: #else
241: typedef int PetscInt;
242: #define MPIU_INT MPI_INT
243: #endif
244: #if defined(PETSC_HAVE_MPI_INT64_T)
245: #  define MPIU_INT64 MPI_INT64_T
246: #else
247: #  define MPIU_INT64 MPI_LONG_LONG_INT
248: #endif


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

254:     Level: intermediate

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

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

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

268:         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

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

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

275: .seealso: PetscMPIInt, PetscInt

277: M*/
278: #if defined(PETSC_HAVE_64BIT_BLAS_INDICES)
279: typedef Petsc64bitInt PetscBLASInt;
280: #else
281: typedef int PetscBLASInt;
282: #endif

284: /*EC

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

288:     Level: advanced

290: .seealso: PetscObjectSetPrecision()
291: E*/
292: typedef enum { PETSC_PRECISION_SINGLE=4,PETSC_PRECISION_DOUBLE=8 } PetscPrecision;
293: PETSC_EXTERN const char *PetscPrecisions[];

295: /*
296:     For the rare cases when one needs to send a size_t object with MPI
297: */
298: #if (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_INT)
299: #define MPIU_SIZE_T MPI_UNSIGNED
300: #elif  (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_LONG)
301: #define MPIU_SIZE_T MPI_UNSIGNED_LONG
302: #elif  (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_LONG_LONG)
303: #define MPIU_SIZE_T MPI_UNSIGNED_LONG_LONG
304: #else
305: #error "Unknown size for size_t! Send us a bugreport at petsc-maint@mcs.anl.gov"
306: #endif


309: /*
310:       You can use PETSC_STDOUT as a replacement of stdout. You can also change
311:     the value of PETSC_STDOUT to redirect all standard output elsewhere
312: */
313: PETSC_EXTERN FILE* PETSC_STDOUT;

315: /*
316:       You can use PETSC_STDERR as a replacement of stderr. You can also change
317:     the value of PETSC_STDERR to redirect all standard error elsewhere
318: */
319: PETSC_EXTERN FILE* PETSC_STDERR;

321: /*MC
322:     PetscUnlikely - hints the compiler that the given condition is usually FALSE

324:     Synopsis:
325:     #include <petscsys.h>
326:     PetscBool  PetscUnlikely(PetscBool  cond)

328:     Not Collective

330:     Input Parameters:
331: .   cond - condition or expression

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

336:     Level: advanced

338: .seealso: PetscLikely(), CHKERRQ
339: M*/

341: /*MC
342:     PetscLikely - hints the compiler that the given condition is usually TRUE

344:     Synopsis:
345:     #include <petscsys.h>
346:     PetscBool  PetscUnlikely(PetscBool  cond)

348:     Not Collective

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

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

356:     Level: advanced

358: .seealso: PetscUnlikely()
359: M*/
360: #if defined(PETSC_HAVE_BUILTIN_EXPECT)
361: #  define PetscUnlikely(cond)   __builtin_expect(!!(cond),0)
362: #  define PetscLikely(cond)     __builtin_expect(!!(cond),1)
363: #else
364: #  define PetscUnlikely(cond)   (cond)
365: #  define PetscLikely(cond)     (cond)
366: #endif

368: /*
369:     Declare extern C stuff after including external header files
370: */


373: /*
374:        Basic PETSc constants
375: */

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

380:    Level: beginner

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

385: E*/
386: typedef enum { PETSC_FALSE,PETSC_TRUE } PetscBool;
387: PETSC_EXTERN const char *const PetscBools[];
388: PETSC_EXTERN MPI_Datatype MPIU_BOOL PetscAttrMPITypeTag(PetscBool);

390: /*
391:     Defines some elementary mathematics functions and constants.
392: */
393: #include <petscmath.h>

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

398:    Level: beginner

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

406: E*/
407: typedef enum { PETSC_COPY_VALUES, PETSC_OWN_POINTER, PETSC_USE_POINTER} PetscCopyMode;
408: PETSC_EXTERN const char *const PetscCopyModes[];

410: /*MC
411:     PETSC_FALSE - False value of PetscBool

413:     Level: beginner

415:     Note: Zero integer

417: .seealso: PetscBool , PETSC_TRUE
418: M*/

420: /*MC
421:     PETSC_TRUE - True value of PetscBool

423:     Level: beginner

425:     Note: Nonzero integer

427: .seealso: PetscBool , PETSC_FALSE
428: M*/

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

433:    Level: beginner

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

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

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

443: M*/
444: #define PETSC_NULL           NULL

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

449:    Level: beginner

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

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

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

459: M*/
460: #define PETSC_IGNORE         NULL

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

466:    Level: beginner

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

470: M*/
471: #define PETSC_DECIDE  -1

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

477:    Level: beginner


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

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

485: M*/
486: #define PETSC_DETERMINE PETSC_DECIDE

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

492:    Level: beginner

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

496: .seealso: PETSC_DECIDE, PETSC_IGNORE, PETSC_DETERMINE

498: M*/
499: #define PETSC_DEFAULT  -2

501: /*MC
502:     PETSC_COMM_WORLD - the equivalent of the MPI_COMM_WORLD communicator which represents
503:            all the processs that PETSc knows about.

505:    Level: beginner

507:    Notes: By default PETSC_COMM_WORLD and MPI_COMM_WORLD are identical unless you wish to
508:           run PETSc on ONLY a subset of MPI_COMM_WORLD. In that case create your new (smaller)
509:           communicator, call it, say comm, and set PETSC_COMM_WORLD = comm BEFORE calling
510:           PetscInitialize(), but after MPI_Init() has been called.

512:           The value of PETSC_COMM_WORLD should never be USED/accessed before PetscInitialize()
513:           is called because it may not have a valid value yet.

515: .seealso: PETSC_COMM_SELF

517: M*/
518: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;

520: /*MC
521:     PETSC_COMM_SELF - This is always MPI_COMM_SELF

523:    Level: beginner

525:    Notes: Do not USE/access or set this variable before PetscInitialize() has been called.

527: .seealso: PETSC_COMM_WORLD

529: M*/
530: #define PETSC_COMM_SELF MPI_COMM_SELF

532: PETSC_EXTERN PetscBool PetscBeganMPI;
533: PETSC_EXTERN PetscBool PetscInitializeCalled;
534: PETSC_EXTERN PetscBool PetscFinalizeCalled;
535: PETSC_EXTERN PetscBool PetscCUSPSynchronize;
536: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;

538: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm));
539: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*);
540: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm*);

542: /*MC
543:    PetscMalloc - Allocates memory, One should use PetscMalloc1() or PetscCalloc1() usually instead of this

545:    Synopsis:
546:     #include <petscsys.h>
547:    PetscErrorCode PetscMalloc(size_t m,void **result)

549:    Not Collective

551:    Input Parameter:
552: .  m - number of bytes to allocate

554:    Output Parameter:
555: .  result - memory allocated

557:    Level: beginner

559:    Notes:
560:    Memory is always allocated at least double aligned

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

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

566:   Concepts: memory allocation

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

571: /*MC
572:    PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment

574:    Synopsis:
575:     #include <petscsys.h>
576:    void *PetscAddrAlign(void *addr)

578:    Not Collective

580:    Input Parameters:
581: .  addr - address to align (any pointer type)

583:    Level: developer

585: .seealso: PetscMallocAlign()

587:   Concepts: memory allocation
588: M*/
589: #define PetscAddrAlign(a) (void*)((((PETSC_UINTPTR_T)(a))+(PETSC_MEMALIGN-1)) & ~(PETSC_MEMALIGN-1))

591: /*MC
592:    PetscMalloc1 - Allocates an array of memory aligned to PETSC_MEMALIGN

594:    Synopsis:
595:     #include <petscsys.h>
596:    PetscErrorCode PetscMalloc1(size_t m1,type **r1)

598:    Not Collective

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

603:    Output Parameter:
604: .  r1 - memory allocated in first chunk

606:    Level: developer

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

610:   Concepts: memory allocation

612: M*/
613: #define PetscMalloc1(m1,r1) PetscMalloc((m1)*sizeof(**(r1)),r1)

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

618:    Synopsis:
619:     #include <petscsys.h>
620:    PetscErrorCode PetscCalloc1(size_t m1,type **r1)

622:    Not Collective

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

627:    Output Parameter:
628: .  r1 - memory allocated in first chunk

630:    Level: developer

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

634:   Concepts: memory allocation

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

639: /*MC
640:    PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN

642:    Synopsis:
643:     #include <petscsys.h>
644:    PetscErrorCode PetscMalloc2(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(), PetscMalloc1(), PetscCalloc2()

660:   Concepts: memory allocation

662: M*/
663: #if !defined(PETSC_USE_MALLOC_COALESCED)
664: #define PetscMalloc2(m1,r1,m2,r2) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)))
665: #else
666: #define PetscMalloc2(m1,r1,m2,r2) ((((m1)+(m2)) ? (*(r2) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(PETSC_MEMALIGN-1),r1)) : 0) \
667:                                    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),0) \
668:                                    || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0))
669: #endif

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

674:    Synopsis:
675:     #include <petscsys.h>
676:    PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)

678:    Not Collective

680:    Input Parameter:
681: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
682: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

684:    Output Parameter:
685: +  r1 - memory allocated in first chunk
686: -  r2 - memory allocated in second chunk

688:    Level: developer

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

692:   Concepts: memory allocation
693: M*/
694: #define PetscCalloc2(m1,r1,m2,r2) (PetscMalloc2((m1),(r1),(m2),(r2)) || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))))

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

699:    Synopsis:
700:     #include <petscsys.h>
701:    PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

703:    Not Collective

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

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

715:    Level: developer

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

719:   Concepts: memory allocation

721: M*/
722: #if !defined(PETSC_USE_MALLOC_COALESCED)
723: #define PetscMalloc3(m1,r1,m2,r2,m3,r3) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)) || PetscMalloc1((m3),(r3)))
724: #else
725: #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) \
726:                                          || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),0) \
727:                                          || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0))
728: #endif

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

733:    Synopsis:
734:     #include <petscsys.h>
735:    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

737:    Not Collective

739:    Input Parameter:
740: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
741: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
742: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

744:    Output Parameter:
745: +  r1 - memory allocated in first chunk
746: .  r2 - memory allocated in second chunk
747: -  r3 - memory allocated in third chunk

749:    Level: developer

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

753:   Concepts: memory allocation
754: M*/
755: #define PetscCalloc3(m1,r1,m2,r2,m3,r3)                                 \
756:   (PetscMalloc3((m1),(r1),(m2),(r2),(m3),(r3))                          \
757:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))))

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

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

766:    Not Collective

768:    Input Parameter:
769: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
770: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
771: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
772: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

774:    Output Parameter:
775: +  r1 - memory allocated in first chunk
776: .  r2 - memory allocated in second chunk
777: .  r3 - memory allocated in third chunk
778: -  r4 - memory allocated in fourth chunk

780:    Level: developer

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

784:   Concepts: memory allocation

786: M*/
787: #if !defined(PETSC_USE_MALLOC_COALESCED)
788: #define PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)) || PetscMalloc1((m3),(r3)) || PetscMalloc1((m4),(r4)))
789: #else
790: #define PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4)                           \
791:   ((((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) \
792:    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),0) \
793:    || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0) || (!(m4) ? (*(r4) = 0,0) : 0))
794: #endif

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

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

803:    Not Collective

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

811:    Output Parameter:
812: +  r1 - memory allocated in first chunk
813: .  r2 - memory allocated in second chunk
814: .  r3 - memory allocated in third chunk
815: -  r4 - memory allocated in fourth chunk

817:    Level: developer

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

821:   Concepts: memory allocation

823: M*/
824: #define PetscCalloc4(m1,r1,m2,r2,m3,r3,m4,r4)                           \
825:   (PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4)                                \
826:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
827:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))))

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

832:    Synopsis:
833:     #include <petscsys.h>
834:    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)

836:    Not Collective

838:    Input Parameter:
839: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
840: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
841: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
842: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
843: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

845:    Output Parameter:
846: +  r1 - memory allocated in first chunk
847: .  r2 - memory allocated in second chunk
848: .  r3 - memory allocated in third chunk
849: .  r4 - memory allocated in fourth chunk
850: -  r5 - memory allocated in fifth chunk

852:    Level: developer

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

856:   Concepts: memory allocation

858: M*/
859: #if !defined(PETSC_USE_MALLOC_COALESCED)
860: #define PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)) || PetscMalloc1((m3),(r3)) || PetscMalloc1((m4),(r4)) || PetscMalloc1((m5),(r5)))
861: #else
862: #define PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5)      \
863:   ((((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) \
864:    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),*(void**)(r5) = PetscAddrAlign(*(r4)+(m4)),0) \
865:    || (!(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))
866: #endif

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

871:    Synopsis:
872:     #include <petscsys.h>
873:    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)

875:    Not Collective

877:    Input Parameter:
878: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
879: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
880: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
881: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
882: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

884:    Output Parameter:
885: +  r1 - memory allocated in first chunk
886: .  r2 - memory allocated in second chunk
887: .  r3 - memory allocated in third chunk
888: .  r4 - memory allocated in fourth chunk
889: -  r5 - memory allocated in fifth chunk

891:    Level: developer

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

895:   Concepts: memory allocation

897: M*/
898: #define PetscCalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5)                     \
899:   (PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5)                          \
900:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
901:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))))

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

906:    Synopsis:
907:     #include <petscsys.h>
908:    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)

910:    Not Collective

912:    Input Parameter:
913: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
914: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
915: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
916: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
917: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
918: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

920:    Output Parameter:
921: +  r1 - memory allocated in first chunk
922: .  r2 - memory allocated in second chunk
923: .  r3 - memory allocated in third chunk
924: .  r4 - memory allocated in fourth chunk
925: .  r5 - memory allocated in fifth chunk
926: -  r6 - memory allocated in sixth chunk

928:    Level: developer

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

932:   Concepts: memory allocation

934: M*/
935: #if !defined(PETSC_USE_MALLOC_COALESCED)
936: #define PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)) || PetscMalloc1((m3),(r3)) || PetscMalloc1((m4),(r4)) || PetscMalloc1((m5),(r5)) || PetscMalloc1((m6),(r6)))
937: #else
938: #define PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) \
939:   ((((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) \
940:    || (*(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) \
941:    || (!(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))
942: #endif

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

947:    Synopsis:
948:     #include <petscsys.h>
949:    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)

951:    Not Collective

953:    Input Parameter:
954: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
955: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
956: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
957: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
958: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
959: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

961:    Output Parameter:
962: +  r1 - memory allocated in first chunk
963: .  r2 - memory allocated in second chunk
964: .  r3 - memory allocated in third chunk
965: .  r4 - memory allocated in fourth chunk
966: .  r5 - memory allocated in fifth chunk
967: -  r6 - memory allocated in sixth chunk

969:    Level: developer

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

973:   Concepts: memory allocation
974: M*/
975: #define PetscCalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6)               \
976:   (PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6)                    \
977:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
978:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))) || PetscMemzero(*(r6),(m6)*sizeof(**(r6))))

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

983:    Synopsis:
984:     #include <petscsys.h>
985:    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)

987:    Not Collective

989:    Input Parameter:
990: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
991: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
992: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
993: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
994: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
995: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
996: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

998:    Output Parameter:
999: +  r1 - memory allocated in first chunk
1000: .  r2 - memory allocated in second chunk
1001: .  r3 - memory allocated in third chunk
1002: .  r4 - memory allocated in fourth chunk
1003: .  r5 - memory allocated in fifth chunk
1004: .  r6 - memory allocated in sixth chunk
1005: -  r7 - memory allocated in seventh chunk

1007:    Level: developer

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

1011:   Concepts: memory allocation

1013: M*/
1014: #if !defined(PETSC_USE_MALLOC_COALESCED)
1015: #define PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)) || PetscMalloc1((m3),(r3)) || PetscMalloc1((m4),(r4)) || PetscMalloc1((m5),(r5)) || PetscMalloc1((m6),(r6)) || PetscMalloc1((m7),(r7)))
1016: #else
1017: #define PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) \
1018:   ((((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) \
1019:    || (*(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) \
1020:    || (!(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))
1021: #endif

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

1026:    Synopsis:
1027:     #include <petscsys.h>
1028:    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)

1030:    Not Collective

1032:    Input Parameter:
1033: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
1034: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
1035: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
1036: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
1037: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
1038: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
1039: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

1041:    Output Parameter:
1042: +  r1 - memory allocated in first chunk
1043: .  r2 - memory allocated in second chunk
1044: .  r3 - memory allocated in third chunk
1045: .  r4 - memory allocated in fourth chunk
1046: .  r5 - memory allocated in fifth chunk
1047: .  r6 - memory allocated in sixth chunk
1048: -  r7 - memory allocated in seventh chunk

1050:    Level: developer

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

1054:   Concepts: memory allocation
1055: M*/
1056: #define PetscCalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7)         \
1057:   (PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7)              \
1058:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
1059:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))) || PetscMemzero(*(r6),(m6)*sizeof(**(r6))) \
1060:    || PetscMemzero(*(r7),(m7)*sizeof(**(r7))))

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

1065:    Synopsis:
1066:     #include <petscsys.h>
1067:    PetscErrorCode PetscNew(type **result)

1069:    Not Collective

1071:    Output Parameter:
1072: .  result - memory allocated, sized to match pointer type

1074:    Level: beginner

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

1078:   Concepts: memory allocation

1080: M*/
1081: #define PetscNew(b)      PetscCalloc1(1,(b))

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

1087:    Synopsis:
1088:     #include <petscsys.h>
1089:    PetscErrorCode PetscNewLog(PetscObject obj,type **result)

1091:    Not Collective

1093:    Input Parameter:
1094: .  obj - object memory is logged to

1096:    Output Parameter:
1097: .  result - memory allocated, sized to match pointer type

1099:    Level: developer

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

1103:   Concepts: memory allocation

1105: M*/
1106: #define PetscNewLog(o,b) (PetscNew((b)) || PetscLogObjectMemory((PetscObject)o,sizeof(**(b))))

1108: /*MC
1109:    PetscFree - Frees memory

1111:    Synopsis:
1112:     #include <petscsys.h>
1113:    PetscErrorCode PetscFree(void *memory)

1115:    Not Collective

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

1120:    Level: beginner

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

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

1128:   Concepts: memory allocation

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

1133: /*MC
1134:    PetscFreeVoid - Frees memory

1136:    Synopsis:
1137:     #include <petscsys.h>
1138:    void PetscFreeVoid(void *memory)

1140:    Not Collective

1142:    Input Parameter:
1143: .   memory - memory to free

1145:    Level: beginner

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

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

1151:   Concepts: memory allocation

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


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

1160:    Synopsis:
1161:     #include <petscsys.h>
1162:    PetscErrorCode PetscFree2(void *memory1,void *memory2)

1164:    Not Collective

1166:    Input Parameter:
1167: +   memory1 - memory to free
1168: -   memory2 - 2nd memory to free

1170:    Level: developer

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

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

1176:   Concepts: memory allocation

1178: M*/
1179: #if !defined(PETSC_USE_MALLOC_COALESCED)
1180: #define PetscFree2(m1,m2)   (PetscFree(m2) || PetscFree(m1))
1181: #else
1182: #define PetscFree2(m1,m2)   ((m1) ? ((m2)=0,PetscFree(m1)) : ((m1)=0,PetscFree(m2)))
1183: #endif

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

1188:    Synopsis:
1189:     #include <petscsys.h>
1190:    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)

1192:    Not Collective

1194:    Input Parameter:
1195: +   memory1 - memory to free
1196: .   memory2 - 2nd memory to free
1197: -   memory3 - 3rd memory to free

1199:    Level: developer

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

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

1205:   Concepts: memory allocation

1207: M*/
1208: #if !defined(PETSC_USE_MALLOC_COALESCED)
1209: #define PetscFree3(m1,m2,m3)   (PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1210: #else
1211: #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))))
1212: #endif

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

1217:    Synopsis:
1218:     #include <petscsys.h>
1219:    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)

1221:    Not Collective

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

1229:    Level: developer

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

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

1235:   Concepts: memory allocation

1237: M*/
1238: #if !defined(PETSC_USE_MALLOC_COALESCED)
1239: #define PetscFree4(m1,m2,m3,m4)   (PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1240: #else
1241: #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)))))
1242: #endif

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

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

1251:    Not Collective

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

1260:    Level: developer

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

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

1266:   Concepts: memory allocation

1268: M*/
1269: #if !defined(PETSC_USE_MALLOC_COALESCED)
1270: #define PetscFree5(m1,m2,m3,m4,m5)   (PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1271: #else
1272: #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)) : \
1273:                                      ((m4) ? ((m5)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m4)) : ((m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m5))))))
1274: #endif


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

1280:    Synopsis:
1281:     #include <petscsys.h>
1282:    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)

1284:    Not Collective

1286:    Input Parameter:
1287: +   m1 - memory to free
1288: .   m2 - 2nd memory to free
1289: .   m3 - 3rd memory to free
1290: .   m4 - 4th memory to free
1291: .   m5 - 5th memory to free
1292: -   m6 - 6th memory to free


1295:    Level: developer

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

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

1301:   Concepts: memory allocation

1303: M*/
1304: #if !defined(PETSC_USE_MALLOC_COALESCED)
1305: #define PetscFree6(m1,m2,m3,m4,m5,m6)   (PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1306: #else
1307: #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)) : \
1308:                                         ((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)) : \
1309:                                         ((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)))))))
1310: #endif

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

1315:    Synopsis:
1316:     #include <petscsys.h>
1317:    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)

1319:    Not Collective

1321:    Input Parameter:
1322: +   m1 - memory to free
1323: .   m2 - 2nd memory to free
1324: .   m3 - 3rd memory to free
1325: .   m4 - 4th memory to free
1326: .   m5 - 5th memory to free
1327: .   m6 - 6th memory to free
1328: -   m7 - 7th memory to free


1331:    Level: developer

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

1335: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1336:           PetscMalloc7()

1338:   Concepts: memory allocation

1340: M*/
1341: #if !defined(PETSC_USE_MALLOC_COALESCED)
1342: #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   (PetscFree(m7) || PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1343: #else
1344: #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)) : \
1345:                                            ((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)) : \
1346:                                            ((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)) : \
1347:                                                    ((m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m7))))))))
1348: #endif

1350: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,int,const char[],const char[],void**);
1351: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1352: PETSC_EXTERN PetscErrorCode PetscMallocSet(PetscErrorCode (*)(size_t,int,const char[],const char[],void**),PetscErrorCode (*)(void*,int,const char[],const char[]));
1353: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);

1355: /*
1356:     PetscLogDouble variables are used to contain double precision numbers
1357:   that are not used in the numerical computations, but rather in logging,
1358:   timing etc.
1359: */
1360: typedef double PetscLogDouble;
1361: #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE

1363: /*
1364:    Routines for tracing memory corruption/bleeding with default PETSc  memory allocation
1365: */
1366: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1367: PETSC_EXTERN PetscErrorCode PetscMallocDumpLog(FILE *);
1368: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1369: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1370: PETSC_EXTERN PetscErrorCode PetscMallocDebug(PetscBool);
1371: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*);
1372: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1373: PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLog(void);
1374: PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLogThreshold(PetscLogDouble);
1375: PETSC_EXTERN PetscErrorCode PetscMallocGetDumpLog(PetscBool*);

1377: /*E
1378:     PetscDataType - Used for handling different basic data types.

1380:    Level: beginner

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

1384: .seealso: PetscBinaryRead(), PetscBinaryWrite(), PetscDataTypeToMPIDataType(),
1385:           PetscDataTypeGetSize()

1387: E*/
1388: typedef enum {PETSC_INT = 0,PETSC_DOUBLE = 1,PETSC_COMPLEX = 2, PETSC_LONG = 3 ,PETSC_SHORT = 4,PETSC_FLOAT = 5,
1389:               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;
1390: PETSC_EXTERN const char *const PetscDataTypes[];

1392: #if defined(PETSC_USE_COMPLEX)
1393: #define  PETSC_SCALAR  PETSC_COMPLEX
1394: #else
1395: #if defined(PETSC_USE_REAL_SINGLE)
1396: #define  PETSC_SCALAR  PETSC_FLOAT
1397: #elif defined(PETSC_USE_REAL___FLOAT128)
1398: #define  PETSC_SCALAR  PETSC___FLOAT128
1399: #else
1400: #define  PETSC_SCALAR  PETSC_DOUBLE
1401: #endif
1402: #endif
1403: #if defined(PETSC_USE_REAL_SINGLE)
1404: #define  PETSC_REAL  PETSC_FLOAT
1405: #elif defined(PETSC_USE_REAL___FLOAT128)
1406: #define  PETSC_REAL  PETSC___FLOAT128
1407: #else
1408: #define  PETSC_REAL  PETSC_DOUBLE
1409: #endif
1410: #define  PETSC_FORTRANADDR  PETSC_LONG

1412: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1413: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1414: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1415: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);

1417: /*
1418:     Basic memory and string operations. These are usually simple wrappers
1419:    around the basic Unix system calls, but a few of them have additional
1420:    functionality and/or error checking.
1421: */
1422: PETSC_EXTERN PetscErrorCode PetscBitMemcpy(void*,PetscInt,const void*,PetscInt,PetscInt,PetscDataType);
1423: PETSC_EXTERN PetscErrorCode PetscMemmove(void*,void *,size_t);
1424: PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1425: PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1426: PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***);
1427: PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1428: PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool  *);
1429: PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool  *);
1430: PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1431: PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1432: PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1433: PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1434: PETSC_EXTERN PetscErrorCode PetscStrncat(char[],const char[],size_t);
1435: PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1436: PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1437: PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1438: PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1439: PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1440: PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1441: PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1442: PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1443: PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*);
1444: PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1445: PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1446: PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1447: PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1448: PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt,const char *const*,char***);
1449: PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt,char***);
1450: PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);

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

1454: /*S
1455:     PetscToken - 'Token' used for managing tokenizing strings

1457:   Level: intermediate

1459: .seealso: PetscTokenCreate(), PetscTokenFind(), PetscTokenDestroy()
1460: S*/
1461: typedef struct _p_PetscToken* PetscToken;

1463: PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1464: PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1465: PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);

1467: PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1468: PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);

1470: /*
1471:    These are  MPI operations for MPI_Allreduce() etc
1472: */
1473: PETSC_EXTERN MPI_Op PetscMaxSum_Op;
1474: #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128)
1475: PETSC_EXTERN MPI_Op MPIU_SUM;
1476: #else
1477: #define MPIU_SUM MPI_SUM
1478: #endif
1479: #if defined(PETSC_USE_REAL___FLOAT128)
1480: PETSC_EXTERN MPI_Op MPIU_MAX;
1481: PETSC_EXTERN MPI_Op MPIU_MIN;
1482: #else
1483: #define MPIU_MAX MPI_MAX
1484: #define MPIU_MIN MPI_MIN
1485: #endif
1486: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);

1488: PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1489: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);

1491: /*S
1492:      PetscObject - any PETSc object, PetscViewer, Mat, Vec, KSP etc

1494:    Level: beginner

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

1498: .seealso:  PetscObjectDestroy(), PetscObjectView(), PetscObjectGetName(), PetscObjectSetName(), PetscObjectReference(), PetscObjectDereference()
1499: S*/
1500: typedef struct _p_PetscObject* PetscObject;

1502: /*MC
1503:     PetscObjectId - unique integer Id for a PetscObject

1505:     Level: developer

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

1509: .seealso: PetscObjectState, PetscObjectGetId()
1510: M*/
1511: typedef Petsc64bitInt PetscObjectId;

1513: /*MC
1514:     PetscObjectState - integer state for a PetscObject

1516:     Level: developer

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

1522: .seealso: PetscObjectId, PetscObjectStateGet(), PetscObjectStateIncrease(), PetscObjectStateSet()
1523: M*/
1524: typedef Petsc64bitInt PetscObjectState;

1526: /*S
1527:      PetscFunctionList - Linked list of functions, possibly stored in dynamic libraries, accessed
1528:       by string name

1530:    Level: advanced

1532: .seealso:  PetscFunctionListAdd(), PetscFunctionListDestroy(), PetscOpFlist
1533: S*/
1534: typedef struct _n_PetscFunctionList *PetscFunctionList;

1536: /*E
1537:   PetscFileMode - Access mode for a file.

1539:   Level: beginner

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

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

1545:   FILE_MODE_APPEND - open a file at end for writing

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

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

1551: .seealso: PetscViewerFileSetMode()
1552: E*/
1553: typedef enum {FILE_MODE_READ, FILE_MODE_WRITE, FILE_MODE_APPEND, FILE_MODE_UPDATE, FILE_MODE_APPEND_UPDATE} PetscFileMode;
1554: extern const char *const PetscFileModes[];

1556: /*
1557:     Defines PETSc error handling.
1558: */
1559: #include <petscerror.h>

1561: #define PETSC_SMALLEST_CLASSID  1211211
1562: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1563: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1564: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);

1566: /*
1567:    Routines that get memory usage information from the OS
1568: */
1569: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1570: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1571: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1572: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);

1574: PETSC_EXTERN PetscErrorCode PetscInfoAllow(PetscBool ,const char []);
1575: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);

1577: /*
1578:    Initialization of PETSc
1579: */
1580: PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1581: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1582: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1583: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1584: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1585: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1586: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1587: PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1588: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1589: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);

1591: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1592: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);

1594: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1595: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1596: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1597: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);

1599: /*
1600:      These are so that in extern C code we can caste function pointers to non-extern C
1601:    function pointers. Since the regular C++ code expects its function pointers to be C++
1602: */
1603: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1604: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1605: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);

1607: /*
1608:     Functions that can act on any PETSc object.
1609: */
1610: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1611: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1612: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1613: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1614: PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1615: PETSC_EXTERN PetscErrorCode PetscObjectSetPrecision(PetscObject,PetscPrecision);
1616: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1617: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1618: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1619: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1620: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1621: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1622: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1623: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1624: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1625: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt *);
1626: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1627: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1628: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject *);
1629: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void));
1630: #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d))
1631: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1632: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1633: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt *);
1634: PETSC_EXTERN PetscErrorCode PetscObjectAddOptionsHandler(PetscObject,PetscErrorCode (*)(PetscObject,void*),PetscErrorCode (*)(PetscObject,void*),void*);
1635: PETSC_EXTERN PetscErrorCode PetscObjectProcessOptionsHandlers(PetscObject);
1636: PETSC_EXTERN PetscErrorCode PetscObjectDestroyOptionsHandlers(PetscObject);
1637: PETSC_EXTERN PetscErrorCode PetscObjectsListGetGlobalNumbering(MPI_Comm,PetscInt,PetscObject*,PetscInt*,PetscInt*);

1639: #include <petscviewertypes.h>
1640: #include <petscoptions.h>

1642: PETSC_EXTERN PetscErrorCode PetscMemoryShowUsage(PetscViewer,const char[]);
1643: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1644: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1645: #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1646: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1647: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1648: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1649: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1650: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1651: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1652: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1653: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1654: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,PetscObject,const char[]);
1655: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1656: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1657: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1658: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1659: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);

1661: #if defined(PETSC_HAVE_SAWS)
1662: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1663: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1664: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1665: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1666: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1667: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1668: PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1669: PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1670: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1671: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);

1673: #else
1674: #define PetscSAWsBlock()                        0
1675: #define PetscObjectSAWsViewOff(obj)             0
1676: #define PetscObjectSAWsSetBlock(obj,flg)        0
1677: #define PetscObjectSAWsBlock(obj)               0
1678: #define PetscObjectSAWsGrantAccess(obj)         0
1679: #define PetscObjectSAWsTakeAccess(obj)          0
1680: #define PetscStackViewSAWs()                    0
1681: #define PetscStackSAWsViewOff()                 0
1682: #define PetscStackSAWsTakeAccess()
1683: #define PetscStackSAWsGrantAccess()

1685: #endif

1687: typedef void* PetscDLHandle;
1688: typedef enum {PETSC_DL_DECIDE=0,PETSC_DL_NOW=1,PETSC_DL_LOCAL=2} PetscDLMode;
1689: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1690: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1691: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);


1694: #if defined(PETSC_USE_DEBUG)
1695: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1696: #endif
1697: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);

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

1702:    Level: developer

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

1706: .seealso:  PetscObjectListAdd(), PetscObjectListDestroy(), PetscObjectListFind(), PetscObjectCompose(), PetscObjectQuery(), PetscFunctionList
1707: S*/
1708: typedef struct _n_PetscObjectList *PetscObjectList;

1710: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1711: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1712: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1713: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1714: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1715: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);

1717: /*
1718:     Dynamic library lists. Lists of names of routines in objects or in dynamic
1719:   link libraries that will be loaded as needed.
1720: */

1722: #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1723: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1724: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1725: #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1726: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1727: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[]);
1728: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1729: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1730: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);

1732: /*S
1733:      PetscDLLibrary - Linked list of dynamics libraries to search for functions

1735:    Level: advanced

1737: .seealso:  PetscDLLibraryOpen()
1738: S*/
1739: typedef struct _n_PetscDLLibrary *PetscDLLibrary;
1740: PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1741: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1742: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1743: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1744: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1745: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1746: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1747: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);

1749: /*
1750:      Useful utility routines
1751: */
1752: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1753: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1754: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1755: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1756: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1757: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);

1759: /*
1760:     PetscNot - negates a logical type value and returns result as a PetscBool

1762:     Notes: This is useful in cases like
1763: $     int        *a;
1764: $     PetscBool  flag = PetscNot(a)
1765:      where !a does not return a PetscBool  because we cannot provide a cast from int to PetscBool  in C.
1766: */
1767: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)

1769: /*MC
1770:     PetscHelpPrintf - Prints help messages.

1772:    Synopsis:
1773:     #include <petscsys.h>
1774:      PetscErrorCode (*PetscHelpPrintf)(const char format[],...);

1776:     Not Collective

1778:     Input Parameters:
1779: .   format - the usual printf() format string

1781:    Level: developer

1783:     Fortran Note:
1784:     This routine is not supported in Fortran.

1786:     Concepts: help messages^printing
1787:     Concepts: printing^help messages

1789: .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1790: M*/
1791: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);

1793: /*
1794:      Defines PETSc profiling.
1795: */
1796: #include <petsclog.h>



1800: /*
1801:       Simple PETSc parallel IO for ASCII printing
1802: */
1803: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1804: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1805: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1806: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1807: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1808: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1809: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);

1811: /* These are used internally by PETSc ASCII IO routines*/
1812: #include <stdarg.h>
1813: PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list);
1814: PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE*,const char[],va_list);
1815: PETSC_EXTERN PetscErrorCode PetscVFPrintfDefault(FILE*,const char[],va_list);

1817: #if defined(PETSC_HAVE_MATLAB_ENGINE)
1818: PETSC_EXTERN PetscErrorCode PetscVFPrintf_Matlab(FILE*,const char[],va_list);
1819: #endif

1821: #if defined(PETSC_HAVE_CLOSURES)
1822: PETSC_EXTERN PetscErrorCode PetscVFPrintfSetClosure(int (^)(const char*));
1823: #endif

1825: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1826: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1827: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);

1829: #if defined(PETSC_HAVE_POPEN)
1830: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1831: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*,int*);
1832: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1833: #endif

1835: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1836: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1837: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*);
1838: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1839: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1840: PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1841: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);

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

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

1848:    Level: advanced

1850: .seealso:  PetscObject, PetscContainerCreate()
1851: S*/
1852: PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1853: typedef struct _p_PetscContainer*  PetscContainer;
1854: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void **);
1855: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void *);
1856: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1857: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer *);
1858: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void*));

1860: /*
1861:    For use in debuggers
1862: */
1863: PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1864: PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1865: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1866: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1867: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);

1869: #include <stddef.h>
1870: #include <string.h>             /* for memcpy, memset */
1871: #if defined(PETSC_HAVE_STDLIB_H)
1872: #include <stdlib.h>
1873: #endif

1875: #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1876: #include <xmmintrin.h>
1877: #endif

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

1886:    Not Collective

1888:    Input Parameters:
1889: +  b - pointer to initial memory space
1890: -  n - length (in bytes) of space to copy

1892:    Output Parameter:
1893: .  a - pointer to copy space

1895:    Level: intermediate

1897:    Compile Option:
1898:     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1899:                                   for memory copies on double precision values.
1900:     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1901:                                   for memory copies on double precision values.
1902:     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1903:                                   for memory copies on double precision values.

1905:    Note:
1906:    This routine is analogous to memcpy().

1908:    Developer Note: this is inlined for fastest performance

1910:   Concepts: memory^copying
1911:   Concepts: copying^memory

1913: .seealso: PetscMemmove()

1915: @*/
1916: PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n)
1917: {
1918: #if defined(PETSC_USE_DEBUG)
1919:   unsigned long al = (unsigned long) a,bl = (unsigned long) b;
1920:   unsigned long nl = (unsigned long) n;
1922:   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1923:   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1924: #else
1926: #endif
1927:   if (a != b && n > 0) {
1928: #if defined(PETSC_USE_DEBUG)
1929:     if ((al > bl && (al - bl) < nl) || (bl - al) < nl)  SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1930:               or make sure your copy regions and lengths are correct. \n\
1931:               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1932: #endif
1933: #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1934:    if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1935:       size_t len = n/sizeof(PetscScalar);
1936: #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1937:       PetscBLASInt   one = 1,blen;
1939:       PetscBLASIntCast(len,&blen);
1940:       PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one));
1941: #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1942:       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1943: #else
1944:       size_t      i;
1945:       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1946:       for (i=0; i<len; i++) y[i] = x[i];
1947: #endif
1948:     } else {
1949:       memcpy((char*)(a),(char*)(b),n);
1950:     }
1951: #else
1952:     memcpy((char*)(a),(char*)(b),n);
1953: #endif
1954:   }
1955:   return(0);
1956: }

1958: /*@C
1959:    PetscMemzero - Zeros the specified memory.

1961:    Not Collective

1963:    Input Parameters:
1964: +  a - pointer to beginning memory location
1965: -  n - length (in bytes) of memory to initialize

1967:    Level: intermediate

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

1973:    Developer Note: this is inlined for fastest performance

1975:    Concepts: memory^zeroing
1976:    Concepts: zeroing^memory

1978: .seealso: PetscMemcpy()
1979: @*/
1980: PETSC_STATIC_INLINE PetscErrorCode  PetscMemzero(void *a,size_t n)
1981: {
1982:   if (n > 0) {
1983: #if defined(PETSC_USE_DEBUG)
1984:     if (!a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer");
1985: #endif
1986: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1987:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1988:       size_t      i,len = n/sizeof(PetscScalar);
1989:       PetscScalar *x = (PetscScalar*)a;
1990:       for (i=0; i<len; i++) x[i] = 0.0;
1991:     } else {
1992: #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1993:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1994:       PetscInt len = n/sizeof(PetscScalar);
1995:       fortranzero_(&len,(PetscScalar*)a);
1996:     } else {
1997: #endif
1998: #if defined(PETSC_PREFER_BZERO)
1999:       bzero((char *)a,n);
2000: #else
2001:       memset((char*)a,0,n);
2002: #endif
2003: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
2004:     }
2005: #endif
2006:   }
2007:   return 0;
2008: }

2010: /*MC
2011:    PetscPrefetchBlock - Prefetches a block of memory

2013:    Synopsis:
2014:     #include <petscsys.h>
2015:     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)

2017:    Not Collective

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

2025:    Level: developer

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

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

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

2040:    Concepts: memory
2041: M*/
2042: #define PetscPrefetchBlock(a,n,rw,t) do {                               \
2043:     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
2044:     for ( ; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
2045:   } while (0)

2047: /*
2048:       Determine if some of the kernel computation routines use
2049:    Fortran (rather than C) for the numerical calculations. On some machines
2050:    and compilers (like complex numbers) the Fortran version of the routines
2051:    is faster than the C/C++ versions. The flag --with-fortran-kernels
2052:    should be used with ./configure to turn these on.
2053: */
2054: #if defined(PETSC_USE_FORTRAN_KERNELS)

2056: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
2057: #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
2058: #endif

2060: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
2061: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
2062: #endif

2064: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
2065: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
2066: #endif

2068: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
2069: #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
2070: #endif

2072: #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
2073: #define PETSC_USE_FORTRAN_KERNEL_NORM
2074: #endif

2076: #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
2077: #define PETSC_USE_FORTRAN_KERNEL_MAXPY
2078: #endif

2080: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
2081: #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
2082: #endif

2084: #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
2085: #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
2086: #endif

2088: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
2089: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
2090: #endif

2092: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
2093: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
2094: #endif

2096: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
2097: #define PETSC_USE_FORTRAN_KERNEL_MDOT
2098: #endif

2100: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
2101: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
2102: #endif

2104: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
2105: #define PETSC_USE_FORTRAN_KERNEL_AYPX
2106: #endif

2108: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
2109: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
2110: #endif

2112: #endif

2114: /*
2115:     Macros for indicating code that should be compiled with a C interface,
2116:    rather than a C++ interface. Any routines that are dynamically loaded
2117:    (such as the PCCreate_XXX() routines) must be wrapped so that the name
2118:    mangler does not change the functions symbol name. This just hides the
2119:    ugly extern "C" {} wrappers.
2120: */
2121: #if defined(__cplusplus)
2122: #define EXTERN_C_BEGIN extern "C" {
2123: #define EXTERN_C_END }
2124: #else
2125: #define EXTERN_C_BEGIN
2126: #define EXTERN_C_END
2127: #endif

2129: /* --------------------------------------------------------------------*/

2131: /*MC
2132:     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
2133:         communication

2135:    Level: beginner

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

2139: .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
2140: M*/

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


2148:    Level: beginner

2150: .seealso: PetscReal, PassiveReal, PassiveScalar, MPIU_SCALAR, PetscInt
2151: M*/

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

2156:    Synopsis:
2157:    #include <petscsys.h>
2158:    PetscComplex number = 1. + 2.*PETSC_i;

2160:    Level: beginner

2162:    Note:
2163:    Complex numbers are automatically available if PETSc was able to find a working complex implementation

2165: .seealso: PetscReal, PetscComplex, MPIU_COMPLEX, PetscInt, PETSC_i
2166: M*/

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

2171:    Level: beginner

2173: .seealso: PetscScalar, PassiveReal, PassiveScalar
2174: M*/

2176: /*MC
2177:     PassiveScalar - PETSc type that represents a PetscScalar
2178:    Level: beginner

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

2183: .seealso: PetscReal, PassiveReal, PetscScalar
2184: M*/

2186: /*MC
2187:     PassiveReal - PETSc type that represents a PetscReal

2189:    Level: beginner

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

2194: .seealso: PetscScalar, PetscReal, PassiveScalar
2195: M*/

2197: /*MC
2198:     MPIU_SCALAR - MPI datatype corresponding to PetscScalar

2200:    Level: beginner

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

2205: .seealso: PetscReal, PassiveReal, PassiveScalar, PetscScalar, MPIU_INT
2206: M*/

2208: #if defined(PETSC_HAVE_MPIIO)
2209: #if !defined(PETSC_WORDS_BIGENDIAN)
2210: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2211: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2212: #else
2213: #define MPIU_File_write_all(a,b,c,d,e) MPI_File_write_all(a,b,c,d,e)
2214: #define MPIU_File_read_all(a,b,c,d,e) MPI_File_read_all(a,b,c,d,e)
2215: #endif
2216: #endif

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

2220: /* Limit MPI to 32-bits */
2221: #define PETSC_MPI_INT_MAX  2147483647
2222: #define PETSC_MPI_INT_MIN -2147483647
2223: /* Limit BLAS to 32-bits */
2224: #define PETSC_BLAS_INT_MAX  2147483647
2225: #define PETSC_BLAS_INT_MIN -2147483647

2229: /*@C
2230:     PetscBLASIntCast - casts a PetscInt (which may be 64 bits in size) to a PetscBLASInt (which may be 32 bits in size), generates an
2231:          error if the PetscBLASInt is not large enough to hold the number.

2233:    Not Collective

2235:    Input Parameter:
2236: .     a - the PetscInt value

2238:    Output Parameter:
2239: .     b - the resulting PetscBLASInt value

2241:    Level: advanced

2243: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast()
2244: @*/
2245: PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
2246: {
2248:   *b =  (PetscBLASInt)(a);
2249: #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
2250:   if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK");
2251: #endif
2252:   return(0);
2253: }

2257: /*@C
2258:     PetscMPIIntCast - casts a PetscInt (which may be 64 bits in size) to a PetscMPIInt (which may be 32 bits in size), generates an
2259:          error if the PetscMPIInt is not large enough to hold the number.

2261:    Not Collective

2263:    Input Parameter:
2264: .     a - the PetscInt value

2266:    Output Parameter:
2267: .     b - the resulting PetscMPIInt value

2269:    Level: advanced

2271: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast()
2272: @*/
2273: PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
2274: {
2276:   *b =  (PetscMPIInt)(a);
2277: #if defined(PETSC_USE_64BIT_INDICES)
2278:   if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for MPI");
2279: #endif
2280:   return(0);
2281: }


2284: /*
2285:      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2286: */
2287: #if defined(hz)
2288: #undef hz
2289: #endif

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


2294: #if defined(PETSC_HAVE_LIMITS_H)
2295: #include <limits.h>
2296: #endif
2297: #if defined(PETSC_HAVE_SYS_PARAM_H)
2298: #include <sys/param.h>
2299: #endif
2300: #if defined(PETSC_HAVE_SYS_TYPES_H)
2301: #include <sys/types.h>
2302: #endif
2303: #if defined(MAXPATHLEN)
2304: #  define PETSC_MAX_PATH_LEN     MAXPATHLEN
2305: #elif defined(MAX_PATH)
2306: #  define PETSC_MAX_PATH_LEN     MAX_PATH
2307: #elif defined(_MAX_PATH)
2308: #  define PETSC_MAX_PATH_LEN     _MAX_PATH
2309: #else
2310: #  define PETSC_MAX_PATH_LEN     4096
2311: #endif

2313: /*MC

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

2317: $    1) classic Fortran 77 style
2318: $#include "petsc/finclude/petscXXX.h" to work with material from the XXX component of PETSc
2319: $       XXX variablename
2320: $      You cannot use this approach if you wish to use the Fortran 90 specific PETSc routines
2321: $      which end in F90; such as VecGetArrayF90()
2322: $
2323: $    2) classic Fortran 90 style
2324: $#include "petsc/finclude/petscXXX.h"
2325: $#include "petsc/finclude/petscXXX.h90" to work with material from the XXX component of PETSc
2326: $       XXX variablename
2327: $
2328: $    3) Using Fortran modules
2329: $#include "petsc/finclude/petscXXXdef.h"
2330: $         use petscXXXX
2331: $       XXX variablename
2332: $
2333: $    4) Use Fortran modules and Fortran data types for PETSc types
2334: $#include "petsc/finclude/petscXXXdef.h"
2335: $         use petscXXXX
2336: $       type(XXX) variablename
2337: $      To use this approach you must ./configure PETSc with the additional
2338: $      option --with-fortran-datatypes You cannot use the type(XXX) declaration approach without using Fortran modules

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

2342: $    3a) skip the #include BUT you cannot use any PETSc data type names like Vec, Mat, PetscInt, PetscErrorCode etc
2343: $        and you must declare the variables as integer, for example
2344: $        integer variablename
2345: $
2346: $    4a) skip the #include, you use the object types like type(Vec) type(Mat) but cannot use the data type
2347: $        names like PetscErrorCode, PetscInt etc. again for those you must use integer

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

2352:    Fortran type checking with interfaces is strick, this means you cannot pass a scalar value when an array value
2353: is expected (even though it is legal Fortran). For example when setting a single value in a matrix with MatSetValues()
2354: you cannot have something like
2355: $      PetscInt row,col
2356: $      PetscScalar val
2357: $        ...
2358: $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)
2359: You must instead have
2360: $      PetscInt row(1),col(1)
2361: $      PetscScalar val(1)
2362: $        ...
2363: $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)


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

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

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

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

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

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

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

2388:     Level: beginner

2390: M*/

2392: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2393: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2394: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2395: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2396: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2397: PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2398: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);

2400: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2401: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2402: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2403: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2404: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2405: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2406: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt*,PetscInt*,PetscInt*);
2407: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2408: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2409: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2410: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2411: PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt,PetscInt[],void*,size_t,void*);
2412: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2413: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2414: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2415: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2416: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2417: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2418: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt*,const PetscInt*,PetscInt,const PetscInt*,const PetscInt*,PetscInt*,PetscInt**,PetscInt**);
2419: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt,const PetscInt*,PetscInt,const PetscInt*,PetscInt*,PetscInt**);

2421: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2422: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);

2424: /*J
2425:     PetscRandomType - String with the name of a PETSc randomizer

2427:    Level: beginner

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

2432: .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2433: J*/
2434: typedef const char* PetscRandomType;
2435: #define PETSCRAND       "rand"
2436: #define PETSCRAND48     "rand48"
2437: #define PETSCSPRNG      "sprng"

2439: /* Logging support */
2440: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;

2442: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);

2444: /*S
2445:      PetscRandom - Abstract PETSc object that manages generating random numbers

2447:    Level: intermediate

2449:   Concepts: random numbers

2451: .seealso:  PetscRandomCreate(), PetscRandomGetValue(), PetscRandomType
2452: S*/
2453: typedef struct _p_PetscRandom*   PetscRandom;

2455: /* Dynamic creation and loading functions */
2456: PETSC_EXTERN PetscFunctionList PetscRandomList;

2458: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2459: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2460: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2461: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType*);
2462: PETSC_STATIC_INLINE PetscErrorCode PetscRandomViewFromOptions(PetscRandom A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
2463: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);

2465: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2466: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2467: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2468: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2469: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2470: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2471: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2472: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2473: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);

2475: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2476: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2477: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2478: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2479: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2480: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2481: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);

2483: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscDataType);
2484: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscDataType);
2485: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,void*,PetscInt,PetscDataType,PetscBool );
2486: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,void*,PetscInt,PetscDataType,PetscBool );
2487: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2488: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2489: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2490: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2491: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2492: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2493: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2494: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);

2496: /*
2497:    In binary files variables are stored using the following lengths,
2498:   regardless of how they are stored in memory on any one particular
2499:   machine. Use these rather then sizeof() in computing sizes for
2500:   PetscBinarySeek().
2501: */
2502: #define PETSC_BINARY_INT_SIZE   (32/8)
2503: #define PETSC_BINARY_FLOAT_SIZE  (32/8)
2504: #define PETSC_BINARY_CHAR_SIZE  (8/8)
2505: #define PETSC_BINARY_SHORT_SIZE  (16/8)
2506: #define PETSC_BINARY_DOUBLE_SIZE  (64/8)
2507: #define PETSC_BINARY_SCALAR_SIZE  sizeof(PetscScalar)

2509: /*E
2510:   PetscBinarySeekType - argument to PetscBinarySeek()

2512:   Level: advanced

2514: .seealso: PetscBinarySeek(), PetscBinarySynchronizedSeek()
2515: E*/
2516: typedef enum {PETSC_BINARY_SEEK_SET = 0,PETSC_BINARY_SEEK_CUR = 1,PETSC_BINARY_SEEK_END = 2} PetscBinarySeekType;
2517: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2518: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2519: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);

2521: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2522: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool );
2523: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2524: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2525: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2526: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);

2528: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2529: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2530: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2531: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2532: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2533: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscInt,const PetscMPIInt*,const void*,PetscInt*,PetscMPIInt**,void*) PetscAttrMPIPointerWithType(6,3);

2535: /*E
2536:     PetscBuildTwoSidedType - algorithm for setting up two-sided communication

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

2544:    Level: developer

2546: .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedSetType(), PetscCommBuildTwoSidedGetType()
2547: E*/
2548: typedef enum {
2549:   PETSC_BUILDTWOSIDED_NOTSET = -1,
2550:   PETSC_BUILDTWOSIDED_ALLREDUCE = 0,
2551:   PETSC_BUILDTWOSIDED_IBARRIER = 1
2552:   /* Updates here must be accompanied by updates in petsc/finclude/petscsys.h and the string array in mpits.c */
2553: } PetscBuildTwoSidedType;
2554: PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2555: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2556: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);

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

2560: /*E
2561:   InsertMode - Whether entries are inserted or added into vectors or matrices

2563:   Level: beginner

2565: .seealso: VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2566:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(),
2567:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd()
2568: E*/
2569:  typedef enum {NOT_SET_VALUES, INSERT_VALUES, ADD_VALUES, MAX_VALUES, INSERT_ALL_VALUES, ADD_ALL_VALUES, INSERT_BC_VALUES, ADD_BC_VALUES} InsertMode;

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

2574:     Level: beginner

2576: .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2577:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), ADD_VALUES,
2578:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES

2580: M*/

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

2586:     Level: beginner

2588: .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2589:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), INSERT_VALUES,
2590:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES

2592: M*/

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

2597:     Level: beginner

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

2601: M*/

2603: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);

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

2608: /*S
2609:    PetscSubcomm - A decomposition of an MPI communicator into subcommunicators

2611:    Notes: After a call to PetscSubcommSetType(), PetscSubcommSetTypeGeneral(), or PetscSubcommSetFromOptions() one may call
2612: $     PetscSubcommChild() returns the associated subcommunicator on this process
2613: $     PetscSubcommContiguousParent() returns a parent communitor but with all child of the same subcommunicator having contiquous rank

2615:    Sample Usage:
2616:        PetscSubcommCreate()
2617:        PetscSubcommSetNumber()
2618:        PetscSubcommSetType(PETSC_SUBCOMM_INTERLACED);
2619:        ccomm = PetscSubcommChild()
2620:        PetscSubcommDestroy()

2622:    Level: advanced

2624:    Concepts: communicator, create

2626:    Notes:
2627: $   PETSC_SUBCOMM_GENERAL - similar to MPI_Comm_split() each process sets the new communicator (color) they will belong to and the order within that communicator
2628: $   PETSC_SUBCOMM_CONTIGUOUS - each new communicator contains a set of process with contiquous ranks in the original MPI communicator
2629: $   PETSC_SUBCOMM_INTERLACED - each new communictor contains a set of processes equally far apart in rank from the others in that new communicator

2631:    Examaple: Consider a communicator with six processes split into 3 subcommunicators.
2632: $     PETSC_SUBCOMM_CONTIGUOUS - the first communicator contains rank 0,1  the second rank 2,3 and the third rank 4,5 in the original ordering of the original communicator
2633: $     PETSC_SUBCOMM_INTERLACED - the first communicator contains rank 0,3, the second 1,4 and the third 2,5

2635:    Developer Notes: This is used in objects such as PCREDUNDANT() to manage the subcommunicators on which the redundant computations
2636:       are performed.


2639: .seealso: PetscSubcommCreate(), PetscSubcommSetNumber(), PetscSubcommSetType(), PetscSubcommView(), PetscSubcommSetFromOptions()

2641: S*/
2642: typedef struct _n_PetscSubcomm* PetscSubcomm;

2644: struct _n_PetscSubcomm {
2645:   MPI_Comm         parent;           /* parent communicator */
2646:   MPI_Comm         dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2647:   MPI_Comm         child;            /* the sub-communicator */
2648:   PetscMPIInt      n;                /* num of subcommunicators under the parent communicator */
2649:   PetscMPIInt      color;            /* color of processors belong to this communicator */
2650:   PetscMPIInt      *subsize;         /* size of subcommunicator[color] */
2651:   PetscSubcommType type;
2652: };

2654: PETSC_STATIC_INLINE MPI_Comm PetscSubcommChild(PetscSubcomm scomm) {return scomm->child;}
2655: PETSC_STATIC_INLINE MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm) {return scomm->dupparent;}
2656: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2657: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2658: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2659: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2660: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2661: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2662: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);

2664: /*S
2665:    PetscSegBuffer - a segmented extendable buffer

2667:    Level: developer

2669: .seealso: PetscSegBufferCreate(), PetscSegBufferGet(), PetscSegBufferExtract(), PetscSegBufferDestroy()
2670: S*/
2671: typedef struct _n_PetscSegBuffer *PetscSegBuffer;
2672: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2673: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2674: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2675: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2676: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2677: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2678: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2679: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);

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

2686: PETSC_EXTERN PetscSegBuffer PetscCitationsList;
2689: /*@C
2690:       PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code.

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

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

2698:    Level: intermediate

2700:      Options Database:
2701: .     -citations [filenmae]   - print out the bibtex entries for the given computation
2702: @*/
2703: PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2704: {
2705:   size_t         len;
2706:   char           *vstring;

2710:   if (set && *set) return(0);
2711:   PetscStrlen(cit,&len);
2712:   PetscSegBufferGet(PetscCitationsList,len,&vstring);
2713:   PetscMemcpy(vstring,cit,len);
2714:   if (set) *set = PETSC_TRUE;
2715:   return(0);
2716: }

2718: PETSC_EXTERN PetscErrorCode PetscURLShorten(const char[],char[],size_t);
2719: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm,char[],char[],size_t);
2720: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm,const char[],char[],size_t);
2721: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm,const char[],const char []);

2723: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm,char[],char[],size_t);
2724: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm,const char[],char[],char[],size_t);

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

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

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

2735: #endif