Actual source code: petscsys.h

petsc-3.7.5 2017-01-01
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:    For --prefix installs the ${PETSC_ARCH}/ does not exist and petscconf.h is in the same
 12:    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: /* ========================================================================== */
 47: /*
 48:    Since PETSc manages its own extern "C" handling users should never include PETSc include
 49:    files within extern "C". This will generate a compiler error if a user does put the include 
 50:    file within an extern "C".
 51: */
 52: #if defined(__cplusplus)
 53: void assert_never_put_petsc_headers_inside_an_extern_c(int); void assert_never_put_petsc_headers_inside_an_extern_c(double);
 54: #endif

 56: #if defined(__cplusplus)
 57: #  define PETSC_RESTRICT PETSC_CXX_RESTRICT
 58: #else
 59: #  define PETSC_RESTRICT PETSC_C_RESTRICT
 60: #endif

 62: #if defined(__cplusplus)
 63: #  define PETSC_STATIC_INLINE PETSC_CXX_STATIC_INLINE
 64: #else
 65: #  define PETSC_STATIC_INLINE PETSC_C_STATIC_INLINE
 66: #endif

 68: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES) /* For Win32 shared libraries */
 69: #  define  __declspec(dllexport)
 70: #  define PETSC_DLLIMPORT __declspec(dllimport)
 71: #  define PETSC_VISIBILITY_INTERNAL
 72: #elif defined(PETSC_USE_VISIBILITY_CXX) && defined(__cplusplus)
 73: #  define  __attribute__((visibility ("default")))
 74: #  define PETSC_DLLIMPORT __attribute__((visibility ("default")))
 75: #  define PETSC_VISIBILITY_INTERNAL __attribute__((visibility ("hidden")))
 76: #elif defined(PETSC_USE_VISIBILITY_C) && !defined(__cplusplus)
 77: #  define  __attribute__((visibility ("default")))
 78: #  define PETSC_DLLIMPORT __attribute__((visibility ("default")))
 79: #  define PETSC_VISIBILITY_INTERNAL __attribute__((visibility ("hidden")))
 80: #else
 81: #  define 
 82: #  define PETSC_DLLIMPORT
 83: #  define PETSC_VISIBILITY_INTERNAL
 84: #endif

 86: #if defined(petsc_EXPORTS)      /* CMake defines this when building the shared library */
 87: #  define PETSC_VISIBILITY_PUBLIC 
 88: #else  /* Win32 users need this to import symbols from petsc.dll */
 89: #  define PETSC_VISIBILITY_PUBLIC PETSC_DLLIMPORT
 90: #endif

 92: /*
 93:     Functions tagged with PETSC_EXTERN in the header files are
 94:   always defined as extern "C" when compiled with C++ so they may be
 95:   used from C and are always visible in the shared libraries
 96: */
 97: #if defined(__cplusplus)
 98: #define PETSC_EXTERN extern "C" PETSC_VISIBILITY_PUBLIC
 99: #define PETSC_EXTERN_TYPEDEF extern "C"
100: #define PETSC_INTERN extern "C" PETSC_VISIBILITY_INTERNAL
101: #else
102: #define PETSC_EXTERN extern PETSC_VISIBILITY_PUBLIC
103: #define PETSC_EXTERN_TYPEDEF
104: #define PETSC_INTERN extern PETSC_VISIBILITY_INTERNAL
105: #endif

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

110: /* ========================================================================== */

112: /*
113:     Defines the interface to MPI allowing the use of all MPI functions.

115:     PETSc does not use the C++ binding of MPI at ALL. The following flag
116:     makes sure the C++ bindings are not included. The C++ bindings REQUIRE
117:     putting mpi.h before ANY C++ include files, we cannot control this
118:     with all PETSc users. Users who want to use the MPI C++ bindings can include
119:     mpicxx.h directly in their code
120: */
121: #if !defined(MPICH_SKIP_MPICXX)
122: #  define MPICH_SKIP_MPICXX 1
123: #endif
124: #if !defined(OMPI_SKIP_MPICXX)
125: #  define OMPI_SKIP_MPICXX 1
126: #endif
127: #if !defined(OMPI_WANT_MPI_INTERFACE_WARNING)
128: #  define OMPI_WANT_MPI_INTERFACE_WARNING 0
129: #endif
130: #include <mpi.h>

132: /*
133:    Perform various sanity checks that the correct mpi.h is being included at compile time.
134:    This usually happens because
135:       * either an unexpected mpi.h is in the default compiler path (i.e. in /usr/include) or
136:       * an extra include path -I/something (which contains the unexpected mpi.h) is being passed to the compiler
137: */
138: #if defined(PETSC_HAVE_MPIUNI)
139: #  if !defined(__MPIUNI_H)
140: #    error "PETSc was configured with --with-mpi=0 but now appears to be compiling using a different mpi.h"
141: #  endif
142: #elif defined(PETSC_HAVE_MPICH_NUMVERSION)
143: #  if !defined(MPICH_NUMVERSION)
144: #    error "PETSc was configured with MPICH but now appears to be compiling using a non-MPICH mpi.h"
145: #  elif MPICH_NUMVERSION != PETSC_HAVE_MPICH_NUMVERSION
146: #    error "PETSc was configured with one MPICH mpi.h version but now appears to be compiling using a different MPICH mpi.h version"
147: #  endif
148: #elif defined(PETSC_HAVE_OMPI_MAJOR_VERSION)
149: #  if !defined(OMPI_MAJOR_VERSION)
150: #    error "PETSc was configured with OpenMPI but now appears to be compiling using a non-OpenMPI mpi.h"
151: #  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)
152: #    error "PETSc was configured with one OpenMPI mpi.h version but now appears to be compiling using a different OpenMPI mpi.h version"
153: #  endif
154: #endif

156: /*
157:     Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler
158:     see the top of mpicxx.h in the MPICH2 distribution.
159: */
160: #include <stdio.h>

162: /* MSMPI on 32bit windows requires this yukky hack - that breaks MPI standard compliance */
163: #if !defined(MPIAPI)
164: #define MPIAPI
165: #endif

167: /*
168:    Support for Clang (>=3.2) matching type tag arguments with void* buffer types.
169:    This allows the compiler to detect cases where the MPI datatype argument passed to a MPI routine
170:    does not match the actual type of the argument being passed in
171: */
172: #if defined(__has_attribute) && defined(works_with_const_which_is_not_true)
173: #  if __has_attribute(argument_with_type_tag) && __has_attribute(pointer_with_type_tag) && __has_attribute(type_tag_for_datatype)
174: #    define PetscAttrMPIPointerWithType(bufno,typeno) __attribute__((pointer_with_type_tag(MPI,bufno,typeno)))
175: #    define PetscAttrMPITypeTag(type)                 __attribute__((type_tag_for_datatype(MPI,type)))
176: #    define PetscAttrMPITypeTagLayoutCompatible(type) __attribute__((type_tag_for_datatype(MPI,type,layout_compatible)))
177: #  endif
178: #endif
179: #if !defined(PetscAttrMPIPointerWithType)
180: #  define PetscAttrMPIPointerWithType(bufno,typeno)
181: #  define PetscAttrMPITypeTag(type)
182: #  define PetscAttrMPITypeTagLayoutCompatible(type)
183: #endif

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

188:     Level: beginner

190: .seealso: CHKERRQ, SETERRQ
191: M*/
192: typedef int PetscErrorCode;

194: /*MC

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

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

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

204:     Level: developer

206: .seealso: PetscClassIdRegister(), PetscLogEventRegister(), PetscHeaderCreate()
207: M*/
208: typedef int PetscClassId;


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

214:     Level: intermediate

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

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

222: .seealso: PetscBLASInt, PetscInt, PetscMPIIntCast()

224: M*/
225: typedef int PetscMPIInt;

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

230:     Level: intermediate

232: .seealso: PetscOptionsGetEnum(), PetscOptionsEnum(), PetscBagRegisterEnum()
233: M*/
234: typedef enum { ENUM_DUMMY } PetscEnum;
235: PETSC_EXTERN MPI_Datatype MPIU_ENUM PetscAttrMPITypeTag(PetscEnum);

237: #if defined(PETSC_HAVE_STDINT_H)
238: #include <stdint.h>
239: #endif
240: #if defined (PETSC_HAVE_INTTYPES_H)
242: #include <inttypes.h>
243: # if !defined(PRId64)
244: # define PRId64 "ld"
245: # endif
246: #endif

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

