Actual source code: petscsys.h

petsc-3.5.4 2015-05-23
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()

512: .seealso: PETSC_COMM_SELF

514: M*/
515: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;

517: /*MC
518:     PETSC_COMM_SELF - This is always MPI_COMM_SELF

520:    Level: beginner

522: .seealso: PETSC_COMM_WORLD

524: M*/
525: #define PETSC_COMM_SELF MPI_COMM_SELF

527: PETSC_EXTERN PetscBool PetscBeganMPI;
528: PETSC_EXTERN PetscBool PetscInitializeCalled;
529: PETSC_EXTERN PetscBool PetscFinalizeCalled;
530: PETSC_EXTERN PetscBool PetscCUSPSynchronize;
531: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;

533: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm));
534: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*);
535: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm*);

537: /*MC
538:    PetscMalloc - Allocates memory

540:    Synopsis:
541:     #include <petscsys.h>
542:    PetscErrorCode PetscMalloc(size_t m,void **result)

544:    Not Collective

546:    Input Parameter:
547: .  m - number of bytes to allocate

549:    Output Parameter:
550: .  result - memory allocated

552:    Level: beginner

554:    Notes:
555:    Memory is always allocated at least double aligned

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

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

561:   Concepts: memory allocation

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

566: /*MC
567:    PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment

569:    Synopsis:
570:     #include <petscsys.h>
571:    void *PetscAddrAlign(void *addr)

573:    Not Collective

575:    Input Parameters:
576: .  addr - address to align (any pointer type)

578:    Level: developer

580: .seealso: PetscMallocAlign()

582:   Concepts: memory allocation
583: M*/
584: #define PetscAddrAlign(a) (void*)((((PETSC_UINTPTR_T)(a))+(PETSC_MEMALIGN-1)) & ~(PETSC_MEMALIGN-1))

586: /*MC
587:    PetscMalloc1 - Allocates an array of memory aligned to PETSC_MEMALIGN

589:    Synopsis:
590:     #include <petscsys.h>
591:    PetscErrorCode PetscMalloc1(size_t m1,type **r1)

593:    Not Collective

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

598:    Output Parameter:
599: .  r1 - memory allocated in first chunk

601:    Level: developer

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

605:   Concepts: memory allocation

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

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

613:    Synopsis:
614:     #include <petscsys.h>
615:    PetscErrorCode PetscCalloc1(size_t m1,type **r1)

617:    Not Collective

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

622:    Output Parameter:
623: .  r1 - memory allocated in first chunk

625:    Level: developer

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

629:   Concepts: memory allocation

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

634: /*MC
635:    PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN

637:    Synopsis:
638:     #include <petscsys.h>
639:    PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)

641:    Not Collective

643:    Input Parameter:
644: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
645: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

647:    Output Parameter:
648: +  r1 - memory allocated in first chunk
649: -  r2 - memory allocated in second chunk

651:    Level: developer

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

655:   Concepts: memory allocation

657: M*/
658: #if !defined(PETSC_USE_MALLOC_COALESCED)
659: #define PetscMalloc2(m1,r1,m2,r2) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)))
660: #else
661: #define PetscMalloc2(m1,r1,m2,r2) ((((m1)+(m2)) ? (*(r2) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(PETSC_MEMALIGN-1),r1)) : 0) \
662:                                    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),0) \
663:                                    || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0))
664: #endif

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

669:    Synopsis:
670:     #include <petscsys.h>
671:    PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)

673:    Not Collective

675:    Input Parameter:
676: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
677: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

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

683:    Level: developer

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

687:   Concepts: memory allocation
688: M*/
689: #define PetscCalloc2(m1,r1,m2,r2) (PetscMalloc2((m1),(r1),(m2),(r2)) || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))))

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

694:    Synopsis:
695:     #include <petscsys.h>
696:    PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

698:    Not Collective

700:    Input Parameter:
701: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
702: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
703: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

705:    Output Parameter:
706: +  r1 - memory allocated in first chunk
707: .  r2 - memory allocated in second chunk
708: -  r3 - memory allocated in third chunk

710:    Level: developer

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

714:   Concepts: memory allocation

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

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

728:    Synopsis:
729:     #include <petscsys.h>
730:    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

732:    Not Collective

734:    Input Parameter:
735: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
736: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
737: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

739:    Output Parameter:
740: +  r1 - memory allocated in first chunk
741: .  r2 - memory allocated in second chunk
742: -  r3 - memory allocated in third chunk

744:    Level: developer

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

748:   Concepts: memory allocation
749: M*/
750: #define PetscCalloc3(m1,r1,m2,r2,m3,r3)                                 \
751:   (PetscMalloc3((m1),(r1),(m2),(r2),(m3),(r3))                          \
752:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))))

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

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

761:    Not Collective

763:    Input Parameter:
764: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
765: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
766: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
767: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

769:    Output Parameter:
770: +  r1 - memory allocated in first chunk
771: .  r2 - memory allocated in second chunk
772: .  r3 - memory allocated in third chunk
773: -  r4 - memory allocated in fourth chunk

775:    Level: developer

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

779:   Concepts: memory allocation

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

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

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

798:    Not Collective

800:    Input Parameter:
801: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
802: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
803: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
804: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

806:    Output Parameter:
807: +  r1 - memory allocated in first chunk
808: .  r2 - memory allocated in second chunk
809: .  r3 - memory allocated in third chunk
810: -  r4 - memory allocated in fourth chunk

812:    Level: developer

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

816:   Concepts: memory allocation

818: M*/
819: #define PetscCalloc4(m1,r1,m2,r2,m3,r3,m4,r4)                           \
820:   (PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4)                                \
821:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
822:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))))

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

827:    Synopsis:
828:     #include <petscsys.h>
829:    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)

831:    Not Collective

833:    Input Parameter:
834: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
835: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
836: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
837: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
838: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

840:    Output Parameter:
841: +  r1 - memory allocated in first chunk
842: .  r2 - memory allocated in second chunk
843: .  r3 - memory allocated in third chunk
844: .  r4 - memory allocated in fourth chunk
845: -  r5 - memory allocated in fifth chunk

847:    Level: developer

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

851:   Concepts: memory allocation

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

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

866:    Synopsis:
867:     #include <petscsys.h>
868:    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)

870:    Not Collective

872:    Input Parameter:
873: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
874: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
875: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
876: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
877: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

879:    Output Parameter:
880: +  r1 - memory allocated in first chunk
881: .  r2 - memory allocated in second chunk
882: .  r3 - memory allocated in third chunk
883: .  r4 - memory allocated in fourth chunk
884: -  r5 - memory allocated in fifth chunk

886:    Level: developer

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

890:   Concepts: memory allocation

892: M*/
893: #define PetscCalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5)                     \
894:   (PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5)                          \
895:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
896:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))))

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

901:    Synopsis:
902:     #include <petscsys.h>
903:    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)

905:    Not Collective

907:    Input Parameter:
908: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
909: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
910: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
911: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
912: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
913: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

915:    Output Parameter:
916: +  r1 - memory allocated in first chunk
917: .  r2 - memory allocated in second chunk
918: .  r3 - memory allocated in third chunk
919: .  r4 - memory allocated in fourth chunk
920: .  r5 - memory allocated in fifth chunk
921: -  r6 - memory allocated in sixth chunk

923:    Level: developer

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

927:   Concepts: memory allocation

929: M*/
930: #if !defined(PETSC_USE_MALLOC_COALESCED)
931: #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)))
932: #else
933: #define PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) \
934:   ((((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) \
935:    || (*(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) \
936:    || (!(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))
937: #endif

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

942:    Synopsis:
943:     #include <petscsys.h>
944:    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)

946:    Not Collective

948:    Input Parameter:
949: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
950: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
951: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
952: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
953: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
954: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

956:    Output Parameter:
957: +  r1 - memory allocated in first chunk
958: .  r2 - memory allocated in second chunk
959: .  r3 - memory allocated in third chunk
960: .  r4 - memory allocated in fourth chunk
961: .  r5 - memory allocated in fifth chunk
962: -  r6 - memory allocated in sixth chunk

964:    Level: developer

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

968:   Concepts: memory allocation
969: M*/
970: #define PetscCalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6)               \
971:   (PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6)                    \
972:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
973:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))) || PetscMemzero(*(r6),(m6)*sizeof(**(r6))))

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

978:    Synopsis:
979:     #include <petscsys.h>
980:    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)

982:    Not Collective

984:    Input Parameter:
985: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
986: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
987: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
988: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
989: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
990: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
991: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

993:    Output Parameter:
994: +  r1 - memory allocated in first chunk
995: .  r2 - memory allocated in second chunk
996: .  r3 - memory allocated in third chunk
997: .  r4 - memory allocated in fourth chunk
998: .  r5 - memory allocated in fifth chunk
999: .  r6 - memory allocated in sixth chunk
1000: -  r7 - memory allocated in seventh chunk

1002:    Level: developer

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

1006:   Concepts: memory allocation

1008: M*/
1009: #if !defined(PETSC_USE_MALLOC_COALESCED)
1010: #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)))
1011: #else
1012: #define PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) \
1013:   ((((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) \
1014:    || (*(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) \
1015:    || (!(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))
1016: #endif

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

1021:    Synopsis:
1022:     #include <petscsys.h>
1023:    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)

1025:    Not Collective

1027:    Input Parameter:
1028: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
1029: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
1030: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
1031: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
1032: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
1033: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
1034: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

1036:    Output Parameter:
1037: +  r1 - memory allocated in first chunk
1038: .  r2 - memory allocated in second chunk
1039: .  r3 - memory allocated in third chunk
1040: .  r4 - memory allocated in fourth chunk
1041: .  r5 - memory allocated in fifth chunk
1042: .  r6 - memory allocated in sixth chunk
1043: -  r7 - memory allocated in seventh chunk

1045:    Level: developer

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

1049:   Concepts: memory allocation
1050: M*/
1051: #define PetscCalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7)         \
1052:   (PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7)              \
1053:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
1054:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))) || PetscMemzero(*(r6),(m6)*sizeof(**(r6))) \
1055:    || PetscMemzero(*(r7),(m7)*sizeof(**(r7))))

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