253:    Level: intermediate

255: .seealso: PetscScalar, PetscBLASInt, PetscMPIInt
256: M*/
257: #if defined(PETSC_HAVE_STDINT_H) && defined(PETSC_HAVE_INTTYPES_H) && defined(PETSC_HAVE_MPI_INT64_T) /* MPI_INT64_T is not guaranteed to be a macro */
258: typedef int64_t Petsc64bitInt;
259: # define MPIU_INT64 MPI_INT64_T
260: # define PetscInt64_FMT PRId64
261: #elif (PETSC_SIZEOF_LONG_LONG == 8)
262: typedef long long Petsc64bitInt;
263: # define MPIU_INT64 MPI_LONG_LONG_INT
264: # define PetscInt64_FMT "lld"
265: #elif defined(PETSC_HAVE___INT64)
266: typedef __int64 Petsc64bitInt;
267: # define MPIU_INT64 MPI_INT64_T
268: # define PetscInt64_FMT "ld"
269: #else
270: #error "cannot determine Petsc64bitInt type"
271: #endif
272: #if defined(PETSC_USE_64BIT_INDICES)
273: typedef Petsc64bitInt PetscInt;
274: #define MPIU_INT MPIU_INT64
275: #define PetscInt_FMT PetscInt64_FMT
276: #else
277: typedef int PetscInt;
278: #define MPIU_INT MPI_INT
279: #define PetscInt_FMT "d"
280: #endif

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

285:     Level: intermediate

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

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

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

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

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

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

306: .seealso: PetscMPIInt, PetscInt, PetscBLASIntCast()

308: M*/
309: #if defined(PETSC_HAVE_64BIT_BLAS_INDICES)
310: typedef Petsc64bitInt PetscBLASInt;
311: #else
312: typedef int PetscBLASInt;
313: #endif

315: /*EC

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

319:     Level: advanced

321: .seealso: PetscObjectSetPrecision()
322: E*/
323: typedef enum { PETSC_PRECISION_SINGLE=4,PETSC_PRECISION_DOUBLE=8 } PetscPrecision;
324: PETSC_EXTERN const char *PetscPrecisions[];

326: /*
327:     For the rare cases when one needs to send a size_t object with MPI
328: */
329: #if (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_INT)
330: #define MPIU_SIZE_T MPI_UNSIGNED
331: #elif  (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_LONG)
332: #define MPIU_SIZE_T MPI_UNSIGNED_LONG
333: #elif  (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_LONG_LONG)
334: #define MPIU_SIZE_T MPI_UNSIGNED_LONG_LONG
335: #else
336: #error "Unknown size for size_t! Send us a bugreport at petsc-maint@mcs.anl.gov"
337: #endif

339: /*
340:       You can use PETSC_STDOUT as a replacement of stdout. You can also change
341:     the value of PETSC_STDOUT to redirect all standard output elsewhere
342: */
343: PETSC_EXTERN FILE* PETSC_STDOUT;

345: /*
346:       You can use PETSC_STDERR as a replacement of stderr. You can also change
347:     the value of PETSC_STDERR to redirect all standard error elsewhere
348: */
349: PETSC_EXTERN FILE* PETSC_STDERR;

351: /*MC
352:     PetscUnlikely - hints the compiler that the given condition is usually FALSE

354:     Synopsis:
355:     #include <petscsys.h>
356:     PetscBool  PetscUnlikely(PetscBool  cond)

358:     Not Collective

360:     Input Parameters:
361: .   cond - condition or expression

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

366:     Level: advanced

368: .seealso: PetscLikely(), CHKERRQ
369: M*/

371: /*MC
372:     PetscLikely - hints the compiler that the given condition is usually TRUE

374:     Synopsis:
375:     #include <petscsys.h>
376:     PetscBool  PetscUnlikely(PetscBool  cond)

378:     Not Collective

380:     Input Parameters:
381: .   cond - condition or expression

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

386:     Level: advanced

388: .seealso: PetscUnlikely()
389: M*/
390: #if defined(PETSC_HAVE_BUILTIN_EXPECT)
391: #  define PetscUnlikely(cond)   __builtin_expect(!!(cond),0)
392: #  define PetscLikely(cond)     __builtin_expect(!!(cond),1)
393: #else
394: #  define PetscUnlikely(cond)   (cond)
395: #  define PetscLikely(cond)     (cond)
396: #endif

398: /*
399:     Declare extern C stuff after including external header files
400: */


403: /*
404:        Basic PETSc constants
405: */

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

410:    Level: beginner

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

415: .seealso: PETSC_TRUE, PETSC_FALSE, PetscNot()
416: E*/
417: typedef enum { PETSC_FALSE,PETSC_TRUE } PetscBool;
418: PETSC_EXTERN const char *const PetscBools[];
419: PETSC_EXTERN MPI_Datatype MPIU_BOOL PetscAttrMPITypeTag(PetscBool);

421: /*
422:     Defines some elementary mathematics functions and constants.
423: */
424: #include <petscmath.h>

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

429:    Level: beginner

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

437: E*/
438: typedef enum { PETSC_COPY_VALUES, PETSC_OWN_POINTER, PETSC_USE_POINTER} PetscCopyMode;
439: PETSC_EXTERN const char *const PetscCopyModes[];

441: /*MC
442:     PETSC_FALSE - False value of PetscBool

444:     Level: beginner

446:     Note: Zero integer

448: .seealso: PetscBool, PETSC_TRUE
449: M*/

451: /*MC
452:     PETSC_TRUE - True value of PetscBool

454:     Level: beginner

456:     Note: Nonzero integer

458: .seealso: PetscBool, PETSC_FALSE
459: M*/

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

464:    Level: beginner

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

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

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

474: M*/
475: #define PETSC_NULL           NULL

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

480:    Level: beginner

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

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

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

490: M*/
491: #define PETSC_IGNORE         NULL

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

497:    Level: beginner

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

501: M*/
502: #define PETSC_DECIDE  -1

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

508:    Level: beginner

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

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

515: M*/
516: #define PETSC_DETERMINE PETSC_DECIDE

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

522:    Level: beginner

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

526: .seealso: PETSC_DECIDE, PETSC_IGNORE, PETSC_DETERMINE

528: M*/
529: #define PETSC_DEFAULT  -2

531: /*MC
532:     PETSC_COMM_WORLD - the equivalent of the MPI_COMM_WORLD communicator which represents
533:            all the processs that PETSc knows about.

535:    Level: beginner

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

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

545: .seealso: PETSC_COMM_SELF

547: M*/
548: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;

550: /*MC
551:     PETSC_COMM_SELF - This is always MPI_COMM_SELF

553:    Level: beginner

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

557: .seealso: PETSC_COMM_WORLD

559: M*/
560: #define PETSC_COMM_SELF MPI_COMM_SELF

562: PETSC_EXTERN PetscBool PetscBeganMPI;
563: PETSC_EXTERN PetscBool PetscInitializeCalled;
564: PETSC_EXTERN PetscBool PetscFinalizeCalled;
565: PETSC_EXTERN PetscBool PetscCUSPSynchronize;
566: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;
567: PETSC_EXTERN PetscBool PetscCUDASynchronize;

569: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm));
570: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*);
571: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm*);

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

576:    Synopsis:
577:     #include <petscsys.h>
578:    PetscErrorCode PetscMalloc(size_t m,void **result)

580:    Not Collective

582:    Input Parameter:
583: .  m - number of bytes to allocate

585:    Output Parameter:
586: .  result - memory allocated

588:    Level: beginner

590:    Notes:
591:    Memory is always allocated at least double aligned

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

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

597:   Concepts: memory allocation

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

602: /*MC
603:    PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment

605:    Synopsis:
606:     #include <petscsys.h>
607:    void *PetscAddrAlign(void *addr)

609:    Not Collective

611:    Input Parameters:
612: .  addr - address to align (any pointer type)

614:    Level: developer

616: .seealso: PetscMallocAlign()

618:   Concepts: memory allocation
619: M*/
620: #define PetscAddrAlign(a) (void*)((((PETSC_UINTPTR_T)(a))+(PETSC_MEMALIGN-1)) & ~(PETSC_MEMALIGN-1))

622: /*MC
623:    PetscMalloc1 - Allocates an array of memory aligned to PETSC_MEMALIGN

625:    Synopsis:
626:     #include <petscsys.h>
627:    PetscErrorCode PetscMalloc1(size_t m1,type **r1)

629:    Not Collective

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

634:    Output Parameter:
635: .  r1 - memory allocated in first chunk

637:    Level: developer

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

641:   Concepts: memory allocation

643: M*/
644: #define PetscMalloc1(m1,r1) PetscMalloc((m1)*sizeof(**(r1)),r1)

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