1060:    Synopsis:
1061:     #include <petscsys.h>
1062:    PetscErrorCode PetscNew(type **result)

1064:    Not Collective

1066:    Output Parameter:
1067: .  result - memory allocated, sized to match pointer type

1069:    Level: beginner

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

1073:   Concepts: memory allocation

1075: M*/
1076: #define PetscNew(b)      PetscCalloc1(1,(b))

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

1082:    Synopsis:
1083:     #include <petscsys.h>
1084:    PetscErrorCode PetscNewLog(PetscObject obj,type **result)

1086:    Not Collective

1088:    Input Parameter:
1089: .  obj - object memory is logged to

1091:    Output Parameter:
1092: .  result - memory allocated, sized to match pointer type

1094:    Level: developer

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

1098:   Concepts: memory allocation

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

1103: /*MC
1104:    PetscFree - Frees memory

1106:    Synopsis:
1107:     #include <petscsys.h>
1108:    PetscErrorCode PetscFree(void *memory)

1110:    Not Collective

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

1115:    Level: beginner

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

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

1123:   Concepts: memory allocation

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

1128: /*MC
1129:    PetscFreeVoid - Frees memory

1131:    Synopsis:
1132:     #include <petscsys.h>
1133:    void PetscFreeVoid(void *memory)

1135:    Not Collective

1137:    Input Parameter:
1138: .   memory - memory to free

1140:    Level: beginner

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

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

1146:   Concepts: memory allocation

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


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

1155:    Synopsis:
1156:     #include <petscsys.h>
1157:    PetscErrorCode PetscFree2(void *memory1,void *memory2)

1159:    Not Collective

1161:    Input Parameter:
1162: +   memory1 - memory to free
1163: -   memory2 - 2nd memory to free

1165:    Level: developer

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

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

1171:   Concepts: memory allocation

1173: M*/
1174: #if !defined(PETSC_USE_MALLOC_COALESCED)
1175: #define PetscFree2(m1,m2)   (PetscFree(m2) || PetscFree(m1))
1176: #else
1177: #define PetscFree2(m1,m2)   ((m1) ? ((m2)=0,PetscFree(m1)) : ((m1)=0,PetscFree(m2)))
1178: #endif

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

1183:    Synopsis:
1184:     #include <petscsys.h>
1185:    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)

1187:    Not Collective

1189:    Input Parameter:
1190: +   memory1 - memory to free
1191: .   memory2 - 2nd memory to free
1192: -   memory3 - 3rd memory to free

1194:    Level: developer

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

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

1200:   Concepts: memory allocation

1202: M*/
1203: #if !defined(PETSC_USE_MALLOC_COALESCED)
1204: #define PetscFree3(m1,m2,m3)   (PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1205: #else
1206: #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))))
1207: #endif

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

1212:    Synopsis:
1213:     #include <petscsys.h>
1214:    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)

1216:    Not Collective

1218:    Input Parameter:
1219: +   m1 - memory to free
1220: .   m2 - 2nd memory to free
1221: .   m3 - 3rd memory to free
1222: -   m4 - 4th memory to free

1224:    Level: developer

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

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

1230:   Concepts: memory allocation

1232: M*/
1233: #if !defined(PETSC_USE_MALLOC_COALESCED)
1234: #define PetscFree4(m1,m2,m3,m4)   (PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1235: #else
1236: #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)))))
1237: #endif

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

1242:    Synopsis:
1243:     #include <petscsys.h>
1244:    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)

1246:    Not Collective

1248:    Input Parameter:
1249: +   m1 - memory to free
1250: .   m2 - 2nd memory to free
1251: .   m3 - 3rd memory to free
1252: .   m4 - 4th memory to free
1253: -   m5 - 5th memory to free

1255:    Level: developer

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

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

1261:   Concepts: memory allocation

1263: M*/
1264: #if !defined(PETSC_USE_MALLOC_COALESCED)
1265: #define PetscFree5(m1,m2,m3,m4,m5)   (PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1266: #else
1267: #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)) : \
1268:                                      ((m4) ? ((m5)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m4)) : ((m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m5))))))
1269: #endif


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

1275:    Synopsis:
1276:     #include <petscsys.h>
1277:    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)

1279:    Not Collective