649:    Synopsis:
650:     #include <petscsys.h>
651:    PetscErrorCode PetscCalloc1(size_t m1,type **r1)

653:    Not Collective

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

658:    Output Parameter:
659: .  r1 - memory allocated in first chunk

661:    Level: developer

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

665:   Concepts: memory allocation

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

670: /*MC
671:    PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN

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

677:    Not Collective

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

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

687:    Level: developer

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

691:   Concepts: memory allocation

693: M*/
694: #if !defined(PETSC_USE_MALLOC_COALESCED)
695: #define PetscMalloc2(m1,r1,m2,r2) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)))
696: #else
697: #define PetscMalloc2(m1,r1,m2,r2) ((((m1)+(m2)) ? (*(r2) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(PETSC_MEMALIGN-1),r1)) : 0) \
698:                                    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),0) \
699:                                    || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0))
700: #endif

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

705:    Synopsis:
706:     #include <petscsys.h>
707:    PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)

709:    Not Collective

711:    Input Parameter:
712: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
713: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

715:    Output Parameter:
716: +  r1 - memory allocated in first chunk
717: -  r2 - memory allocated in second chunk

719:    Level: developer

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

723:   Concepts: memory allocation
724: M*/
725: #define PetscCalloc2(m1,r1,m2,r2) (PetscMalloc2((m1),(r1),(m2),(r2)) || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))))

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

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

734:    Not Collective

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

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

746:    Level: developer

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

750:   Concepts: memory allocation

752: M*/
753: #if !defined(PETSC_USE_MALLOC_COALESCED)
754: #define PetscMalloc3(m1,r1,m2,r2,m3,r3) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)) || PetscMalloc1((m3),(r3)))
755: #else
756: #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) \
757:                                          || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),0) \
758:                                          || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0))
759: #endif

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

764:    Synopsis:
765:     #include <petscsys.h>
766:    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

768:    Not Collective

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

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

780:    Level: developer

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

784:   Concepts: memory allocation
785: M*/
786: #define PetscCalloc3(m1,r1,m2,r2,m3,r3)                                 \
787:   (PetscMalloc3((m1),(r1),(m2),(r2),(m3),(r3))                          \
788:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))))

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

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

797:    Not Collective

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

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

811:    Level: developer

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

815:   Concepts: memory allocation

817: M*/
818: #if !defined(PETSC_USE_MALLOC_COALESCED)
819: #define PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)) || PetscMalloc1((m3),(r3)) || PetscMalloc1((m4),(r4)))
820: #else
821: #define PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4)                           \
822:   ((((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) \
823:    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),0) \
824:    || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0) || (!(m4) ? (*(r4) = 0,0) : 0))
825: #endif

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

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

834:    Not Collective

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

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

848:    Level: developer

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

852:   Concepts: memory allocation

854: M*/
855: #define PetscCalloc4(m1,r1,m2,r2,m3,r3,m4,r4)                           \
856:   (PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4)                                \
857:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
858:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))))

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

863:    Synopsis:
864:     #include <petscsys.h>
865:    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)

867:    Not Collective

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

876:    Output Parameter:
877: +  r1 - memory allocated in first chunk
878: .  r2 - memory allocated in second chunk
879: .  r3 - memory allocated in third chunk
880: .  r4 - memory allocated in fourth chunk
881: -  r5 - memory allocated in fifth chunk

883:    Level: developer

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

887:   Concepts: memory allocation

889: M*/
890: #if !defined(PETSC_USE_MALLOC_COALESCED)
891: #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)))
892: #else
893: #define PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5)      \
894:   ((((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) \
895:    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),*(void**)(r5) = PetscAddrAlign(*(r4)+(m4)),0) \
896:    || (!(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))
897: #endif

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

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

906:    Not Collective

908:    Input Parameter:
909: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
910: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
911: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
912: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
913: -  m5 - number of elements to allocate in 5th 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

922:    Level: developer

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

926:   Concepts: memory allocation

928: M*/
929: #define PetscCalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5)                     \
930:   (PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5)                          \
931:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
932:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))))

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

937:    Synopsis:
938:     #include <petscsys.h>
939:    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)

941:    Not Collective

943:    Input Parameter:
944: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
945: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
946: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
947: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
948: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
949: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

951:    Output Parameter:
952: +  r1 - memory allocated in first chunk
953: .  r2 - memory allocated in second chunk
954: .  r3 - memory allocated in third chunk
955: .  r4 - memory allocated in fourth chunk
956: .  r5 - memory allocated in fifth chunk
957: -  r6 - memory allocated in sixth chunk

959:    Level: developer

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

963:   Concepts: memory allocation

965: M*/
966: #if !defined(PETSC_USE_MALLOC_COALESCED)
967: #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)))
968: #else
969: #define PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) \
970:   ((((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) \
971:    || (*(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) \
972:    || (!(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))
973: #endif

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

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

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)

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

1000:    Level: developer

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

1004:   Concepts: memory allocation
1005: M*/
1006: #define PetscCalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6)               \
1007:   (PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6)                    \
1008:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
1009:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))) || PetscMemzero(*(r6),(m6)*sizeof(**(r6))))

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

1014:    Synopsis:
1015:     #include <petscsys.h>
1016:    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)

1018:    Not Collective

1020:    Input Parameter:
1021: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
1022: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
1023: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
1024: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
1025: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
1026: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
1027: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

1029:    Output Parameter:
1030: +  r1 - memory allocated in first chunk
1031: .  r2 - memory allocated in second chunk
1032: .  r3 - memory allocated in third chunk
1033: .  r4 - memory allocated in fourth chunk
1034: .  r5 - memory allocated in fifth chunk
1035: .  r6 - memory allocated in sixth chunk
1036: -  r7 - memory allocated in seventh chunk

1038:    Level: developer

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

1042:   Concepts: memory allocation

1044: M*/
1045: #if !defined(PETSC_USE_MALLOC_COALESCED)
1046: #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)))
1047: #else
1048: #define PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) \
1049:   ((((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) \
1050:    || (*(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) \
1051:    || (!(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))
1052: #endif

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

1057:    Synopsis:
1058:     #include <petscsys.h>
1059:    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)

1061:    Not Collective

1063:    Input Parameter:
1064: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
1065: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
1066: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
1067: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
1068: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
1069: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
1070: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

1072:    Output Parameter:
1073: +  r1 - memory allocated in first chunk
1074: .  r2 - memory allocated in second chunk
1075: .  r3 - memory allocated in third chunk
1076: .  r4 - memory allocated in fourth chunk
1077: .  r5 - memory allocated in fifth chunk
1078: .  r6 - memory allocated in sixth chunk
1079: -  r7 - memory allocated in seventh chunk

1081:    Level: developer

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

1085:   Concepts: memory allocation
1086: M*/
1087: #define PetscCalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7)         \
1088:   (PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7)              \
1089:    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
1090:    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))) || PetscMemzero(*(r6),(m6)*sizeof(**(r6))) \
1091:    || PetscMemzero(*(r7),(m7)*sizeof(**(r7))))

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

1096:    Synopsis:
1097:     #include <petscsys.h>
1098:    PetscErrorCode PetscNew(type **result)

1100:    Not Collective

1102:    Output Parameter:
1103: .  result - memory allocated, sized to match pointer type

1105:    Level: beginner

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

1109:   Concepts: memory allocation

1111: M*/
1112: #define PetscNew(b)      PetscCalloc1(1,(b))

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

1118:    Synopsis:
1119:     #include <petscsys.h>
1120:    PetscErrorCode PetscNewLog(PetscObject obj,type **result)

1122:    Not Collective

1124:    Input Parameter:
1125: .  obj - object memory is logged to

1127:    Output Parameter:
1128: .  result - memory allocated, sized to match pointer type

1130:    Level: developer

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

1134:   Concepts: memory allocation

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