1281:    Input Parameter:
1282: +   m1 - memory to free
1283: .   m2 - 2nd memory to free
1284: .   m3 - 3rd memory to free
1285: .   m4 - 4th memory to free
1286: .   m5 - 5th memory to free
1287: -   m6 - 6th memory to free


1290:    Level: developer

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

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

1296:   Concepts: memory allocation

1298: M*/
1299: #if !defined(PETSC_USE_MALLOC_COALESCED)
1300: #define PetscFree6(m1,m2,m3,m4,m5,m6)   (PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1301: #else
1302: #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)) : \
1303:                                         ((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)) : \
1304:                                         ((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)))))))
1305: #endif

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

1310:    Synopsis:
1311:     #include <petscsys.h>
1312:    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)

1314:    Not Collective

1316:    Input Parameter:
1317: +   m1 - memory to free
1318: .   m2 - 2nd memory to free
1319: .   m3 - 3rd memory to free
1320: .   m4 - 4th memory to free
1321: .   m5 - 5th memory to free
1322: .   m6 - 6th memory to free
1323: -   m7 - 7th memory to free


1326:    Level: developer

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

1330: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1331:           PetscMalloc7()

1333:   Concepts: memory allocation

1335: M*/
1336: #if !defined(PETSC_USE_MALLOC_COALESCED)
1337: #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   (PetscFree(m7) || PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1338: #else
1339: #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)) : \
1340:                                            ((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)) : \
1341:                                            ((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)) : \
1342:                                                    ((m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m7))))))))
1343: #endif

1345: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,int,const char[],const char[],void**);
1346: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1347: PETSC_EXTERN PetscErrorCode PetscMallocSet(PetscErrorCode (*)(size_t,int,const char[],const char[],void**),PetscErrorCode (*)(void*,int,const char[],const char[]));
1348: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);

1350: /*
1351:     PetscLogDouble variables are used to contain double precision numbers
1352:   that are not used in the numerical computations, but rather in logging,
1353:   timing etc.
1354: */
1355: typedef double PetscLogDouble;
1356: #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE

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

1372: /*E
1373:     PetscDataType - Used for handling different basic data types.

1375:    Level: beginner

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

1379: .seealso: PetscBinaryRead(), PetscBinaryWrite(), PetscDataTypeToMPIDataType(),
1380:           PetscDataTypeGetSize()

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

1387: #if defined(PETSC_USE_COMPLEX)
1388: #define  PETSC_SCALAR  PETSC_COMPLEX
1389: #else
1390: #if defined(PETSC_USE_REAL_SINGLE)
1391: #define  PETSC_SCALAR  PETSC_FLOAT
1392: #elif defined(PETSC_USE_REAL___FLOAT128)
1393: #define  PETSC_SCALAR  PETSC___FLOAT128
1394: #else
1395: #define  PETSC_SCALAR  PETSC_DOUBLE
1396: #endif
1397: #endif
1398: #if defined(PETSC_USE_REAL_SINGLE)
1399: #define  PETSC_REAL  PETSC_FLOAT
1400: #elif defined(PETSC_USE_REAL___FLOAT128)
1401: #define  PETSC_REAL  PETSC___FLOAT128
1402: #else
1403: #define  PETSC_REAL  PETSC_DOUBLE
1404: #endif
1405: #define  PETSC_FORTRANADDR  PETSC_LONG

1407: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1408: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1409: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1410: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);

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

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

1447: /*S
1448:     PetscToken - 'Token' used for managing tokenizing strings

1450:   Level: intermediate

1452: .seealso: PetscTokenCreate(), PetscTokenFind(), PetscTokenDestroy()
1453: S*/
1454: typedef struct _p_PetscToken* PetscToken;

1456: PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1457: PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1458: PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);

1460: PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1461: PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);

1463: /*
1464:    These are  MPI operations for MPI_Allreduce() etc
1465: */
1466: PETSC_EXTERN MPI_Op PetscMaxSum_Op;
1467: #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128)
1468: PETSC_EXTERN MPI_Op MPIU_SUM;
1469: #else
1470: #define MPIU_SUM MPI_SUM
1471: #endif
1472: #if defined(PETSC_USE_REAL___FLOAT128)
1473: PETSC_EXTERN MPI_Op MPIU_MAX;
1474: PETSC_EXTERN MPI_Op MPIU_MIN;
1475: #else
1476: #define MPIU_MAX MPI_MAX
1477: #define MPIU_MIN MPI_MIN
1478: #endif
1479: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);

1481: PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1482: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);

1484: /*S
1485:      PetscObject - any PETSc object, PetscViewer, Mat, Vec, KSP etc

1487:    Level: beginner

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

1491: .seealso:  PetscObjectDestroy(), PetscObjectView(), PetscObjectGetName(), PetscObjectSetName(), PetscObjectReference(), PetscObjectDereferenc()
1492: S*/
1493: typedef struct _p_PetscObject* PetscObject;