1139: /*MC
1140:    PetscFree - Frees memory

1142:    Synopsis:
1143:     #include <petscsys.h>
1144:    PetscErrorCode PetscFree(void *memory)

1146:    Not Collective

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

1151:    Level: beginner

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

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

1159:   Concepts: memory allocation

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

1164: /*MC
1165:    PetscFreeVoid - Frees memory

1167:    Synopsis:
1168:     #include <petscsys.h>
1169:    void PetscFreeVoid(void *memory)

1171:    Not Collective

1173:    Input Parameter:
1174: .   memory - memory to free

1176:    Level: beginner

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

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

1182:   Concepts: memory allocation

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


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

1191:    Synopsis:
1192:     #include <petscsys.h>
1193:    PetscErrorCode PetscFree2(void *memory1,void *memory2)

1195:    Not Collective

1197:    Input Parameter:
1198: +   memory1 - memory to free
1199: -   memory2 - 2nd memory to free

1201:    Level: developer

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

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

1207:   Concepts: memory allocation

1209: M*/
1210: #if !defined(PETSC_USE_MALLOC_COALESCED)
1211: #define PetscFree2(m1,m2)   (PetscFree(m2) || PetscFree(m1))
1212: #else
1213: #define PetscFree2(m1,m2)   ((m1) ? ((m2)=0,PetscFree(m1)) : ((m1)=0,PetscFree(m2)))
1214: #endif

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

1219:    Synopsis:
1220:     #include <petscsys.h>
1221:    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)

1223:    Not Collective

1225:    Input Parameter:
1226: +   memory1 - memory to free
1227: .   memory2 - 2nd memory to free
1228: -   memory3 - 3rd memory to free

1230:    Level: developer

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

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

1236:   Concepts: memory allocation

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

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

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

1252:    Not Collective

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

1260:    Level: developer

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

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

1266:   Concepts: memory allocation

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

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

1278:    Synopsis:
1279:     #include <petscsys.h>
1280:    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)

1282:    Not Collective

1284:    Input Parameter:
1285: +   m1 - memory to free
1286: .   m2 - 2nd memory to free
1287: .   m3 - 3rd memory to free
1288: .   m4 - 4th memory to free
1289: -   m5 - 5th memory to free

1291:    Level: developer

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

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

1297:   Concepts: memory allocation

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


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

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

1315:    Not Collective

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


1326:    Level: developer

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

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

1332:   Concepts: memory allocation

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

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

1346:    Synopsis:
1347:     #include <petscsys.h>
1348:    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)

1350:    Not Collective

1352:    Input Parameter:
1353: +   m1 - memory to free
1354: .   m2 - 2nd memory to free
1355: .   m3 - 3rd memory to free
1356: .   m4 - 4th memory to free
1357: .   m5 - 5th memory to free
1358: .   m6 - 6th memory to free
1359: -   m7 - 7th memory to free


1362:    Level: developer

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

1366: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1367:           PetscMalloc7()

1369:   Concepts: memory allocation

1371: M*/
1372: #if !defined(PETSC_USE_MALLOC_COALESCED)
1373: #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   (PetscFree(m7) || PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1374: #else
1375: #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)) : \
1376:                                            ((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)) : \
1377:                                            ((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)) : \
1378:                                                    ((m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m7))))))))
1379: #endif

1381: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,int,const char[],const char[],void**);
1382: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1383: PETSC_EXTERN PetscErrorCode PetscMallocSet(PetscErrorCode (*)(size_t,int,const char[],const char[],void**),PetscErrorCode (*)(void*,int,const char[],const char[]));
1384: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);

1386: /*
1387:     PetscLogDouble variables are used to contain double precision numbers
1388:   that are not used in the numerical computations, but rather in logging,
1389:   timing etc.
1390: */
1391: typedef double PetscLogDouble;
1392: #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE

1394: /*
1395:    Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1396: */
1397: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1398: PETSC_EXTERN PetscErrorCode PetscMallocDumpLog(FILE *);
1399: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1400: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1401: PETSC_EXTERN PetscErrorCode PetscMallocDebug(PetscBool);
1402: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*);
1403: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1404: PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLog(void);
1405: PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLogThreshold(PetscLogDouble);
1406: PETSC_EXTERN PetscErrorCode PetscMallocGetDumpLog(PetscBool*);

1408: /*E
1409:     PetscDataType - Used for handling different basic data types.

1411:    Level: beginner

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

1415: .seealso: PetscBinaryRead(), PetscBinaryWrite(), PetscDataTypeToMPIDataType(),
1416:           PetscDataTypeGetSize()

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

1423: #if defined(PETSC_USE_COMPLEX)
1424: #define  PETSC_SCALAR  PETSC_COMPLEX
1425: #else
1426: #if defined(PETSC_USE_REAL_SINGLE)
1427: #define  PETSC_SCALAR  PETSC_FLOAT
1428: #elif defined(PETSC_USE_REAL___FLOAT128)
1429: #define  PETSC_SCALAR  PETSC___FLOAT128
1430: #else
1431: #define  PETSC_SCALAR  PETSC_DOUBLE
1432: #endif
1433: #endif
1434: #if defined(PETSC_USE_REAL_SINGLE)
1435: #define  PETSC_REAL  PETSC_FLOAT
1436: #elif defined(PETSC_USE_REAL___FLOAT128)
1437: #define  PETSC_REAL  PETSC___FLOAT128
1438: #else
1439: #define  PETSC_REAL  PETSC_DOUBLE
1440: #endif
1441: #define  PETSC_FORTRANADDR  PETSC_LONG

1443: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1444: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1445: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1446: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);

1448: /*
1449:     Basic memory and string operations. These are usually simple wrappers
1450:    around the basic Unix system calls, but a few of them have additional
1451:    functionality and/or error checking.
1452: */
1453: PETSC_EXTERN PetscErrorCode PetscBitMemcpy(void*,PetscInt,const void*,PetscInt,PetscInt,PetscDataType);
1454: PETSC_EXTERN PetscErrorCode PetscMemmove(void*,void *,size_t);
1455: PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1456: PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1457: PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***);
1458: PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1459: PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool  *);
1460: PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool  *);
1461: PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1462: PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1463: PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1464: PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1465: PETSC_EXTERN PetscErrorCode PetscStrncat(char[],const char[],size_t);
1466: PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1467: PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1468: PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1469: PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1470: PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1471: PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1472: PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1473: PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1474: PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*);
1475: PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1476: PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1477: PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1478: PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1479: PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt,const char *const*,char***);
1480: PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt,char***);
1481: PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);

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

1485: /*S
1486:     PetscToken - 'Token' used for managing tokenizing strings

1488:   Level: intermediate

1490: .seealso: PetscTokenCreate(), PetscTokenFind(), PetscTokenDestroy()
1491: S*/
1492: typedef struct _p_PetscToken* PetscToken;

1494: PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1495: PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1496: PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);

1498: PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1499: PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);

1501: /*
1502:    These are MPI operations for MPI_Allreduce() etc
1503: */
1504: PETSC_EXTERN MPI_Op PetscMaxSum_Op;
1505: #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128)
1506: PETSC_EXTERN MPI_Op MPIU_SUM;
1507: #else
1508: #define MPIU_SUM MPI_SUM
1509: #endif
1510: #if defined(PETSC_USE_REAL___FLOAT128)
1511: PETSC_EXTERN MPI_Op MPIU_MAX;
1512: PETSC_EXTERN MPI_Op MPIU_MIN;
1513: #else
1514: #define MPIU_MAX MPI_MAX
1515: #define MPIU_MIN MPI_MIN
1516: #endif
1517: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);

1519: PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1520: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);

1522: /*S
1523:      PetscObject - any PETSc object, PetscViewer, Mat, Vec, KSP etc

1525:    Level: beginner

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

1529: .seealso:  PetscObjectDestroy(), PetscObjectView(), PetscObjectGetName(), PetscObjectSetName(), PetscObjectReference(), PetscObjectDereference()
1530: S*/
1531: typedef struct _p_PetscObject* PetscObject;

1533: /*MC
1534:     PetscObjectId - unique integer Id for a PetscObject

1536:     Level: developer

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

1540: .seealso: PetscObjectState, PetscObjectGetId()
1541: M*/
1542: typedef Petsc64bitInt PetscObjectId;

1544: /*MC
1545:     PetscObjectState - integer state for a PetscObject

1547:     Level: developer

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

1553: .seealso: PetscObjectId, PetscObjectStateGet(), PetscObjectStateIncrease(), PetscObjectStateSet()
1554: M*/
1555: typedef Petsc64bitInt PetscObjectState;

1557: /*S
1558:      PetscFunctionList - Linked list of functions, possibly stored in dynamic libraries, accessed
1559:       by string name

1561:    Level: advanced

1563: .seealso:  PetscFunctionListAdd(), PetscFunctionListDestroy(), PetscOpFlist
1564: S*/
1565: typedef struct _n_PetscFunctionList *PetscFunctionList;

1567: /*E
1568:   PetscFileMode - Access mode for a file.

1570:   Level: beginner

1572: $  FILE_MODE_READ - open a file at its beginning for reading
1573: $  FILE_MODE_WRITE - open a file at its beginning for writing (will create if the file does not exist)
1574: $  FILE_MODE_APPEND - open a file at end for writing
1575: $  FILE_MODE_UPDATE - open a file for updating, meaning for reading and writing
1576: $  FILE_MODE_APPEND_UPDATE - open a file for updating, meaning for reading and writing, at the end

1578: .seealso: PetscViewerFileSetMode()
1579: E*/
1580: typedef enum {FILE_MODE_READ, FILE_MODE_WRITE, FILE_MODE_APPEND, FILE_MODE_UPDATE, FILE_MODE_APPEND_UPDATE} PetscFileMode;
1581: extern const char *const PetscFileModes[];

1583: /*
1584:     Defines PETSc error handling.
1585: */
1586: #include <petscerror.h>

1588: #define PETSC_SMALLEST_CLASSID  1211211
1589: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1590: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1591: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);

1593: /*
1594:    Routines that get memory usage information from the OS
1595: */
1596: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1597: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1598: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1599: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);

1601: PETSC_EXTERN PetscErrorCode PetscInfoAllow(PetscBool ,const char []);
1602: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);

1604: /*
1605:    Initialization of PETSc
1606: */
1607: PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1608: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1609: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1610: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1611: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1612: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1613: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1614: PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1615: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1616: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);

1618: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1619: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);

1621: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1622: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1623: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1624: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);

1626: /*
1627:      These are so that in extern C code we can caste function pointers to non-extern C
1628:    function pointers. Since the regular C++ code expects its function pointers to be C++
1629: */
1630: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1631: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1632: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);

1634: /*
1635:     Functions that can act on any PETSc object.
1636: */
1637: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1638: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1639: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1640: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1641: PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1642: PETSC_EXTERN PetscErrorCode PetscObjectSetPrecision(PetscObject,PetscPrecision);
1643: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1644: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1645: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1646: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1647: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1648: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1649: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1650: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1651: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1652: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt *);
1653: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1654: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1655: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject *);
1656: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void));
1657: #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d))
1658: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1659: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1660: PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1661: PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject,PetscObject);
1662: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt *);

1664: #include <petscviewertypes.h>
1665: #include <petscoptions.h>

1667: PETSC_EXTERN PetscErrorCode PetscObjectsListGetGlobalNumbering(MPI_Comm,PetscInt,PetscObject*,PetscInt*,PetscInt*);

1669: PETSC_EXTERN PetscErrorCode PetscMemoryShowUsage(PetscViewer,const char[]);
1670: PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer,const char[]);
1671: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1672: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1673: #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1674: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1675: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1676: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1677: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1678: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1679: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1680: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1681: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1682: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,PetscObject,const char[]);
1683: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1684: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1685: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1686: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1687: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);

1689: #if defined(PETSC_HAVE_SAWS)
1690: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1691: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1692: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1693: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1694: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1695: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1696: PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1697: PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1698: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1699: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);

1701: #else
1702: #define PetscSAWsBlock()                        0
1703: #define PetscObjectSAWsViewOff(obj)             0
1704: #define PetscObjectSAWsSetBlock(obj,flg)        0
1705: #define PetscObjectSAWsBlock(obj)               0
1706: #define PetscObjectSAWsGrantAccess(obj)         0
1707: #define PetscObjectSAWsTakeAccess(obj)          0
1708: #define PetscStackViewSAWs()                    0
1709: #define PetscStackSAWsViewOff()                 0
1710: #define PetscStackSAWsTakeAccess()
1711: #define PetscStackSAWsGrantAccess()

1713: #endif

1715: typedef void* PetscDLHandle;
1716: typedef enum {PETSC_DL_DECIDE=0,PETSC_DL_NOW=1,PETSC_DL_LOCAL=2} PetscDLMode;
1717: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1718: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1719: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);

1721: #if defined(PETSC_USE_DEBUG)
1722: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1723: #endif
1724: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);

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

1729:    Level: developer

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

1733: .seealso:  PetscObjectListAdd(), PetscObjectListDestroy(), PetscObjectListFind(), PetscObjectCompose(), PetscObjectQuery(), PetscFunctionList
1734: S*/
1735: typedef struct _n_PetscObjectList *PetscObjectList;

1737: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1738: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1739: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1740: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1741: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1742: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);

1744: /*
1745:     Dynamic library lists. Lists of names of routines in objects or in dynamic
1746:   link libraries that will be loaded as needed.
1747: */

1749: #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1750: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1751: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1752: #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1753: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1754: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[]);
1755: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1756: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1757: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);

1759: /*S
1760:      PetscDLLibrary - Linked list of dynamics libraries to search for functions

1762:    Level: advanced

1764: .seealso:  PetscDLLibraryOpen()
1765: S*/
1766: typedef struct _n_PetscDLLibrary *PetscDLLibrary;
1767: PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1768: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1769: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1770: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1771: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1772: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1773: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1774: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);

1776: /*
1777:      Useful utility routines
1778: */
1779: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1780: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1781: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1782: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1783: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1784: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);

1786: /*
1787:     PetscNot - negates a logical type value and returns result as a PetscBool

1789:     Notes: This is useful in cases like
1790: $     int        *a;
1791: $     PetscBool  flag = PetscNot(a)
1792:      where !a would not return a PetscBool because we cannot provide a cast from int to PetscBool in C.
1793: */
1794: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)

1796: /*MC
1797:     PetscHelpPrintf - Prints help messages.

1799:    Synopsis:
1800:     #include <petscsys.h>
1801:      PetscErrorCode (*PetscHelpPrintf)(const char format[],...);

1803:     Not Collective

1805:     Input Parameters:
1806: .   format - the usual printf() format string

1808:    Level: developer

1810:     Fortran Note:
1811:     This routine is not supported in Fortran.

1813:     Concepts: help messages^printing
1814:     Concepts: printing^help messages

1816: .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1817: M*/
1818: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);

1820: /*
1821:      Defines PETSc profiling.
1822: */
1823: #include <petsclog.h>

1825: /*
1826:       Simple PETSc parallel IO for ASCII printing
1827: */
1828: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1829: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1830: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1831: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1832: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1833: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1834: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);

1836: /* These are used internally by PETSc ASCII IO routines*/
1837: #include <stdarg.h>
1838: PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list);
1839: PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE*,const char[],va_list);
1840: PETSC_EXTERN PetscErrorCode PetscVFPrintfDefault(FILE*,const char[],va_list);

1842: #if defined(PETSC_HAVE_MATLAB_ENGINE)
1843: PETSC_EXTERN PetscErrorCode PetscVFPrintf_Matlab(FILE*,const char[],va_list);
1844: #endif

1846: #if defined(PETSC_HAVE_CLOSURES)
1847: PETSC_EXTERN PetscErrorCode PetscVFPrintfSetClosure(int (^)(const char*));
1848: #endif

1850: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1851: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1852: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);

1854: #if defined(PETSC_HAVE_POPEN)
1855: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1856: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*,int*);
1857: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1858: #endif

1860: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1861: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1862: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*);
1863: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1864: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1865: PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1866: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);

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

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

1873:    Level: advanced

1875: .seealso:  PetscObject, PetscContainerCreate()
1876: S*/
1877: PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1878: typedef struct _p_PetscContainer*  PetscContainer;
1879: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void **);
1880: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void *);
1881: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1882: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer *);
1883: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void*));

1885: /*
1886:    For use in debuggers
1887: */
1888: PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1889: PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1890: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1891: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1892: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);

1894: #include <stddef.h>
1895: #include <string.h>             /* for memcpy, memset */
1896: #if defined(PETSC_HAVE_STDLIB_H)
1897: #include <stdlib.h>
1898: #endif

1900: #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1901: #include <xmmintrin.h>
1902: #endif

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

1911:    Not Collective

1913:    Input Parameters:
1914: +  b - pointer to initial memory space
1915: -  n - length (in bytes) of space to copy

1917:    Output Parameter:
1918: .  a - pointer to copy space

1920:    Level: intermediate

1922:    Compile Option:
1923:     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1924:                                   for memory copies on double precision values.
1925:     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1926:                                   for memory copies on double precision values.
1927:     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1928:                                   for memory copies on double precision values.