1495: /*MC
1496:     PetscObjectId - unique integer Id for a PetscObject

1498:     Level: developer

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

1502: .seealso: PetscObjectState, PetscObjectGetId()
1503: M*/
1504: typedef Petsc64bitInt PetscObjectId;

1506: /*MC
1507:     PetscObjectState - integer state for a PetscObject

1509:     Level: developer

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

1515: .seealso: PetscObjectId, PetscObjectStateGet(), PetscObjectStateIncrease(), PetscObjectStateSet()
1516: M*/
1517: typedef Petsc64bitInt PetscObjectState;

1519: /*S
1520:      PetscFunctionList - Linked list of functions, possibly stored in dynamic libraries, accessed
1521:       by string name

1523:    Level: advanced

1525: .seealso:  PetscFunctionListAdd(), PetscFunctionListDestroy(), PetscOpFlist
1526: S*/
1527: typedef struct _n_PetscFunctionList *PetscFunctionList;

1529: /*E
1530:   PetscFileMode - Access mode for a file.

1532:   Level: beginner

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

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

1538:   FILE_MODE_APPEND - open a file at end for writing

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

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

1544: .seealso: PetscViewerFileSetMode()
1545: E*/
1546: typedef enum {FILE_MODE_READ, FILE_MODE_WRITE, FILE_MODE_APPEND, FILE_MODE_UPDATE, FILE_MODE_APPEND_UPDATE} PetscFileMode;
1547: extern const char *const PetscFileModes[];

1549: /*
1550:     Defines PETSc error handling.
1551: */
1552: #include <petscerror.h>

1554: #define PETSC_SMALLEST_CLASSID  1211211
1555: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1556: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1557: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);

1559: /*
1560:    Routines that get memory usage information from the OS
1561: */
1562: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1563: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1564: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1565: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);

1567: PETSC_EXTERN PetscErrorCode PetscInfoAllow(PetscBool ,const char []);
1568: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);

1570: /*
1571:    Initialization of PETSc
1572: */
1573: PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1574: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1575: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1576: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1577: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1578: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1579: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1580: PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1581: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1582: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);

1584: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1585: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);

1587: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1588: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1589: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1590: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);

1592: /*
1593:      These are so that in extern C code we can caste function pointers to non-extern C
1594:    function pointers. Since the regular C++ code expects its function pointers to be C++
1595: */
1596: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1597: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1598: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);

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

1632: #include <petscviewertypes.h>
1633: #include <petscoptions.h>

1635: PETSC_EXTERN PetscErrorCode PetscMemoryShowUsage(PetscViewer,const char[]);
1636: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1637: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1638: #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1639: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1640: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1641: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1642: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1643: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1644: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1645: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1646: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1647: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,const char[],const char[]);
1648: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1649: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1650: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1651: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1652: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);

1654: #if defined(PETSC_HAVE_SAWS)
1655: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1656: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1657: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1658: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1659: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1660: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1661: PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1662: PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1663: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1664: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);

1666: #else
1667: #define PetscSAWsBlock()                        0
1668: #define PetscObjectSAWsViewOff(obj)             0
1669: #define PetscObjectSAWsSetBlock(obj,flg)        0
1670: #define PetscObjectSAWsBlock(obj)               0
1671: #define PetscObjectSAWsGrantAccess(obj)         0
1672: #define PetscObjectSAWsTakeAccess(obj)          0
1673: #define PetscStackViewSAWs()                    0
1674: #define PetscStackSAWsViewOff()                 0
1675: #define PetscStackSAWsTakeAccess()
1676: #define PetscStackSAWsGrantAccess()

1678: #endif

1680: typedef void* PetscDLHandle;
1681: typedef enum {PETSC_DL_DECIDE=0,PETSC_DL_NOW=1,PETSC_DL_LOCAL=2} PetscDLMode;
1682: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1683: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1684: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);


1687: #if defined(PETSC_USE_DEBUG)
1688: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1689: #endif
1690: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);

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

1695:    Level: developer

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

1699: .seealso:  PetscObjectListAdd(), PetscObjectListDestroy(), PetscObjectListFind(), PetscObjectCompose(), PetscObjectQuery(), PetscFunctionList
1700: S*/
1701: typedef struct _n_PetscObjectList *PetscObjectList;

1703: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1704: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1705: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1706: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1707: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1708: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);

1710: /*
1711:     Dynamic library lists. Lists of names of routines in objects or in dynamic
1712:   link libraries that will be loaded as needed.
1713: */

1715: #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1716: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1717: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1718: #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1719: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1720: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[]);
1721: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1722: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1723: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);

1725: /*S
1726:      PetscDLLibrary - Linked list of dynamics libraries to search for functions

1728:    Level: advanced

1730: .seealso:  PetscDLLibraryOpen()
1731: S*/
1732: typedef struct _n_PetscDLLibrary *PetscDLLibrary;
1733: PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1734: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1735: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1736: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1737: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1738: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1739: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1740: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);