1930:    Note:
1931:    This routine is analogous to memcpy().

1933:    Developer Note: this is inlined for fastest performance

1935:   Concepts: memory^copying
1936:   Concepts: copying^memory

1938: .seealso: PetscMemmove()

1940: @*/
1941: PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n)
1942: {
1943: #if defined(PETSC_USE_DEBUG)
1944:   size_t al = (size_t) a,bl = (size_t) b;
1945:   size_t nl = (size_t) n;
1947:   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1948:   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1949: #else
1951: #endif
1952:   if (a != b && n > 0) {
1953: #if defined(PETSC_USE_DEBUG)
1954:     if ((al > bl && (al - bl) < nl) || (bl - al) < nl)  SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1955:               or make sure your copy regions and lengths are correct. \n\
1956:               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1957: #endif
1958: #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1959:    if (!(a % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1960:       size_t len = n/sizeof(PetscScalar);
1961: #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1962:       PetscBLASInt   one = 1,blen;
1964:       PetscBLASIntCast(len,&blen);
1965:       PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one));
1966: #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1967:       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1968: #else
1969:       size_t      i;
1970:       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1971:       for (i=0; i<len; i++) y[i] = x[i];
1972: #endif
1973:     } else {
1974:       memcpy((char*)(a),(char*)(b),n);
1975:     }
1976: #else
1977:     memcpy((char*)(a),(char*)(b),n);
1978: #endif
1979:   }
1980:   return(0);
1981: }

1983: /*@C
1984:    PetscMemzero - Zeros the specified memory.

1986:    Not Collective

1988:    Input Parameters:
1989: +  a - pointer to beginning memory location
1990: -  n - length (in bytes) of memory to initialize

1992:    Level: intermediate

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

1998:    Developer Note: this is inlined for fastest performance

2000:    Concepts: memory^zeroing
2001:    Concepts: zeroing^memory

2003: .seealso: PetscMemcpy()
2004: @*/
2005: PETSC_STATIC_INLINE PetscErrorCode  PetscMemzero(void *a,size_t n)
2006: {
2007:   if (n > 0) {
2008: #if defined(PETSC_USE_DEBUG)
2009:     if (!a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer");
2010: #endif
2011: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
2012:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
2013:       size_t      i,len = n/sizeof(PetscScalar);
2014:       PetscScalar *x = (PetscScalar*)a;
2015:       for (i=0; i<len; i++) x[i] = 0.0;
2016:     } else {
2017: #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
2018:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
2019:       PetscInt len = n/sizeof(PetscScalar);
2020:       fortranzero_(&len,(PetscScalar*)a);
2021:     } else {
2022: #endif
2023: #if defined(PETSC_PREFER_BZERO)
2024:       bzero((char *)a,n);
2025: #else
2026:       memset((char*)a,0,n);
2027: #endif
2028: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
2029:     }
2030: #endif
2031:   }
2032:   return 0;
2033: }

2035: /*MC
2036:    PetscPrefetchBlock - Prefetches a block of memory

2038:    Synopsis:
2039:     #include <petscsys.h>
2040:     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)

2042:    Not Collective

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

2050:    Level: developer

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

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

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

2065:    Concepts: memory
2066: M*/
2067: #define PetscPrefetchBlock(a,n,rw,t) do {                               \
2068:     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
2069:     for ( ; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
2070:   } while (0)

2072: /*
2073:       Determine if some of the kernel computation routines use
2074:    Fortran (rather than C) for the numerical calculations. On some machines
2075:    and compilers (like complex numbers) the Fortran version of the routines
2076:    is faster than the C/C++ versions. The flag --with-fortran-kernels
2077:    should be used with ./configure to turn these on.
2078: */
2079: #if defined(PETSC_USE_FORTRAN_KERNELS)

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

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

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

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

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

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

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

2109: #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
2110: #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
2111: #endif

2113: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
2114: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
2115: #endif

2117: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
2118: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
2119: #endif

2121: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
2122: #define PETSC_USE_FORTRAN_KERNEL_MDOT
2123: #endif

2125: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
2126: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
2127: #endif

2129: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
2130: #define PETSC_USE_FORTRAN_KERNEL_AYPX
2131: #endif

2133: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
2134: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
2135: #endif

2137: #endif

2139: /*
2140:     Macros for indicating code that should be compiled with a C interface,
2141:    rather than a C++ interface. Any routines that are dynamically loaded
2142:    (such as the PCCreate_XXX() routines) must be wrapped so that the name
2143:    mangler does not change the functions symbol name. This just hides the
2144:    ugly extern "C" {} wrappers.
2145: */
2146: #if defined(__cplusplus)
2147: #define EXTERN_C_BEGIN extern "C" {
2148: #define EXTERN_C_END }
2149: #else
2150: #define EXTERN_C_BEGIN
2151: #define EXTERN_C_END
2152: #endif

2154: /* --------------------------------------------------------------------*/

2156: /*MC
2157:     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
2158:         communication

2160:    Level: beginner

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

2164: .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
2165: M*/

2167: /*MC
2168:     PetscScalar - PETSc type that represents either a double precision real number, a double precision
2169:        complex number, a single precision real number, a long double or an int - if the code is configured
2170:        with --with-scalar-type=real,complex --with-precision=single,double,__float128

2172:    Level: beginner

2174: .seealso: PetscReal, MPIU_SCALAR, PetscInt, MPIU_REAL
2175: M*/

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

2180:    Synopsis:
2181:    #include <petscsys.h>
2182:    PetscComplex number = 1. + 2.*PETSC_i;

2184:    Level: beginner

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

2189: .seealso: PetscReal, PetscComplex, MPIU_COMPLEX, PetscInt, PETSC_i
2190: M*/

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

2195:    Level: beginner

2197: .seealso: PetscScalar
2198: M*/

2200: /*MC
2201:     MPIU_SCALAR - MPI datatype corresponding to PetscScalar

2203:    Level: beginner

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

2208: .seealso: PetscReal, PetscScalar, MPIU_INT
2209: M*/

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

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

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

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

2236:    Not Collective

2238:    Input Parameter:
2239: .     a - the PetscInt value

2241:    Output Parameter:
2242: .     b - the resulting PetscBLASInt value

2244:    Level: advanced

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

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

2264:    Not Collective

2266:    Input Parameter:
2267: .     a - the PetscInt value

2269:    Output Parameter:
2270: .     b - the resulting PetscMPIInt value

2272:    Level: advanced

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

2286: #define PetscIntMult64bit(a,b)   ((Petsc64bitInt)(a))*((Petsc64bitInt)(b))

2290: /*@C

2292:    PetscRealIntMultTruncate - Computes the product of a positive PetscReal and a positive PetscInt and truncates the value to slightly less than the maximal possible value

2294:    Not Collective

2296:    Input Parameter:
2297: +     a - the PetscReal value
2298: -     b - the second value

2300:    Output Parameter:
2301: .     c - the result as a PetscInt value

2303:    Use PetscIntMult64bit() to compute the product of two PetscInt as a Petsc64bitInt
2304:    Use PetscIntMultTruncate() to compute the product of two positive PetscInt and truncate to fit a PetscInt
2305:    Use PetscIntMultError() to compute the product of two PetscInt if you wish to generate an error if the result will not fit in a PetscInt

2307:    Developers Note: We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks.

2309:    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.

2311:    Level: advanced

2313: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64()
2314: @*/
2315: PETSC_STATIC_INLINE PetscInt PetscRealIntMultTruncate(PetscReal a,PetscInt b)
2316: {
2317:   Petsc64bitInt r;

2319:   r  =  (Petsc64bitInt) (a*(PetscReal)b);
2320:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2321:   return (PetscInt) r;
2322: }

2326: /*@C

2328:    PetscIntMultTruncate - Computes the product of two positive PetscInt and truncates the value to slightly less than the maximal possible value

2330:    Not Collective

2332:    Input Parameter:
2333: +     a - the PetscInt value
2334: -     b - the second value

2336:    Output Parameter:
2337: .     c - the result as a PetscInt value

2339:    Use PetscIntMult64bit() to compute the product of two PetscInt as a Petsc64bitInt
2340:    Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2341:    Use PetscIntMultError() to compute the product of two PetscInt if you wish to generate an error if the result will not fit in a PetscInt

2343:    Developers Note: We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks.

2345:    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.

2347:    Level: advanced

2349: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64()
2350: @*/
2351: PETSC_STATIC_INLINE PetscInt PetscIntMultTruncate(PetscInt a,PetscInt b)
2352: {
2353:   Petsc64bitInt r;

2355:   r  =  PetscIntMult64bit(a,b);
2356:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2357:   return (PetscInt) r;
2358: }