1742: /*
1743:      Useful utility routines
1744: */
1745: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1746: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1747: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1748: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1749: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1750: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);

1752: /*
1753:     PetscNot - negates a logical type value and returns result as a PetscBool

1755:     Notes: This is useful in cases like
1756: $     int        *a;
1757: $     PetscBool  flag = PetscNot(a)
1758:      where !a does not return a PetscBool  because we cannot provide a cast from int to PetscBool  in C.
1759: */
1760: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)

1762: #if defined(PETSC_HAVE_VALGRIND)
1763: #  include <valgrind/valgrind.h>
1764: #  define PETSC_RUNNING_ON_VALGRIND RUNNING_ON_VALGRIND
1765: #else
1766: #  define PETSC_RUNNING_ON_VALGRIND PETSC_FALSE
1767: #endif


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

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

1777:     Not Collective

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

1782:    Level: developer

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

1787:     Concepts: help messages^printing
1788:     Concepts: printing^help messages

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

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



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

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

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

1822: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1823: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1824: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);

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

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

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

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

1845:    Level: advanced

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

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

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

1872: #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1873: #include <xmmintrin.h>
1874: #endif

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

1883:    Not Collective

1885:    Input Parameters:
1886: +  b - pointer to initial memory space
1887: -  n - length (in bytes) of space to copy

1889:    Output Parameter:
1890: .  a - pointer to copy space

1892:    Level: intermediate

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

1902:    Note:
1903:    This routine is analogous to memcpy().

1905:    Developer Note: this is inlined for fastest performance

1907:   Concepts: memory^copying
1908:   Concepts: copying^memory

1910: .seealso: PetscMemmove()

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

1955: /*@C
1956:    PetscMemzero - Zeros the specified memory.

1958:    Not Collective

1960:    Input Parameters:
1961: +  a - pointer to beginning memory location
1962: -  n - length (in bytes) of memory to initialize

1964:    Level: intermediate

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

1970:    Developer Note: this is inlined for fastest performance

1972:    Concepts: memory^zeroing
1973:    Concepts: zeroing^memory

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

2007: /*MC
2008:    PetscPrefetchBlock - Prefetches a block of memory

2010:    Synopsis:
2011:     #include <petscsys.h>
2012:     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)

2014:    Not Collective

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

2022:    Level: developer

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

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

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

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

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

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

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

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

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

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

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

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

2081: #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
2082: #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
2083: #endif

2085: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
2086: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
2087: #endif

2089: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
2090: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
2091: #endif

2093: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
2094: #define PETSC_USE_FORTRAN_KERNEL_MDOT
2095: #endif

2097: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
2098: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
2099: #endif

2101: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
2102: #define PETSC_USE_FORTRAN_KERNEL_AYPX
2103: #endif

2105: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
2106: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
2107: #endif

2109: #endif

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

2126: /* --------------------------------------------------------------------*/

2128: /*MC
2129:     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
2130:         communication

2132:    Level: beginner

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

2136: .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
2137: M*/

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


2145:    Level: beginner

2147: .seealso: PetscReal, PassiveReal, PassiveScalar, MPIU_SCALAR, PetscInt
2148: M*/

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

2153:    Synopsis:
2154:    #define PETSC_DESIRE_COMPLEX
2155:    #include <petscsys.h>
2156:    PetscComplex number = 1. + 2.*PETSC_i;

2158:    Level: beginner

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

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

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

2172:    Level: beginner

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

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

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

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

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

2190:    Level: beginner

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

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

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

2201:    Level: beginner

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

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

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

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

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

2230: PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
2231: {
2233:   *b =  (PetscBLASInt)(a);
2234: #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
2235:   if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK");
2236: #endif
2237:   return(0);
2238: }

2242: PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
2243: {
2245:   *b =  (PetscMPIInt)(a);
2246: #if defined(PETSC_USE_64BIT_INDICES)
2247:   if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for MPI");
2248: #endif
2249:   return(0);
2250: }


2253: /*
2254:      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2255: */
2256: #if defined(hz)
2257: #undef hz
2258: #endif

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


2263: #if defined(PETSC_HAVE_LIMITS_H)
2264: #include <limits.h>
2265: #endif
2266: #if defined(PETSC_HAVE_SYS_PARAM_H)
2267: #include <sys/param.h>
2268: #endif
2269: #if defined(PETSC_HAVE_SYS_TYPES_H)
2270: #include <sys/types.h>
2271: #endif
2272: #if defined(MAXPATHLEN)
2273: #  define PETSC_MAX_PATH_LEN     MAXPATHLEN
2274: #elif defined(MAX_PATH)
2275: #  define PETSC_MAX_PATH_LEN     MAX_PATH
2276: #elif defined(_MAX_PATH)
2277: #  define PETSC_MAX_PATH_LEN     _MAX_PATH
2278: #else
2279: #  define PETSC_MAX_PATH_LEN     4096
2280: #endif