2362: /*@C

2364:    PetscIntSumTruncate - Computes the sum of two positive PetscInt and truncates the value to slightly less than the maximal possible value

2366:    Not Collective

2368:    Input Parameter:
2369: +     a - the PetscInt value
2370: -     b - the second value

2372:    Output Parameter:
2373: .     c - the result as a PetscInt value

2375:    Use PetscIntMult64bit() to compute the product of two PetscInt as a Petsc64bitInt
2376:    Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2377:    Use PetscIntMultError() to compute the product of two PetscInt if you wish to generate an error if the result will not fit in a PetscInt

2379:    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.

2381:    Level: advanced

2383: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64()
2384: @*/
2385: PETSC_STATIC_INLINE PetscInt PetscIntSumTruncate(PetscInt a,PetscInt b)
2386: {
2387:   Petsc64bitInt r;

2389:   r  =  ((Petsc64bitInt)a) + ((Petsc64bitInt)b);
2390:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2391:   return (PetscInt) r;
2392: }

2396: /*@C

2398:    PetscIntMultError - Computes the product of two positive PetscInt and generates an error with overflow.

2400:    Not Collective

2402:    Input Parameter:
2403: +     a - the PetscInt value
2404: -     b - the second value

2406:    Output Parameter:ma
2407: .     c - the result as a PetscInt value

2409:    Use PetscIntMult64bit() to compute the product of two 32 bit PetscInt and store in a Petsc64bitInt
2410:    Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt

2412:    Developers Note: We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks.

2414:    Level: advanced

2416: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64()
2417: @*/
2418: PETSC_STATIC_INLINE PetscErrorCode PetscIntMultError(PetscInt a,PetscInt b,PetscInt *result)
2419: {
2420:   Petsc64bitInt r;

2423:   r  =  PetscIntMult64bit(a,b);
2424: #if !defined(PETSC_USE_64BIT_INDICES)
2425:   if (r > PETSC_MAX_INT) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Product of two integer %d %d overflow, you must ./configure PETSc with --with-64-bit-indices for the case you are running",a,b);
2426: #endif
2427:   if (result) *result = (PetscInt) r;
2428:   return(0);
2429: }

2433: /*@C

2435:    PetscIntSumError - Computes the product of two positive PetscInt and generates an error with overflow.

2437:    Not Collective

2439:    Input Parameter:
2440: +     a - the PetscInt value
2441: -     b - the second value

2443:    Output Parameter:ma
2444: .     c - the result as a PetscInt value

2446:    Use PetscIntMult64bit() to compute the product of two 32 bit PetscInt and store in a Petsc64bitInt
2447:    Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt

2449:    Level: advanced

2451: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64()
2452: @*/
2453: PETSC_STATIC_INLINE PetscErrorCode PetscIntSumError(PetscInt a,PetscInt b,PetscInt *result)
2454: {
2455:   Petsc64bitInt r;

2458:   r  =  ((Petsc64bitInt)a) + ((Petsc64bitInt)b);
2459: #if !defined(PETSC_USE_64BIT_INDICES)
2460:   if (r > PETSC_MAX_INT) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Sum of two integer %d %d overflow, you must ./configure PETSc with --with-64-bit-indices for the case you are running",a,b);
2461: #endif
2462:   if (result) *result = (PetscInt) r;
2463:   return(0);
2464: }
2465: 
2466: /*
2467:      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2468: */
2469: #if defined(hz)
2470: #undef hz
2471: #endif

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


2476: #if defined(PETSC_HAVE_LIMITS_H)
2477: #include <limits.h>
2478: #endif
2479: #if defined(PETSC_HAVE_SYS_PARAM_H)
2480: #include <sys/param.h>
2481: #endif
2482: #if defined(PETSC_HAVE_SYS_TYPES_H)
2483: #include <sys/types.h>
2484: #endif
2485: #if defined(MAXPATHLEN)
2486: #  define PETSC_MAX_PATH_LEN     MAXPATHLEN
2487: #elif defined(MAX_PATH)
2488: #  define PETSC_MAX_PATH_LEN     MAX_PATH
2489: #elif defined(_MAX_PATH)
2490: #  define PETSC_MAX_PATH_LEN     _MAX_PATH
2491: #else
2492: #  define PETSC_MAX_PATH_LEN     4096
2493: #endif

2495: /*MC

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

2499: $    1) classic Fortran 77 style
2500: $#include "petsc/finclude/petscXXX.h" to work with material from the XXX component of PETSc
2501: $       XXX variablename
2502: $      You cannot use this approach if you wish to use the Fortran 90 specific PETSc routines
2503: $      which end in F90; such as VecGetArrayF90()
2504: $
2505: $    2) classic Fortran 90 style
2506: $#include "petsc/finclude/petscXXX.h"
2507: $#include "petsc/finclude/petscXXX.h90" to work with material from the XXX component of PETSc
2508: $       XXX variablename
2509: $
2510: $    3) Using Fortran modules
2511: $#include "petsc/finclude/petscXXXdef.h"
2512: $         use petscXXXX
2513: $       XXX variablename
2514: $
2515: $    4) Use Fortran modules and Fortran data types for PETSc types
2516: $#include "petsc/finclude/petscXXXdef.h"
2517: $         use petscXXXX
2518: $       type(XXX) variablename
2519: $      To use this approach you must ./configure PETSc with the additional
2520: $      option --with-fortran-datatypes You cannot use the type(XXX) declaration approach without using Fortran modules

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

2524: $    3a) skip the #include BUT you cannot use any PETSc data type names like Vec, Mat, PetscInt, PetscErrorCode etc
2525: $        and you must declare the variables as integer, for example
2526: $        integer variablename
2527: $
2528: $    4a) skip the #include, you use the object types like type(Vec) type(Mat) but cannot use the data type
2529: $        names like PetscErrorCode, PetscInt etc. again for those you must use integer

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

2534:    Fortran type checking with interfaces is strick, this means you cannot pass a scalar value when an array value
2535: is expected (even though it is legal Fortran). For example when setting a single value in a matrix with MatSetValues()
2536: you cannot have something like
2537: $      PetscInt row,col
2538: $      PetscScalar val
2539: $        ...
2540: $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)
2541: You must instead have
2542: $      PetscInt row(1),col(1)
2543: $      PetscScalar val(1)
2544: $        ...
2545: $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)


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

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

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

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

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

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

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

2570:     Level: beginner

2572: M*/

2574: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2575: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2576: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2577: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2578: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2579: PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2580: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);

2582: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2583: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2584: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2585: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2586: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2587: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2588: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt*,PetscInt*,PetscInt*);
2589: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2590: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2591: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2592: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2593: PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt,PetscInt[],void*,size_t,void*);
2594: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2595: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2596: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2597: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2598: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2599: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2600: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt*,const PetscInt*,PetscInt,const PetscInt*,const PetscInt*,PetscInt*,PetscInt**,PetscInt**);
2601: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt,const PetscInt*,PetscInt,const PetscInt*,PetscInt*,PetscInt**);
2602: PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt,const PetscMPIInt[],PetscInt,const PetscMPIInt[],PetscInt*,PetscMPIInt**);

2604: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2605: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);

2607: /*J
2608:     PetscRandomType - String with the name of a PETSc randomizer

2610:    Level: beginner

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

2615: .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2616: J*/
2617: typedef const char* PetscRandomType;
2618: #define PETSCRAND       "rand"
2619: #define PETSCRAND48     "rand48"
2620: #define PETSCSPRNG      "sprng"
2621: #define PETSCRANDER48   "rander48"

2623: /* Logging support */
2624: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;

2626: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);

2628: /*S
2629:      PetscRandom - Abstract PETSc object that manages generating random numbers

2631:    Level: intermediate

2633:   Concepts: random numbers

2635: .seealso:  PetscRandomCreate(), PetscRandomGetValue(), PetscRandomType
2636: S*/
2637: typedef struct _p_PetscRandom*   PetscRandom;

2639: /* Dynamic creation and loading functions */
2640: PETSC_EXTERN PetscFunctionList PetscRandomList;

2642: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2643: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2644: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2645: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType*);
2646: PETSC_STATIC_INLINE PetscErrorCode PetscRandomViewFromOptions(PetscRandom A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
2647: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);

2649: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2650: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2651: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2652: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2653: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2654: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2655: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2656: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2657: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);

2659: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2660: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2661: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2662: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2663: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2664: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2665: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);
2666: PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2667: PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);

2669: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscDataType);
2670: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscDataType);
2671: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,void*,PetscInt,PetscDataType,PetscBool );
2672: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,void*,PetscInt,PetscDataType,PetscBool );
2673: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2674: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2675: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2676: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2677: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2678: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2679: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2680: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);

2682: /*
2683:    In binary files variables are stored using the following lengths,
2684:   regardless of how they are stored in memory on any one particular
2685:   machine. Use these rather then sizeof() in computing sizes for
2686:   PetscBinarySeek().
2687: */
2688: #define PETSC_BINARY_INT_SIZE   (32/8)
2689: #define PETSC_BINARY_FLOAT_SIZE  (32/8)
2690: #define PETSC_BINARY_CHAR_SIZE  (8/8)
2691: #define PETSC_BINARY_SHORT_SIZE  (16/8)
2692: #define PETSC_BINARY_DOUBLE_SIZE  (64/8)
2693: #define PETSC_BINARY_SCALAR_SIZE  sizeof(PetscScalar)

2695: /*E
2696:   PetscBinarySeekType - argument to PetscBinarySeek()

2698:   Level: advanced

2700: .seealso: PetscBinarySeek(), PetscBinarySynchronizedSeek()
2701: E*/
2702: typedef enum {PETSC_BINARY_SEEK_SET = 0,PETSC_BINARY_SEEK_CUR = 1,PETSC_BINARY_SEEK_END = 2} PetscBinarySeekType;
2703: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2704: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2705: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);

2707: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2708: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool );
2709: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2710: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2711: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2712: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);

2714: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2715: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2716: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2717: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2718: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2719: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt*,const void*,PetscMPIInt*,PetscMPIInt**,void*)
2720:   PetscAttrMPIPointerWithType(6,3);
2721: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2722:                                                     PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2723:                                                     PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2724:   PetscAttrMPIPointerWithType(6,3);
2725: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2726:                                                        MPI_Request**,MPI_Request**,
2727:                                                        PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2728:                                                        PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2729:   PetscAttrMPIPointerWithType(6,3);

2731: /*E
2732:     PetscBuildTwoSidedType - algorithm for setting up two-sided communication

2734: $  PETSC_BUILDTWOSIDED_ALLREDUCE - classical algorithm using an MPI_Allreduce with
2735: $      a buffer of length equal to the communicator size. Not memory-scalable due to
2736: $      the large reduction size. Requires only MPI-1.
2737: $  PETSC_BUILDTWOSIDED_IBARRIER - nonblocking algorithm based on MPI_Issend and MPI_Ibarrier.
2738: $      Proved communication-optimal in Hoefler, Siebert, and Lumsdaine (2010). Requires MPI-3.
2739: $  PETSC_BUILDTWOSIDED_REDSCATTER - similar to above, but use more optimized function
2740: $      that only communicates the part of the reduction that is necessary.  Requires MPI-2.

2742:    Level: developer

2744: .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedSetType(), PetscCommBuildTwoSidedGetType()
2745: E*/
2746: typedef enum {
2747:   PETSC_BUILDTWOSIDED_NOTSET = -1,
2748:   PETSC_BUILDTWOSIDED_ALLREDUCE = 0,
2749:   PETSC_BUILDTWOSIDED_IBARRIER = 1,
2750:   PETSC_BUILDTWOSIDED_REDSCATTER = 2
2751:   /* Updates here must be accompanied by updates in finclude/petscsys.h and the string array in mpits.c */
2752: } PetscBuildTwoSidedType;
2753: PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2754: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2755: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);

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

2759: /*E
2760:   InsertMode - Whether entries are inserted or added into vectors or matrices

2762:   Level: beginner

2764: .seealso: VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2765:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(),
2766:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd()
2767: E*/
2768:  typedef enum {NOT_SET_VALUES, INSERT_VALUES, ADD_VALUES, MAX_VALUES, INSERT_ALL_VALUES, ADD_ALL_VALUES, INSERT_BC_VALUES, ADD_BC_VALUES} InsertMode;

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

2773:     Level: beginner

2775: .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2776:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), ADD_VALUES,
2777:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES

2779: M*/

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

2785:     Level: beginner

2787: .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2788:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), INSERT_VALUES,
2789:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES

2791: M*/

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

2796:     Level: beginner

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

2800: M*/

2802: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);

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

2807: /*S
2808:    PetscSubcomm - A decomposition of an MPI communicator into subcommunicators

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

2814:    Sample Usage:
2815:        PetscSubcommCreate()
2816:        PetscSubcommSetNumber()
2817:        PetscSubcommSetType(PETSC_SUBCOMM_INTERLACED);
2818:        ccomm = PetscSubcommChild()
2819:        PetscSubcommDestroy()

2821:    Level: advanced

2823:    Concepts: communicator, create

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

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

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


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

2840: S*/
2841: typedef struct _n_PetscSubcomm* PetscSubcomm;

2843: struct _n_PetscSubcomm {
2844:   MPI_Comm         parent;           /* parent communicator */
2845:   MPI_Comm         dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2846:   MPI_Comm         child;            /* the sub-communicator */
2847:   PetscMPIInt      n;                /* num of subcommunicators under the parent communicator */
2848:   PetscMPIInt      color;            /* color of processors belong to this communicator */
2849:   PetscMPIInt      *subsize;         /* size of subcommunicator[color] */
2850:   PetscSubcommType type;
2851:   char             *subcommprefix;
2852: };

2854: PETSC_STATIC_INLINE MPI_Comm PetscSubcommParent(PetscSubcomm scomm) {return scomm->parent;}
2855: PETSC_STATIC_INLINE MPI_Comm PetscSubcommChild(PetscSubcomm scomm) {return scomm->child;}
2856: PETSC_STATIC_INLINE MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm) {return scomm->dupparent;}
2857: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2858: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2859: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2860: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2861: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2862: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2863: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2864: PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm,const char[]);

2866: /*S
2867:    PetscSegBuffer - a segmented extendable buffer

2869:    Level: developer

2871: .seealso: PetscSegBufferCreate(), PetscSegBufferGet(), PetscSegBufferExtract(), PetscSegBufferDestroy()
2872: S*/
2873: typedef struct _n_PetscSegBuffer *PetscSegBuffer;
2874: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2875: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2876: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2877: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2878: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2879: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2880: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2881: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);

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

2888: typedef struct _n_PetscOptionsHelpPrinted *PetscOptionsHelpPrinted;
2889: extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2890: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted*);
2891: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted*);
2892: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted,const char*,const char*,PetscBool*);

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

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

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

2906:    Level: intermediate

2908:      Options Database:
2909: .     -citations [filenmae]   - print out the bibtex entries for the given computation
2910: @*/
2911: PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2912: {
2913:   size_t         len;
2914:   char           *vstring;

2918:   if (set && *set) return(0);
2919:   PetscStrlen(cit,&len);
2920:   PetscSegBufferGet(PetscCitationsList,len,&vstring);
2921:   PetscMemcpy(vstring,cit,len);
2922:   if (set) *set = PETSC_TRUE;
2923:   return(0);
2924: }

2926: PETSC_EXTERN PetscErrorCode PetscURLShorten(const char[],char[],size_t);
2927: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm,char[],char[],size_t);
2928: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm,const char[],char[],size_t);
2929: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm,const char[],const char []);

2931: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm,char[],char[],size_t);
2932: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm,const char[],char[],char[],size_t);

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

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


2940: #if defined(PETSC_USE_DEBUG)
2941: /*
2942:    Verify that all processes in the communicator have called this from the same line of code
2943:  */
2944: PETSC_EXTERN PetscErrorCode PetscAllreduceBarrierCheck(MPI_Comm,PetscMPIInt,int,const char*,const char *);
2945: #define MPIU_Allreduce(a,b,c,d,e,fcomm) (PetscAllreduceBarrierCheck(fcomm,c,__LINE__,__FUNCT__,__FILE__) || MPI_Allreduce(a,b,c,d,e,fcomm))
2946: #else
2947: #define MPIU_Allreduce(a,b,c,d,e,fcomm) MPI_Allreduce(a,b,c,d,e,fcomm)
2948: #endif

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

2954: #endif