2282: /*MC

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

2286: $    1) classic Fortran 77 style
2287: $#include "finclude/petscXXX.h" to work with material from the XXX component of PETSc
2288: $       XXX variablename
2289: $      You cannot use this approach if you wish to use the Fortran 90 specific PETSc routines
2290: $      which end in F90; such as VecGetArrayF90()
2291: $
2292: $    2) classic Fortran 90 style
2293: $#include "finclude/petscXXX.h"
2294: $#include "finclude/petscXXX.h90" to work with material from the XXX component of PETSc
2295: $       XXX variablename
2296: $
2297: $    3) Using Fortran modules
2298: $#include "finclude/petscXXXdef.h"
2299: $         use petscXXXX
2300: $       XXX variablename
2301: $
2302: $    4) Use Fortran modules and Fortran data types for PETSc types
2303: $#include "finclude/petscXXXdef.h"
2304: $         use petscXXXX
2305: $       type(XXX) variablename
2306: $      To use this approach you must ./configure PETSc with the additional
2307: $      option --with-fortran-datatypes You cannot use the type(XXX) declaration approach without using Fortran modules

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

2311: $    3a) skip the #include BUT you cannot use any PETSc data type names like Vec, Mat, PetscInt, PetscErrorCode etc
2312: $        and you must declare the variables as integer, for example
2313: $        integer variablename
2314: $
2315: $    4a) skip the #include, you use the object types like type(Vec) type(Mat) but cannot use the data type
2316: $        names like PetscErrorCode, PetscInt etc. again for those you must use integer

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

2321:    Fortran type checking with interfaces is strick, this means you cannot pass a scalar value when an array value
2322: is expected (even though it is legal Fortran). For example when setting a single value in a matrix with MatSetValues()
2323: you cannot have something like
2324: $      PetscInt row,col
2325: $      PetscScalar val
2326: $        ...
2327: $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)
2328: You must instead have
2329: $      PetscInt row(1),col(1)
2330: $      PetscScalar val(1)
2331: $        ...
2332: $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)


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

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

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

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

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

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

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

2357:     Level: beginner

2359: M*/

2361: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2362: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2363: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2364: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2365: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2366: PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2367: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);

2369: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2370: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2371: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2372: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2373: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2374: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2375: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt*,PetscInt*,PetscInt*);
2376: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2377: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2378: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2379: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2380: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2381: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2382: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2383: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2384: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2385: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2386: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt*,const PetscInt*,PetscInt,const PetscInt*,const PetscInt*,PetscInt*,PetscInt**,PetscInt**);
2387: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt,const PetscInt*,PetscInt,const PetscInt*,PetscInt*,PetscInt**);

2389: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2390: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);

2392: /*J
2393:     PetscRandomType - String with the name of a PETSc randomizer

2395:    Level: beginner

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

2400: .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2401: J*/
2402: typedef const char* PetscRandomType;
2403: #define PETSCRAND       "rand"
2404: #define PETSCRAND48     "rand48"
2405: #define PETSCSPRNG      "sprng"

2407: /* Logging support */
2408: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;

2410: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);

2412: /*S
2413:      PetscRandom - Abstract PETSc object that manages generating random numbers

2415:    Level: intermediate

2417:   Concepts: random numbers

2419: .seealso:  PetscRandomCreate(), PetscRandomGetValue(), PetscRandomType
2420: S*/
2421: typedef struct _p_PetscRandom*   PetscRandom;

2423: /* Dynamic creation and loading functions */
2424: PETSC_EXTERN PetscFunctionList PetscRandomList;
2425: PETSC_EXTERN PetscBool         PetscRandomRegisterAllCalled;

2427: PETSC_EXTERN PetscErrorCode PetscRandomRegisterAll(void);
2428: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2429: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2430: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2431: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType*);
2432: PETSC_STATIC_INLINE PetscErrorCode PetscRandomViewFromOptions(PetscRandom A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
2433: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);

2435: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2436: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2437: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2438: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2439: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2440: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2441: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2442: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2443: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);

2445: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2446: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2447: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2448: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2449: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2450: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2451: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);

2453: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscDataType);
2454: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscDataType);
2455: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,void*,PetscInt,PetscDataType,PetscBool );
2456: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,void*,PetscInt,PetscDataType,PetscBool );
2457: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2458: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2459: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2460: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2461: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2462: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2463: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2464: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);
2465: PETSC_EXTERN PetscErrorCode PetscWebServe(MPI_Comm,int);

2467: /*
2468:    In binary files variables are stored using the following lengths,
2469:   regardless of how they are stored in memory on any one particular
2470:   machine. Use these rather then sizeof() in computing sizes for
2471:   PetscBinarySeek().
2472: */
2473: #define PETSC_BINARY_INT_SIZE   (32/8)
2474: #define PETSC_BINARY_FLOAT_SIZE  (32/8)
2475: #define PETSC_BINARY_CHAR_SIZE  (8/8)
2476: #define PETSC_BINARY_SHORT_SIZE  (16/8)
2477: #define PETSC_BINARY_DOUBLE_SIZE  (64/8)
2478: #define PETSC_BINARY_SCALAR_SIZE  sizeof(PetscScalar)

2480: /*E
2481:   PetscBinarySeekType - argument to PetscBinarySeek()

2483:   Level: advanced

2485: .seealso: PetscBinarySeek(), PetscBinarySynchronizedSeek()
2486: E*/
2487: typedef enum {PETSC_BINARY_SEEK_SET = 0,PETSC_BINARY_SEEK_CUR = 1,PETSC_BINARY_SEEK_END = 2} PetscBinarySeekType;
2488: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2489: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2490: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);

2492: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2493: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool );
2494: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2495: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2496: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2497: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);

2499: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2500: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2501: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2502: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2503: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2504: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscInt,const PetscMPIInt*,const void*,PetscInt*,PetscMPIInt**,void*) PetscAttrMPIPointerWithType(6,3);

2506: /*E
2507:     PetscBuildTwoSidedType - algorithm for setting up two-sided communication

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

2515:    Level: developer

2517: .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedSetType(), PetscCommBuildTwoSidedGetType()
2518: E*/
2519: typedef enum {
2520:   PETSC_BUILDTWOSIDED_NOTSET = -1,
2521:   PETSC_BUILDTWOSIDED_ALLREDUCE = 0,
2522:   PETSC_BUILDTWOSIDED_IBARRIER = 1
2523:   /* Updates here must be accompanied by updates in finclude/petscsys.h and the string array in mpits.c */
2524: } PetscBuildTwoSidedType;
2525: PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2526: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2527: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);

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

2531: /*E
2532:   InsertMode - Whether entries are inserted or added into vectors or matrices

2534:   Level: beginner

2536: .seealso: VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2537:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(),
2538:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd()
2539: E*/
2540:  typedef enum {NOT_SET_VALUES, INSERT_VALUES, ADD_VALUES, MAX_VALUES, INSERT_ALL_VALUES, ADD_ALL_VALUES, INSERT_BC_VALUES, ADD_BC_VALUES} InsertMode;

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

2545:     Level: beginner

2547: .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2548:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), ADD_VALUES,
2549:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES

2551: M*/

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

2557:     Level: beginner

2559: .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2560:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), INSERT_VALUES,
2561:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES

2563: M*/

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

2568:     Level: beginner

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

2572: M*/

2574: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);

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

2579: /*S
2580:    PetscSubcomm - Context of MPI subcommunicators, used by PCREDUNDANT

2582:    Level: advanced

2584:    Concepts: communicator, create
2585: S*/
2586: typedef struct _n_PetscSubcomm* PetscSubcomm;

2588: struct _n_PetscSubcomm {
2589:   MPI_Comm    parent;           /* parent communicator */
2590:   MPI_Comm    dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2591:   MPI_Comm    comm;             /* this communicator */
2592:   PetscMPIInt n;                /* num of subcommunicators under the parent communicator */
2593:   PetscMPIInt color;            /* color of processors belong to this communicator */
2594:   PetscMPIInt *subsize;         /* size of subcommunicator[color] */
2595:   PetscSubcommType type;
2596: };

2598: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2599: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2600: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2601: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2602: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2603: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2604: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);

2606: /*S
2607:    PetscSegBuffer - a segmented extendable buffer

2609:    Level: developer

2611: .seealso: PetscSegBufferCreate(), PetscSegBufferGet(), PetscSegBufferExtract(), PetscSegBufferDestroy()
2612: S*/
2613: typedef struct _n_PetscSegBuffer *PetscSegBuffer;
2614: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2615: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2616: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2617: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2618: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2619: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2620: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2621: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);

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

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

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

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

2640:    Level: intermediate

2642:      Options Database:
2643: .     -citations [filenmae]   - print out the bibtex entries for the given computation
2644: @*/
2645: PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2646: {
2647:   size_t         len;
2648:   char           *vstring;

2652:   if (set && *set) return(0);
2653:   PetscStrlen(cit,&len);
2654:   PetscSegBufferGet(PetscCitationsList,(PetscInt)len,&vstring);
2655:   PetscMemcpy(vstring,cit,len);
2656:   if (set) *set = PETSC_TRUE;
2657:   return(0);
2658: }

2660: PETSC_EXTERN PetscErrorCode PetscURLShorten(const char[],char[],size_t);
2661: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm,char[],char[],size_t);
2662: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm,const char[],char[],size_t);
2663: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm,const char[],const char []);

2665: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm,char[],char[],size_t);
2666: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm,const char[],char[],char[],size_t);

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

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

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

2677: #endif