Actual source code: petscsys.h

  1: /*
  2:    This is the main PETSc include file (for C and C++).  It is included by all
  3:    other PETSc include files, so it almost never has to be specifically included.
  4: */
  5: #if !defined(PETSCSYS_H)
  6: #define PETSCSYS_H

  8: /* ========================================================================== */
  9: /*
 10:    petscconf.h is contained in ${PETSC_ARCH}/include/petscconf.h it is
 11:    found automatically by the compiler due to the -I${PETSC_DIR}/${PETSC_ARCH}/include that
 12:    PETSc's makefiles add to the compiler rules.
 13:    For --prefix installs the directory ${PETSC_ARCH} does not exist and petscconf.h is in the same
 14:    directory as the other PETSc include files.
 15: */
 16: #include <petscconf.h>
 17: #include <petscconf_poison.h>
 18: #include <petscfix.h>
 19: #include <petscmacros.h>

 21: /* placeholder defines */
 22: #if PetscHasAttribute(format)
 23: #  define PETSC_ATTRIBUTE_FORMAT(strIdx,vaArgIdx)
 24: #else
 25: #  define PETSC_ATTRIBUTE_FORMAT(strIdx,vaArgIdx)
 26: #endif

 28: #if defined(PETSC_DESIRE_FEATURE_TEST_MACROS)
 29: /*
 30:    Feature test macros must be included before headers defined by IEEE Std 1003.1-2001
 31:    We only turn these in PETSc source files that require them by setting PETSC_DESIRE_FEATURE_TEST_MACROS
 32: */
 33: #  if defined(PETSC__POSIX_C_SOURCE_200112L) && !defined(_POSIX_C_SOURCE)
 34: #    define _POSIX_C_SOURCE 200112L
 35: #  endif
 36: #  if defined(PETSC__BSD_SOURCE) && !defined(_BSD_SOURCE)
 37: #    define _BSD_SOURCE
 38: #  endif
 39: #  if defined(PETSC__DEFAULT_SOURCE) && !defined(_DEFAULT_SOURCE)
 40: #    define _DEFAULT_SOURCE
 41: #  endif
 42: #  if defined(PETSC__GNU_SOURCE) && !defined(_GNU_SOURCE)
 43: #    define _GNU_SOURCE
 44: #  endif
 45: #endif

 47: #include <petscsystypes.h>

 49: /* ========================================================================== */

 51: /*
 52:     Defines the interface to MPI allowing the use of all MPI functions.

 54:     PETSc does not use the C++ binding of MPI at ALL. The following flag
 55:     makes sure the C++ bindings are not included. The C++ bindings REQUIRE
 56:     putting mpi.h before ANY C++ include files, we cannot control this
 57:     with all PETSc users. Users who want to use the MPI C++ bindings can include
 58:     mpicxx.h directly in their code
 59: */
 60: #if !defined(MPICH_SKIP_MPICXX)
 61: #  define MPICH_SKIP_MPICXX 1
 62: #endif
 63: #if !defined(OMPI_SKIP_MPICXX)
 64: #  define OMPI_SKIP_MPICXX 1
 65: #endif
 66: #if defined(PETSC_HAVE_MPIUNI)
 67: #  include <petsc/mpiuni/mpi.h>
 68: #else
 69: #  include <mpi.h>
 70: #endif

 72: /*
 73:    Perform various sanity checks that the correct mpi.h is being included at compile time.
 74:    This usually happens because
 75:       * either an unexpected mpi.h is in the default compiler path (i.e. in /usr/include) or
 76:       * an extra include path -I/something (which contains the unexpected mpi.h) is being passed to the compiler
 77: */
 78: #if defined(PETSC_HAVE_MPIUNI)
 79: #  if !defined(MPIUNI_H)
 80: #    error "PETSc was configured with --with-mpi=0 but now appears to be compiling using a different mpi.h"
 81: #  endif
 82: #elif defined(PETSC_HAVE_I_MPI_NUMVERSION)
 83: #  if !defined(I_MPI_NUMVERSION)
 84: #    error "PETSc was configured with I_MPI but now appears to be compiling using a non-I_MPI mpi.h"
 85: #  elif I_MPI_NUMVERSION != PETSC_HAVE_I_MPI_NUMVERSION
 86: #    error "PETSc was configured with one I_MPI mpi.h version but now appears to be compiling using a different I_MPI mpi.h version"
 87: #  endif
 88: #elif defined(PETSC_HAVE_MVAPICH2_NUMVERSION)
 89: #  if !defined(MVAPICH2_NUMVERSION)
 90: #    error "PETSc was configured with MVAPICH2 but now appears to be compiling using a non-MVAPICH2 mpi.h"
 91: #  elif MVAPICH2_NUMVERSION != PETSC_HAVE_MVAPICH2_NUMVERSION
 92: #    error "PETSc was configured with one MVAPICH2 mpi.h version but now appears to be compiling using a different MVAPICH2 mpi.h version"
 93: #  endif
 94: #elif defined(PETSC_HAVE_MPICH_NUMVERSION)
 95: #  if !defined(MPICH_NUMVERSION) || defined(MVAPICH2_NUMVERSION) || defined(I_MPI_NUMVERSION)
 96: #    error "PETSc was configured with MPICH but now appears to be compiling using a non-MPICH mpi.h"
 97: #  elif (MPICH_NUMVERSION/100000000 != PETSC_HAVE_MPICH_NUMVERSION/100000000) || (MPICH_NUMVERSION/100000 < PETSC_HAVE_MPICH_NUMVERSION/100000) || (MPICH_NUMVERSION/100000 == PETSC_HAVE_MPICH_NUMVERSION/100000 && MPICH_NUMVERSION%100000/1000 < PETSC_HAVE_MPICH_NUMVERSION%100000/1000)
 98: #    error "PETSc was configured with one MPICH mpi.h version but now appears to be compiling using a different MPICH mpi.h version"
 99: #  endif
100: #elif defined(PETSC_HAVE_OMPI_MAJOR_VERSION)
101: #  if !defined(OMPI_MAJOR_VERSION)
102: #    error "PETSc was configured with OpenMPI but now appears to be compiling using a non-OpenMPI mpi.h"
103: #  elif (OMPI_MAJOR_VERSION != PETSC_HAVE_OMPI_MAJOR_VERSION) || (OMPI_MINOR_VERSION < PETSC_HAVE_OMPI_MINOR_VERSION) || (OMPI_MINOR_VERSION == PETSC_HAVE_OMPI_MINOR_VERSION && OMPI_RELEASE_VERSION < PETSC_HAVE_OMPI_RELEASE_VERSION)
104: #    error "PETSc was configured with one OpenMPI mpi.h version but now appears to be compiling using a different OpenMPI mpi.h version"
105: #  endif
106: #  define PETSC_MPI_COMM_FMT "p"
107: #  define PETSC_MPI_WIN_FMT  "p"
108: #elif defined(PETSC_HAVE_MSMPI_VERSION)
109: #  if !defined(MSMPI_VER)
110: #    error "PETSc was configured with MSMPI but now appears to be compiling using a non-MSMPI mpi.h"
111: #  elif (MSMPI_VER != PETSC_HAVE_MSMPI_VERSION)
112: #    error "PETSc was configured with one MSMPI mpi.h version but now appears to be compiling using a different MSMPI mpi.h version"
113: #  endif
114: #elif defined(OMPI_MAJOR_VERSION) || defined(MPICH_NUMVERSION) || defined(MSMPI_VER)
115: #  error "PETSc was configured with undetermined MPI - but now appears to be compiling using any of OpenMPI, MS-MPI or a MPICH variant"
116: #endif

118: /* Format specifier for printing MPI_Comm (most implementations use 'int' as type) */
119: #if !defined(PETSC_MPI_COMM_FMT)
120: #  define PETSC_MPI_COMM_FMT "d"
121: #endif

123: #if !defined(PETSC_MPI_WIN_FMT)
124: #  define PETSC_MPI_WIN_FMT "d"
125: #endif

127: /*
128:     Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler
129:     see the top of mpicxx.h in the MPICH2 distribution.
130: */
131: #include <stdio.h>

133: /* MSMPI on 32bit windows requires this yukky hack - that breaks MPI standard compliance */
134: #if !defined(MPIAPI)
135: #define MPIAPI
136: #endif

138: /*
139:    Support for Clang (>=3.2) matching type tag arguments with void* buffer types.
140:    This allows the compiler to detect cases where the MPI datatype argument passed to a MPI routine
141:    does not match the actual type of the argument being passed in
142: */
143: #if defined(__has_attribute) && defined(works_with_const_which_is_not_true)
144: #  if __has_attribute(argument_with_type_tag) && __has_attribute(pointer_with_type_tag) && __has_attribute(type_tag_for_datatype)
145: #    define PetscAttrMPIPointerWithType(bufno,typeno) __attribute__((pointer_with_type_tag(MPI,bufno,typeno)))
146: #    define PetscAttrMPITypeTag(type)                 __attribute__((type_tag_for_datatype(MPI,type)))
147: #    define PetscAttrMPITypeTagLayoutCompatible(type) __attribute__((type_tag_for_datatype(MPI,type,layout_compatible)))
148: #  endif
149: #endif
150: #if !defined(PetscAttrMPIPointerWithType)
151: #  define PetscAttrMPIPointerWithType(bufno,typeno)
152: #  define PetscAttrMPITypeTag(type)
153: #  define PetscAttrMPITypeTagLayoutCompatible(type)
154: #endif

156: PETSC_EXTERN MPI_Datatype MPIU_ENUM PetscAttrMPITypeTag(PetscEnum);
157: PETSC_EXTERN MPI_Datatype MPIU_BOOL PetscAttrMPITypeTag(PetscBool);

159: /*MC
160:    MPIU_INT - MPI datatype corresponding to PetscInt

162:    Notes:
163:    In MPI calls that require an MPI datatype that matches a PetscInt or array of PetscInt values, pass this value.

165:    Level: beginner

167: .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_SCALAR, MPIU_COMPLEX
168: M*/

170: PETSC_EXTERN MPI_Datatype MPIU_FORTRANADDR;

172: #if defined(PETSC_USE_64BIT_INDICES)
173: #  define MPIU_INT MPIU_INT64
174: #  define PetscInt_FMT PetscInt64_FMT
175: #else
176: #  define MPIU_INT MPI_INT
177: #  define PetscInt_FMT "d"
178: #endif

180: /*
181:     For the rare cases when one needs to send a size_t object with MPI
182: */
183: PETSC_EXTERN MPI_Datatype MPIU_SIZE_T;

185: /*
186:       You can use PETSC_STDOUT as a replacement of stdout. You can also change
187:     the value of PETSC_STDOUT to redirect all standard output elsewhere
188: */
189: PETSC_EXTERN FILE* PETSC_STDOUT;

191: /*
192:       You can use PETSC_STDERR as a replacement of stderr. You can also change
193:     the value of PETSC_STDERR to redirect all standard error elsewhere
194: */
195: PETSC_EXTERN FILE* PETSC_STDERR;

197: /* PetscPragmaSIMD - from CeedPragmaSIMD */

199: #if defined(__NEC__)
200: #  define PetscPragmaSIMD _Pragma("_NEC ivdep")
201: #elif defined(__INTEL_COMPILER) && !defined(_WIN32)
202: #  define PetscPragmaSIMD _Pragma("vector")
203: #elif defined(__GNUC__) && __GNUC__ >= 5 && !defined(__PGI)
204: #  define PetscPragmaSIMD _Pragma("GCC ivdep")
205: #elif defined(_OPENMP) && _OPENMP >= 201307
206: #  if defined(_MSC_VER)
207: #    define PetscPragmaSIMD __pragma(omp simd)
208: #  else
209: #    define PetscPragmaSIMD _Pragma("omp simd")
210: #  endif
211: #elif defined(PETSC_HAVE_CRAY_VECTOR)
212: #  define PetscPragmaSIMD _Pragma("_CRI ivdep")
213: #else
214: #  define PetscPragmaSIMD
215: #endif

217: /*
218:     Declare extern C stuff after including external header files
219: */

221: PETSC_EXTERN const char *const PetscBools[];

223: PETSC_EXTERN PetscBool PETSC_RUNNING_ON_VALGRIND;
224: /*
225:     Defines elementary mathematics functions and constants.
226: */
227: #include <petscmath.h>

229: PETSC_EXTERN const char *const PetscCopyModes[];

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

234:    Level: beginner

236:    Note:
237:    Accepted by many PETSc functions to not set a parameter and instead use
238:           some default

240:    Fortran Notes:
241:    This macro does not exist in Fortran; you must use PETSC_NULL_INTEGER,
242:           PETSC_NULL_DOUBLE_PRECISION etc

244: .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_DETERMINE

246: M*/
247: #define PETSC_IGNORE NULL

249: /* This is deprecated */
250: #define PETSC_NULL NULL

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

256:    Level: beginner

258: .seealso: PETSC_DEFAULT, PETSC_IGNORE, PETSC_DETERMINE

260: M*/
261: #define PETSC_DECIDE -1

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

267:    Level: beginner

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

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

275: M*/
276: #define PETSC_DETERMINE PETSC_DECIDE

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

282:    Level: beginner

284:    Fortran Notes:
285:    You need to use PETSC_DEFAULT_INTEGER or PETSC_DEFAULT_REAL.

287: .seealso: PETSC_DECIDE, PETSC_IGNORE, PETSC_DETERMINE

289: M*/
290: #define PETSC_DEFAULT -2

292: /*MC
293:     PETSC_COMM_WORLD - the equivalent of the MPI_COMM_WORLD communicator which represents
294:            all the processes that PETSc knows about.

296:    Level: beginner

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

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

307: .seealso: PETSC_COMM_SELF

309: M*/
310: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;

312: /*MC
313:     PETSC_COMM_SELF - This is always MPI_COMM_SELF

315:    Level: beginner

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

320: .seealso: PETSC_COMM_WORLD

322: M*/
323: #define PETSC_COMM_SELF MPI_COMM_SELF

325: /*MC
326:     PETSC_MPI_THREAD_REQUIRED - the required threading support used if PETSc initializes
327:            MPI with MPI_Init_thread.

329:    Level: beginner

331:    Notes:
332:    By default PETSC_MPI_THREAD_REQUIRED equals MPI_THREAD_FUNNELED.

334: .seealso: PetscInitialize()

336: M*/
337: PETSC_EXTERN PetscMPIInt PETSC_MPI_THREAD_REQUIRED;

339: PETSC_EXTERN PetscBool PetscBeganMPI;
340: PETSC_EXTERN PetscBool PetscErrorHandlingInitialized;
341: PETSC_EXTERN PetscBool PetscInitializeCalled;
342: PETSC_EXTERN PetscBool PetscFinalizeCalled;
343: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;

345: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm));
346: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*);
347: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm*);
348: PETSC_EXTERN PetscErrorCode PetscCommGetComm(MPI_Comm,MPI_Comm*);
349: PETSC_EXTERN PetscErrorCode PetscCommRestoreComm(MPI_Comm,MPI_Comm*);

351: #if defined(PETSC_HAVE_KOKKOS)
352: PETSC_EXTERN PetscErrorCode PetscKokkosInitializeCheck(void);  /* Initialize Kokkos if not yet. */
353: #endif

355: #if defined(PETSC_HAVE_NVSHMEM)
356: PETSC_EXTERN PetscBool      PetscBeganNvshmem;
357: PETSC_EXTERN PetscBool      PetscNvshmemInitialized;
358: PETSC_EXTERN PetscErrorCode PetscNvshmemFinalize(void);
359: #endif

361: #if defined(PETSC_HAVE_ELEMENTAL)
362: PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
363: PETSC_EXTERN PetscErrorCode PetscElementalInitialized(PetscBool*);
364: PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
365: #endif

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

370:    Synopsis:
371: #include <petscsys.h>
372:    PetscErrorCode PetscMalloc(size_t m,void **result)

374:    Not Collective

376:    Input Parameter:
377: .  m - number of bytes to allocate

379:    Output Parameter:
380: .  result - memory allocated

382:    Level: beginner

384:    Notes:
385:    Memory is always allocated at least double aligned

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

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

391: M*/
392: #define PetscMalloc(a,b)  ((*PetscTrMalloc)((a),PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(void**)(b)))

394: /*MC
395:    PetscRealloc - Rellocates memory

397:    Synopsis:
398: #include <petscsys.h>
399:    PetscErrorCode PetscRealloc(size_t m,void **result)

401:    Not Collective

403:    Input Parameters:
404: +  m - number of bytes to allocate
405: -  result - previous memory

407:    Output Parameter:
408: .  result - new memory allocated

410:    Level: developer

412:    Notes:
413:    Memory is always allocated at least double aligned

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

417: M*/
418: #define PetscRealloc(a,b)  ((*PetscTrRealloc)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__,(void**)(b)))

420: /*MC
421:    PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment

423:    Synopsis:
424: #include <petscsys.h>
425:    void *PetscAddrAlign(void *addr)

427:    Not Collective

429:    Input Parameters:
430: .  addr - address to align (any pointer type)

432:    Level: developer

434: .seealso: PetscMallocAlign()

436: M*/
437: #define PetscAddrAlign(a) (void*)((((PETSC_UINTPTR_T)(a))+(PETSC_MEMALIGN-1)) & ~(PETSC_MEMALIGN-1))

439: /*MC
440:    PetscCalloc - Allocates a cleared (zeroed) memory region aligned to PETSC_MEMALIGN

442:    Synopsis:
443: #include <petscsys.h>
444:    PetscErrorCode PetscCalloc(size_t m,void **result)

446:    Not Collective

448:    Input Parameter:
449: .  m - number of bytes to allocate

451:    Output Parameter:
452: .  result - memory allocated

454:    Level: beginner

456:    Notes:
457:    Memory is always allocated at least double aligned. This macro is useful in allocating memory pointed by void pointers

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

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

463: M*/
464: #define PetscCalloc(m,result)  PetscMallocA(1,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)((size_t)m),(result))

466: /*MC
467:    PetscMalloc1 - Allocates an array of memory aligned to PETSC_MEMALIGN

469:    Synopsis:
470: #include <petscsys.h>
471:    PetscErrorCode PetscMalloc1(size_t m1,type **r1)

473:    Not Collective

475:    Input Parameter:
476: .  m1 - number of elements to allocate  (may be zero)

478:    Output Parameter:
479: .  r1 - memory allocated

481:    Note:
482:    This uses the sizeof() of the memory type requested to determine the total memory to be allocated, therefore you should not
483:          multiply the number of elements requested by the sizeof() the type. For example use
484: $  PetscInt *id;
485: $  PetscMalloc1(10,&id);
486:         not
487: $  PetscInt *id;
488: $  PetscMalloc1(10*sizeof(PetscInt),&id);

490:         Does not zero the memory allocated, use PetscCalloc1() to obtain memory that has been zeroed.

492:    Level: beginner

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

496: M*/
497: #define PetscMalloc1(m1,r1) PetscMallocA(1,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)((size_t)m1)*sizeof(**(r1)),(r1))

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

502:    Synopsis:
503: #include <petscsys.h>
504:    PetscErrorCode PetscCalloc1(size_t m1,type **r1)

506:    Not Collective

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

511:    Output Parameter:
512: .  r1 - memory allocated

514:    Notes:
515:    See PetsMalloc1() for more details on usage.

517:    Level: beginner

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

521: M*/
522: #define PetscCalloc1(m1,r1) PetscMallocA(1,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)((size_t)m1)*sizeof(**(r1)),(r1))

524: /*MC
525:    PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN

527:    Synopsis:
528: #include <petscsys.h>
529:    PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)

531:    Not Collective

533:    Input Parameters:
534: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
535: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

537:    Output Parameters:
538: +  r1 - memory allocated in first chunk
539: -  r2 - memory allocated in second chunk

541:    Level: developer

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

545: M*/
546: #define PetscMalloc2(m1,r1,m2,r2) PetscMallocA(2,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)((size_t)m1)*sizeof(**(r1)),(r1),(size_t)((size_t)m2)*sizeof(**(r2)),(r2))

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

551:    Synopsis:
552: #include <petscsys.h>
553:    PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)

555:    Not Collective

557:    Input Parameters:
558: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
559: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

561:    Output Parameters:
562: +  r1 - memory allocated in first chunk
563: -  r2 - memory allocated in second chunk

565:    Level: developer

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

569: M*/
570: #define PetscCalloc2(m1,r1,m2,r2) PetscMallocA(2,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)((size_t)m1)*sizeof(**(r1)),(r1),(size_t)((size_t)m2)*sizeof(**(r2)),(r2))

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

575:    Synopsis:
576: #include <petscsys.h>
577:    PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

579:    Not Collective

581:    Input Parameters:
582: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
583: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
584: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

586:    Output Parameters:
587: +  r1 - memory allocated in first chunk
588: .  r2 - memory allocated in second chunk
589: -  r3 - memory allocated in third chunk

591:    Level: developer

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

595: M*/
596: #define PetscMalloc3(m1,r1,m2,r2,m3,r3) PetscMallocA(3,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)((size_t)m1)*sizeof(**(r1)),(r1),(size_t)((size_t)m2)*sizeof(**(r2)),(r2),(size_t)((size_t)m3)*sizeof(**(r3)),(r3))

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

601:    Synopsis:
602: #include <petscsys.h>
603:    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

605:    Not Collective

607:    Input Parameters:
608: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
609: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
610: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

612:    Output Parameters:
613: +  r1 - memory allocated in first chunk
614: .  r2 - memory allocated in second chunk
615: -  r3 - memory allocated in third chunk

617:    Level: developer

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

621: M*/
622: #define PetscCalloc3(m1,r1,m2,r2,m3,r3) PetscMallocA(3,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)((size_t)m1)*sizeof(**(r1)),(r1),(size_t)((size_t)m2)*sizeof(**(r2)),(r2),(size_t)((size_t)m3)*sizeof(**(r3)),(r3))

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

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

631:    Not Collective

633:    Input Parameters:
634: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
635: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
636: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
637: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

639:    Output Parameters:
640: +  r1 - memory allocated in first chunk
641: .  r2 - memory allocated in second chunk
642: .  r3 - memory allocated in third chunk
643: -  r4 - memory allocated in fourth chunk

645:    Level: developer

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

649: M*/
650: #define PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4) PetscMallocA(4,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)((size_t)m1)*sizeof(**(r1)),(r1),(size_t)((size_t)m2)*sizeof(**(r2)),(r2),(size_t)((size_t)m3)*sizeof(**(r3)),(r3),(size_t)((size_t)m4)*sizeof(**(r4)),(r4))

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

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

659:    Not Collective

661:    Input Parameters:
662: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
663: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
664: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
665: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

667:    Output Parameters:
668: +  r1 - memory allocated in first chunk
669: .  r2 - memory allocated in second chunk
670: .  r3 - memory allocated in third chunk
671: -  r4 - memory allocated in fourth chunk

673:    Level: developer

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

677: M*/
678: #define PetscCalloc4(m1,r1,m2,r2,m3,r3,m4,r4) PetscMallocA(4,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)((size_t)m1)*sizeof(**(r1)),(r1),(size_t)((size_t)m2)*sizeof(**(r2)),(r2),(size_t)((size_t)m3)*sizeof(**(r3)),(r3),(size_t)((size_t)m4)*sizeof(**(r4)),(r4))

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

683:    Synopsis:
684: #include <petscsys.h>
685:    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)

687:    Not Collective

689:    Input Parameters:
690: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
691: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
692: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
693: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
694: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

696:    Output Parameters:
697: +  r1 - memory allocated in first chunk
698: .  r2 - memory allocated in second chunk
699: .  r3 - memory allocated in third chunk
700: .  r4 - memory allocated in fourth chunk
701: -  r5 - memory allocated in fifth chunk

703:    Level: developer

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

707: M*/
708: #define PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5) PetscMallocA(5,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)((size_t)m1)*sizeof(**(r1)),(r1),(size_t)((size_t)m2)*sizeof(**(r2)),(r2),(size_t)((size_t)m3)*sizeof(**(r3)),(r3),(size_t)((size_t)m4)*sizeof(**(r4)),(r4),(size_t)((size_t)m5)*sizeof(**(r5)),(r5))

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

713:    Synopsis:
714: #include <petscsys.h>
715:    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)

717:    Not Collective

719:    Input Parameters:
720: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
721: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
722: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
723: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
724: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

726:    Output Parameters:
727: +  r1 - memory allocated in first chunk
728: .  r2 - memory allocated in second chunk
729: .  r3 - memory allocated in third chunk
730: .  r4 - memory allocated in fourth chunk
731: -  r5 - memory allocated in fifth chunk

733:    Level: developer

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

737: M*/
738: #define PetscCalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5) PetscMallocA(5,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)((size_t)m1)*sizeof(**(r1)),(r1),(size_t)((size_t)m2)*sizeof(**(r2)),(r2),(size_t)((size_t)m3)*sizeof(**(r3)),(r3),(size_t)((size_t)m4)*sizeof(**(r4)),(r4),(size_t)((size_t)m5)*sizeof(**(r5)),(r5))

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

743:    Synopsis:
744: #include <petscsys.h>
745:    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)

747:    Not Collective

749:    Input Parameters:
750: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
751: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
752: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
753: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
754: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
755: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

757:    Output Parameteasr:
758: +  r1 - memory allocated in first chunk
759: .  r2 - memory allocated in second chunk
760: .  r3 - memory allocated in third chunk
761: .  r4 - memory allocated in fourth chunk
762: .  r5 - memory allocated in fifth chunk
763: -  r6 - memory allocated in sixth chunk

765:    Level: developer

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

769: M*/
770: #define PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) PetscMallocA(6,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)((size_t)m1)*sizeof(**(r1)),(r1),(size_t)((size_t)m2)*sizeof(**(r2)),(r2),(size_t)((size_t)m3)*sizeof(**(r3)),(r3),(size_t)((size_t)m4)*sizeof(**(r4)),(r4),(size_t)((size_t)m5)*sizeof(**(r5)),(r5),(size_t)((size_t)m6)*sizeof(**(r6)),(r6))

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

775:    Synopsis:
776: #include <petscsys.h>
777:    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)

779:    Not Collective

781:    Input Parameters:
782: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
783: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
784: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
785: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
786: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
787: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

789:    Output Parameters:
790: +  r1 - memory allocated in first chunk
791: .  r2 - memory allocated in second chunk
792: .  r3 - memory allocated in third chunk
793: .  r4 - memory allocated in fourth chunk
794: .  r5 - memory allocated in fifth chunk
795: -  r6 - memory allocated in sixth chunk

797:    Level: developer

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

801: M*/
802: #define PetscCalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) PetscMallocA(6,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)((size_t)m1)*sizeof(**(r1)),(r1),(size_t)((size_t)m2)*sizeof(**(r2)),(r2),(size_t)((size_t)m3)*sizeof(**(r3)),(r3),(size_t)((size_t)m4)*sizeof(**(r4)),(r4),(size_t)((size_t)m5)*sizeof(**(r5)),(r5),(size_t)((size_t)m6)*sizeof(**(r6)),(r6))

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

807:    Synopsis:
808: #include <petscsys.h>
809:    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)

811:    Not Collective

813:    Input Parameters:
814: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
815: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
816: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
817: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
818: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
819: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
820: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

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

831:    Level: developer

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

835: M*/
836: #define PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) PetscMallocA(7,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)((size_t)m1)*sizeof(**(r1)),(r1),(size_t)((size_t)m2)*sizeof(**(r2)),(r2),(size_t)((size_t)m3)*sizeof(**(r3)),(r3),(size_t)((size_t)m4)*sizeof(**(r4)),(r4),(size_t)((size_t)m5)*sizeof(**(r5)),(r5),(size_t)((size_t)m6)*sizeof(**(r6)),(r6),(size_t)((size_t)m7)*sizeof(**(r7)),(r7))

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

841:    Synopsis:
842: #include <petscsys.h>
843:    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)

845:    Not Collective

847:    Input Parameters:
848: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
849: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
850: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
851: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
852: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
853: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
854: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

856:    Output Parameters:
857: +  r1 - memory allocated in first chunk
858: .  r2 - memory allocated in second chunk
859: .  r3 - memory allocated in third chunk
860: .  r4 - memory allocated in fourth chunk
861: .  r5 - memory allocated in fifth chunk
862: .  r6 - memory allocated in sixth chunk
863: -  r7 - memory allocated in seventh chunk

865:    Level: developer

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

869: M*/
870: #define PetscCalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) PetscMallocA(7,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)((size_t)m1)*sizeof(**(r1)),(r1),(size_t)((size_t)m2)*sizeof(**(r2)),(r2),(size_t)((size_t)m3)*sizeof(**(r3)),(r3),(size_t)((size_t)m4)*sizeof(**(r4)),(r4),(size_t)((size_t)m5)*sizeof(**(r5)),(r5),(size_t)((size_t)m6)*sizeof(**(r6)),(r6),(size_t)((size_t)m7)*sizeof(**(r7)),(r7))

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

875:    Synopsis:
876: #include <petscsys.h>
877:    PetscErrorCode PetscNew(type **result)

879:    Not Collective

881:    Output Parameter:
882: .  result - memory allocated, sized to match pointer type

884:    Level: beginner

886: .seealso: PetscFree(), PetscMalloc(), PetscNewLog(), PetscCalloc1(), PetscMalloc1()

888: M*/
889: #define PetscNew(b)      PetscCalloc1(1,(b))

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

895:    Synopsis:
896: #include <petscsys.h>
897:    PetscErrorCode PetscNewLog(PetscObject obj,type **result)

899:    Not Collective

901:    Input Parameter:
902: .  obj - object memory is logged to

904:    Output Parameter:
905: .  result - memory allocated, sized to match pointer type

907:    Level: developer

909: .seealso: PetscFree(), PetscMalloc(), PetscNew(), PetscLogObjectMemory(), PetscCalloc1(), PetscMalloc1()

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

914: /*MC
915:    PetscFree - Frees memory

917:    Synopsis:
918: #include <petscsys.h>
919:    PetscErrorCode PetscFree(void *memory)

921:    Not Collective

923:    Input Parameter:
924: .   memory - memory to free (the pointer is ALWAYS set to NULL upon success)

926:    Level: beginner

928:    Notes:
929:    Do not free memory obtained with PetscMalloc2(), PetscCalloc2() etc, they must be freed with PetscFree2() etc.

931:    It is safe to call PetscFree() on a NULL pointer.

933: .seealso: PetscNew(), PetscMalloc(), PetscNewLog(), PetscMalloc1(), PetscCalloc1()

935: M*/
936: #define PetscFree(a)   ((*PetscTrFree)((void*)(a),__LINE__,PETSC_FUNCTION_NAME,__FILE__) || ((a) = NULL,0))

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

941:    Synopsis:
942: #include <petscsys.h>
943:    PetscErrorCode PetscFree2(void *memory1,void *memory2)

945:    Not Collective

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

951:    Level: developer

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

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

958: M*/
959: #define PetscFree2(m1,m2)   PetscFreeA(2,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2))

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

964:    Synopsis:
965: #include <petscsys.h>
966:    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)

968:    Not Collective

970:    Input Parameters:
971: +   memory1 - memory to free
972: .   memory2 - 2nd memory to free
973: -   memory3 - 3rd memory to free

975:    Level: developer

977:    Notes:
978:     Memory must have been obtained with PetscMalloc3()

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

982: M*/
983: #define PetscFree3(m1,m2,m3)   PetscFreeA(3,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3))

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

988:    Synopsis:
989: #include <petscsys.h>
990:    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)

992:    Not Collective

994:    Input Parameters:
995: +   m1 - memory to free
996: .   m2 - 2nd memory to free
997: .   m3 - 3rd memory to free
998: -   m4 - 4th memory to free

1000:    Level: developer

1002:    Notes:
1003:     Memory must have been obtained with PetscMalloc4()

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

1007: M*/
1008: #define PetscFree4(m1,m2,m3,m4)   PetscFreeA(4,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4))

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

1013:    Synopsis:
1014: #include <petscsys.h>
1015:    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)

1017:    Not Collective

1019:    Input Parameters:
1020: +   m1 - memory to free
1021: .   m2 - 2nd memory to free
1022: .   m3 - 3rd memory to free
1023: .   m4 - 4th memory to free
1024: -   m5 - 5th memory to free

1026:    Level: developer

1028:    Notes:
1029:     Memory must have been obtained with PetscMalloc5()

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

1033: M*/
1034: #define PetscFree5(m1,m2,m3,m4,m5)   PetscFreeA(5,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5))

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

1039:    Synopsis:
1040: #include <petscsys.h>
1041:    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)

1043:    Not Collective

1045:    Input Parameters:
1046: +   m1 - memory to free
1047: .   m2 - 2nd memory to free
1048: .   m3 - 3rd memory to free
1049: .   m4 - 4th memory to free
1050: .   m5 - 5th memory to free
1051: -   m6 - 6th memory to free

1053:    Level: developer

1055:    Notes:
1056:     Memory must have been obtained with PetscMalloc6()

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

1060: M*/
1061: #define PetscFree6(m1,m2,m3,m4,m5,m6)   PetscFreeA(6,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5),&(m6))

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

1066:    Synopsis:
1067: #include <petscsys.h>
1068:    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)

1070:    Not Collective

1072:    Input Parameters:
1073: +   m1 - memory to free
1074: .   m2 - 2nd memory to free
1075: .   m3 - 3rd memory to free
1076: .   m4 - 4th memory to free
1077: .   m5 - 5th memory to free
1078: .   m6 - 6th memory to free
1079: -   m7 - 7th memory to free

1081:    Level: developer

1083:    Notes:
1084:     Memory must have been obtained with PetscMalloc7()

1086: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1087:           PetscMalloc7()

1089: M*/
1090: #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   PetscFreeA(7,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5),&(m6),&(m7))

1092: PETSC_EXTERN PetscErrorCode PetscMallocA(int,PetscBool,int,const char *,const char *,size_t,void *,...);
1093: PETSC_EXTERN PetscErrorCode PetscFreeA(int,int,const char *,const char *,void *,...);
1094: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,PetscBool,int,const char[],const char[],void**);
1095: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1096: PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t,int,const char[],const char[],void**);
1097: PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1098: PETSC_EXTERN PetscErrorCode PetscMallocSet(PetscErrorCode (*)(size_t,PetscBool,int,const char[],const char[],void**),PetscErrorCode (*)(void*,int,const char[],const char[]),PetscErrorCode (*)(size_t,int,const char[],const char[], void **));
1099: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);

1101: /*
1102:   Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1103: */
1104: PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1105: PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);
1106: #if defined(PETSC_HAVE_CUDA)
1107: PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void);
1108: PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void);
1109: #endif
1110: #if defined(PETSC_HAVE_HIP)
1111: PETSC_EXTERN PetscErrorCode PetscMallocSetHIPHost(void);
1112: PETSC_EXTERN PetscErrorCode PetscMallocResetHIPHost(void);
1113: #endif

1115: #define MPIU_PETSCLOGDOUBLE  MPI_DOUBLE
1116: #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION

1118: /*
1119:    Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1120: */
1121: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1122: PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1123: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1124: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1125: PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1126: PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int,PetscLogDouble*);
1127: PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool,PetscBool);
1128: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*,PetscBool*,PetscBool*);
1129: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1130: PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1131: PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool*);
1132: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool);
1133: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool*);

1135: PETSC_EXTERN const char *const PetscDataTypes[];
1136: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1137: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1138: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1139: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);

1141: /*
1142:     Basic memory and string operations. These are usually simple wrappers
1143:    around the basic Unix system calls, but a few of them have additional
1144:    functionality and/or error checking.
1145: */
1146: PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1147: PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1148: PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***);
1149: PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1150: PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool  *);
1151: PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool  *);
1152: PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1153: PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1154: PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1155: PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1156: PETSC_EXTERN PetscErrorCode PetscStrlcat(char[],const char[],size_t);
1157: PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1158: PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1159: PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1160: PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1161: PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1162: PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1163: PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1164: PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1165: PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*);
1166: PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1167: PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1168: PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1169: PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1170: PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt,const char *const*,char***);
1171: PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt,char***);
1172: PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);

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

1176: PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1177: PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1178: PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);

1180: PETSC_EXTERN PetscErrorCode PetscStrInList(const char[],const char[],char,PetscBool*);
1181: PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1182: PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);

1184: /*
1185:    These are MPI operations for MPI_Allreduce() etc
1186: */
1187: PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1188: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1189: PETSC_EXTERN MPI_Op MPIU_SUM;
1190: PETSC_EXTERN MPI_Op MPIU_MAX;
1191: PETSC_EXTERN MPI_Op MPIU_MIN;
1192: #else
1193: #define MPIU_SUM MPI_SUM
1194: #define MPIU_MAX MPI_MAX
1195: #define MPIU_MIN MPI_MIN
1196: #endif
1197: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);

1199: PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1200: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);

1202: PETSC_EXTERN const char *const PetscFileModes[];

1204: /*
1205:     Defines PETSc error handling.
1206: */
1207: #include <petscerror.h>

1209: #define PETSC_SMALLEST_CLASSID  1211211
1210: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1211: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1212: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);
1213: PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject,PetscObjectId*);
1214: PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject,PetscObjectId,PetscBool*);

1216: /*
1217:    Routines that get memory usage information from the OS
1218: */
1219: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1220: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1221: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1222: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);

1224: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);

1226: /*
1227:    Initialization of PETSc
1228: */
1229: PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1230: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1231: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1232: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1233: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1234: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1235: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1236: PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1237: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1238: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);

1240: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1241: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);

1243: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1244: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1245: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1246: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);

1248: PETSC_EXTERN PetscErrorCode PetscMonitorCompare(PetscErrorCode (*)(void),void *,PetscErrorCode (*)(void**),PetscErrorCode (*)(void),void *,PetscErrorCode (*)(void**),PetscBool *);

1250: /*
1251:      These are so that in extern C code we can caste function pointers to non-extern C
1252:    function pointers. Since the regular C++ code expects its function pointers to be C++
1253: */
1254: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1255: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1256: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);

1258: /*
1259:     Functions that can act on any PETSc object.
1260: */
1261: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1262: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1263: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1264: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1265: PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1266: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1267: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1268: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1269: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1270: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1271: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1272: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1273: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1274: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1275: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt*);
1276: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1277: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1278: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject*);
1279: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void));
1280: #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d))
1281: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1282: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1283: PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1284: PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject,PetscObject);
1285: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt*);

1287: #include <petscviewertypes.h>
1288: #include <petscoptions.h>

1290: PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer,PetscBool,PetscLogDouble);
1291: PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool*);

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

1295: PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer,const char[]);
1296: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1297: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1298: #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1299: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1300: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1301: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1302: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1303: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1304: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1305: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1306: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1307: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,PetscObject,const char[]);
1308: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1309: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1310: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject,const char[],PetscBool *);
1311: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1312: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1313: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1314: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);

1316: #if defined(PETSC_HAVE_SAWS)
1317: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1318: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1319: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1320: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1321: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1322: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1323: PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1324: PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1325: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1326: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);

1328: #else
1329: #define PetscSAWsBlock()                        0
1330: #define PetscObjectSAWsViewOff(obj)             0
1331: #define PetscObjectSAWsSetBlock(obj,flg)        0
1332: #define PetscObjectSAWsBlock(obj)               0
1333: #define PetscObjectSAWsGrantAccess(obj)         0
1334: #define PetscObjectSAWsTakeAccess(obj)          0
1335: #define PetscStackViewSAWs()                    0
1336: #define PetscStackSAWsViewOff()                 0
1337: #define PetscStackSAWsTakeAccess()
1338: #define PetscStackSAWsGrantAccess()

1340: #endif

1342: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1343: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1344: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);
1345: PETSC_EXTERN PetscErrorCode PetscDLAddr(void (*)(void), char **);
1346: #ifdef PETSC_HAVE_CXX
1347: PETSC_EXTERN PetscErrorCode PetscDemangleSymbol(const char *, char **);
1348: #endif

1350: #if defined(PETSC_USE_DEBUG)
1351: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1352: #endif
1353: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);

1355: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1356: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1357: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1358: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1359: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1360: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);

1362: /*
1363:     Dynamic library lists. Lists of names of routines in objects or in dynamic
1364:   link libraries that will be loaded as needed.
1365: */

1367: #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1368: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1369: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1370: #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1371: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1372: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[],const char[]);
1373: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1374: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1375: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);

1377: PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1378: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1379: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1380: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1381: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1382: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1383: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1384: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);

1386: /*
1387:      Useful utility routines
1388: */
1389: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1390: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1391: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm,PetscInt*,PetscInt*);
1392: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1393: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1394: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1395: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);
1396: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm,const PetscInt[2],PetscInt[2]);
1397: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm,const PetscReal[2],PetscReal[2]);

1399: /*MC
1400:     PetscNot - negates a logical type value and returns result as a PetscBool

1402:     Notes:
1403:     This is useful in cases like
1404: $     int        *a;
1405: $     PetscBool  flag = PetscNot(a)
1406:      where !a would not return a PetscBool because we cannot provide a cast from int to PetscBool in C.

1408:     Level: beginner

1410:     .seealso : PetscBool, PETSC_TRUE, PETSC_FALSE
1411: M*/
1412: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)

1414: /*MC
1415:     PetscHelpPrintf - Prints help messages.

1417:    Synopsis:
1418: #include <petscsys.h>
1419:      PetscErrorCode (*PetscHelpPrintf)(MPI_Comm comm, const char format[],args);

1421:     Collective on comm

1423:     Input Parameters:
1424: +  comm - the MPI communicator over which the help message is printed
1425: .  format - the usual printf() format string
1426: -  args - arguments to be printed

1428:    Level: developer

1430:    Fortran Note:
1431:      This routine is not supported in Fortran.

1433:    Note:
1434:      You can change how help messages are printed by replacing the function pointer with a function that does not simply write to stdout.

1436:       To use, write your own function, for example,
1437: $PetscErrorCode mypetschelpprintf(MPI_Comm comm,const char format[],....)
1438: ${
1439: $ return 0;
1440: $}
1441: then do the assigment
1442: $    PetscHelpPrintf = mypetschelpprintf;
1443:    You can do the assignment before PetscInitialize().

1445:   The default routine used is called PetscHelpPrintfDefault().

1447: .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1448: M*/
1449: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...) PETSC_ATTRIBUTE_FORMAT(2,3);

1451: /*
1452:      Defines PETSc profiling.
1453: */
1454: #include <petsclog.h>

1456: /*
1457:       Simple PETSc parallel IO for ASCII printing
1458: */
1459: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1460: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1461: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1462: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...) PETSC_ATTRIBUTE_FORMAT(3,4);
1463: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...) PETSC_ATTRIBUTE_FORMAT(2,3);
1464: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...) PETSC_ATTRIBUTE_FORMAT(3,4);
1465: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...) PETSC_ATTRIBUTE_FORMAT(3,5);
1466: PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[],size_t,const char*,PetscInt,const PetscReal[]);

1468: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...) PETSC_ATTRIBUTE_FORMAT(1,2);
1469: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...) PETSC_ATTRIBUTE_FORMAT(1,2);
1470: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...) PETSC_ATTRIBUTE_FORMAT(2,3);

1472: PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char*,size_t*);
1473: PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char*,char *);

1475: #if defined(PETSC_HAVE_POPEN)
1476: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1477: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*);
1478: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1479: #endif

1481: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...) PETSC_ATTRIBUTE_FORMAT(2,3);
1482: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...) PETSC_ATTRIBUTE_FORMAT(3,4);
1483: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*);
1484: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1485: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1486: PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1487: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);

1489: PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1490: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void**);
1491: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void*);
1492: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1493: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer*);
1494: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer,PetscErrorCode (*)(void*));
1495: PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void*);

1497: /*
1498:    For use in debuggers
1499: */
1500: PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1501: PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1502: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1503: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1504: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);

1506: #include <stddef.h>
1507: #include <string.h>             /* for memcpy, memset */
1508: #include <stdlib.h>

1510: #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1511: #include <xmmintrin.h>
1512: #endif

1514: /*@C
1515:    PetscMemmove - Copies n bytes, beginning at location b, to the space
1516:    beginning at location a. Copying  between regions that overlap will
1517:    take place correctly. Use PetscMemcpy() if the locations do not overlap

1519:    Not Collective

1521:    Input Parameters:
1522: +  b - pointer to initial memory space
1523: .  a - pointer to copy space
1524: -  n - length (in bytes) of space to copy

1526:    Level: intermediate

1528:    Note:
1529:    PetscArraymove() is preferred
1530:    This routine is analogous to memmove().

1532:    Developers Note: This is inlined for performance

1534: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscStrallocpy(),
1535:           PetscArraymove()
1536: @*/
1537: static inline PetscErrorCode PetscMemmove(void *a,const void *b,size_t n)
1538: {
1539:   if (n > 0) {
1542:   }
1543: #if !defined(PETSC_HAVE_MEMMOVE)
1544:   if (a < b) {
1545:     if (a <= b - n) memcpy(a,b,n);
1546:     else {
1547:       memcpy(a,b,(int)(b - a));
1548:       PetscMemmove(b,b + (int)(b - a),n - (int)(b - a));
1549:     }
1550:   } else {
1551:     if (b <= a - n) memcpy(a,b,n);
1552:     else {
1553:       memcpy(b + n,b + (n - (int)(a - b)),(int)(a - b));
1554:       PetscMemmove(a,b,n - (int)(a - b));
1555:     }
1556:   }
1557: #else
1558:   memmove((char*)(a),(char*)(b),n);
1559: #endif
1560:   return 0;
1561: }

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

1568:    Not Collective

1570:    Input Parameters:
1571: +  b - pointer to initial memory space
1572: -  n - length (in bytes) of space to copy

1574:    Output Parameter:
1575: .  a - pointer to copy space

1577:    Level: intermediate

1579:    Compile Option:
1580:     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1581:                                   for memory copies on double precision values.
1582:     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1583:                                   for memory copies on double precision values.
1584:     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1585:                                   for memory copies on double precision values.

1587:    Note:
1588:    Prefer PetscArraycpy()
1589:    This routine is analogous to memcpy().
1590:    Not available from Fortran

1592:    Developer Note: this is inlined for fastest performance

1594: .seealso: PetscMemzero(), PetscMemcmp(), PetscArrayzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy()

1596: @*/
1597: static inline PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n)
1598: {
1599: #if defined(PETSC_USE_DEBUG)
1600:   size_t al = (size_t) a,bl = (size_t) b;
1601:   size_t nl = (size_t) n;
1602:   if (n > 0) {
1605:   }
1606: #else
1607: #endif
1608:   if (a != b && n > 0) {
1609: #if defined(PETSC_USE_DEBUG)
1611:               or make sure your copy regions and lengths are correct. \n\
1612:               Length (bytes) %zu first address %zu second address %zu",nl,al,bl);
1613: #endif
1614: #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1615:    if (!(a % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1616:       size_t len = n/sizeof(PetscScalar);
1617: #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1618:       PetscBLASInt one = 1,blen;
1619:       PetscBLASIntCast(len,&blen);
1620:       PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one));
1621: #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1622:       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1623: #else
1624:       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1625:       for (size_t i=0; i<len; i++) y[i] = x[i];
1626: #endif
1627:     } else {
1628:       memcpy((char*)(a),(char*)(b),n);
1629:     }
1630: #else
1631:     memcpy((char*)(a),(char*)(b),n);
1632: #endif
1633:   }
1634:   return 0;
1635: }

1637: /*@C
1638:    PetscMemzero - Zeros the specified memory.

1640:    Not Collective

1642:    Input Parameters:
1643: +  a - pointer to beginning memory location
1644: -  n - length (in bytes) of memory to initialize

1646:    Level: intermediate

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

1652:    Not available from Fortran
1653:    Prefer PetscArrayzero()

1655:    Developer Note: this is inlined for fastest performance

1657: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy()
1658: @*/
1659: static inline PetscErrorCode PetscMemzero(void *a,size_t n)
1660: {
1661:   if (n > 0) {
1662: #if defined(PETSC_USE_DEBUG)
1664: #endif
1665: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1666:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1667:       size_t      i,len = n/sizeof(PetscScalar);
1668:       PetscScalar *x = (PetscScalar*)a;
1669:       for (i=0; i<len; i++) x[i] = 0.0;
1670:     } else {
1671: #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1672:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1673:       PetscInt len = n/sizeof(PetscScalar);
1674:       fortranzero_(&len,(PetscScalar*)a);
1675:     } else {
1676: #endif
1677: #if defined(PETSC_PREFER_BZERO)
1678:       bzero((char *)a,n);
1679: #else
1680:       memset((char*)a,0,n);
1681: #endif
1682: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1683:     }
1684: #endif
1685:   }
1686:   return 0;
1687: }

1689: /*MC
1690:    PetscArraycmp - Compares two arrays in memory.

1692:    Synopsis:
1693: #include <petscsys.h>
1694:     PetscErrorCode PetscArraycmp(const anytype *str1,const anytype *str2,size_t cnt,PetscBool *e)

1696:    Not Collective

1698:    Input Parameters:
1699: +  str1 - First array
1700: .  str2 - Second array
1701: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1703:    Output Parameters:
1704: .   e - PETSC_TRUE if equal else PETSC_FALSE.

1706:    Level: intermediate

1708:    Note:
1709:    This routine is a preferred replacement to PetscMemcmp()
1710:    The arrays must be of the same type

1712: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy(),
1713:           PetscArraymove()
1714: M*/
1715: #define  PetscArraycmp(str1,str2,cnt,e) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcmp(str1,str2,(size_t)(cnt)*sizeof(*(str1)),e))

1717: /*MC
1718:    PetscArraymove - Copies from one array in memory to another, the arrays may overlap. Use PetscArraycpy() when the arrays
1719:                     do not overlap

1721:    Synopsis:
1722: #include <petscsys.h>
1723:     PetscErrorCode PetscArraymove(anytype *str1,const anytype *str2,size_t cnt)

1725:    Not Collective

1727:    Input Parameters:
1728: +  str1 - First array
1729: .  str2 - Second array
1730: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1732:    Level: intermediate

1734:    Note:
1735:    This routine is a preferred replacement to PetscMemmove()
1736:    The arrays must be of the same type

1738: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycpy(), PetscMemmove(), PetscArraycmp(), PetscStrallocpy()
1739: M*/
1740: #define  PetscArraymove(str1,str2,cnt) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemmove(str1,str2,(size_t)(cnt)*sizeof(*(str1))))

1742: /*MC
1743:    PetscArraycpy - Copies from one array in memory to another

1745:    Synopsis:
1746: #include <petscsys.h>
1747:     PetscErrorCode PetscArraycpy(anytype *str1,const anytype *str2,size_t cnt)

1749:    Not Collective

1751:    Input Parameters:
1752: +  str1 - First array (destination)
1753: .  str2 - Second array (source)
1754: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1756:    Level: intermediate

1758:    Note:
1759:    This routine is a preferred replacement to PetscMemcpy()
1760:    The arrays must be of the same type

1762: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraymove(), PetscMemmove(), PetscArraycmp(), PetscStrallocpy()
1763: M*/
1764: #define  PetscArraycpy(str1,str2,cnt) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcpy(str1,str2,(size_t)(cnt)*sizeof(*(str1))))

1766: /*MC
1767:    PetscArrayzero - Zeros an array in memory.

1769:    Synopsis:
1770: #include <petscsys.h>
1771:     PetscErrorCode PetscArrayzero(anytype *str1,size_t cnt)

1773:    Not Collective

1775:    Input Parameters:
1776: +  str1 - array
1777: -  cnt  - Count of the array, not in bytes, but number of entries in the array

1779:    Level: intermediate

1781:    Note:
1782:    This routine is a preferred replacement to PetscMemzero()

1784: .seealso: PetscMemcpy(), PetscMemcmp(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy(), PetscArraymove()
1785: M*/
1786: #define  PetscArrayzero(str1,cnt) PetscMemzero(str1,(size_t)(cnt)*sizeof(*(str1)))

1788: /*MC
1789:    PetscPrefetchBlock - Prefetches a block of memory

1791:    Synopsis:
1792: #include <petscsys.h>
1793:     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)

1795:    Not Collective

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

1803:    Level: developer

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

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

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

1818: M*/
1819: #define PetscPrefetchBlock(a,n,rw,t) do {                               \
1820:     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
1821:     for (; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
1822:   } while (0)

1824: /*
1825:       Determine if some of the kernel computation routines use
1826:    Fortran (rather than C) for the numerical calculations. On some machines
1827:    and compilers (like complex numbers) the Fortran version of the routines
1828:    is faster than the C/C++ versions. The flag --with-fortran-kernels
1829:    should be used with ./configure to turn these on.
1830: */
1831: #if defined(PETSC_USE_FORTRAN_KERNELS)

1833: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1834: #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1835: #endif

1837: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
1838: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
1839: #endif

1841: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1842: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
1843: #endif

1845: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1846: #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
1847: #endif

1849: #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
1850: #define PETSC_USE_FORTRAN_KERNEL_NORM
1851: #endif

1853: #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
1854: #define PETSC_USE_FORTRAN_KERNEL_MAXPY
1855: #endif

1857: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1858: #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
1859: #endif

1861: #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1862: #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
1863: #endif

1865: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
1866: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
1867: #endif

1869: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1870: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
1871: #endif

1873: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
1874: #define PETSC_USE_FORTRAN_KERNEL_MDOT
1875: #endif

1877: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
1878: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
1879: #endif

1881: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
1882: #define PETSC_USE_FORTRAN_KERNEL_AYPX
1883: #endif

1885: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
1886: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
1887: #endif

1889: #endif

1891: /*
1892:     Macros for indicating code that should be compiled with a C interface,
1893:    rather than a C++ interface. Any routines that are dynamically loaded
1894:    (such as the PCCreate_XXX() routines) must be wrapped so that the name
1895:    mangler does not change the functions symbol name. This just hides the
1896:    ugly extern "C" {} wrappers.
1897: */
1898: #if defined(__cplusplus)
1899: #  define EXTERN_C_BEGIN extern "C" {
1900: #  define EXTERN_C_END }
1901: #else
1902: #  define EXTERN_C_BEGIN
1903: #  define EXTERN_C_END
1904: #endif

1906: /* --------------------------------------------------------------------*/

1908: /*MC
1909:     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
1910:         communication

1912:    Level: beginner

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

1916: .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
1917: M*/

1919: #if defined(PETSC_HAVE_MPIIO)
1920: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1921: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1922: PETSC_EXTERN PetscErrorCode MPIU_File_write_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1923: PETSC_EXTERN PetscErrorCode MPIU_File_read_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1924: PETSC_EXTERN PetscErrorCode MPIU_File_write_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1925: PETSC_EXTERN PetscErrorCode MPIU_File_read_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1926: #endif

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

1930: /* Limit MPI to 32-bits */
1931: #define PETSC_MPI_INT_MAX  2147483647
1932: #define PETSC_MPI_INT_MIN -2147483647
1933: /* Limit BLAS to 32-bits */
1934: #define PETSC_BLAS_INT_MAX  2147483647
1935: #define PETSC_BLAS_INT_MIN -2147483647
1936: #define PETSC_CUBLAS_INT_MAX  2147483647

1938: /*@C
1939:     PetscIntCast - casts a PetscInt64 (which is 64 bits in size) to a PetscInt (which may be 32 bits in size), generates an
1940:          error if the PetscInt is not large enough to hold the number.

1942:    Not Collective

1944:    Input Parameter:
1945: .     a - the PetscInt64 value

1947:    Output Parameter:
1948: .     b - the resulting PetscInt value

1950:    Level: advanced

1952:    Notes: If integers needed for the applications are too large to fit in 32 bit ints you can ./configure using --with-64-bit-indices to make PetscInt use 64 bit ints

1954:    Not available from Fortran

1956: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscBLASIntCast(), PetscIntMultError(), PetscIntSumError()
1957: @*/
1958: static inline PetscErrorCode PetscIntCast(PetscInt64 a,PetscInt *b)
1959: {
1960: #if !defined(PETSC_USE_64BIT_INDICES)
1961:   if (a > PETSC_MAX_INT) {
1962:     *b = 0;
1963:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"%" PetscInt64_FMT " is too big for PetscInt, you may need to ./configure using --with-64-bit-indices",a);
1964:   }
1965: #endif
1966:   *b = (PetscInt)(a);
1967:   return 0;
1968: }

1970: /*@C
1971:     PetscCountCast - casts a PetscCount to a PetscInt (which may be 32 bits in size), generates an
1972:          error if the PetscInt is not large enough to hold the number.

1974:    Not Collective

1976:    Input Parameter:
1977: .     a - the PetscCount value

1979:    Output Parameter:
1980: .     b - the resulting PetscInt value

1982:    Level: advanced

1984:    Notes:
1985:      If integers needed for the applications are too large to fit in 32 bit ints you can ./configure using --with-64-bit-indices to make PetscInt use 64 bit ints

1987:    Not available from Fortran

1989: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscBLASIntCast(), PetscIntMultError(), PetscIntSumError(), PetscIntCast()
1990: @*/
1991: static inline PetscErrorCode PetscCountCast(PetscCount a,PetscInt *b)
1992: {
1993:   if (sizeof(PetscCount) > sizeof(PetscInt) && a > PETSC_MAX_INT) {
1994:     *b = 0;
1995:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"%" PetscCount_FMT " is too big for PetscInt, you may need to ./configure using --with-64-bit-indices",a);
1996:   }
1997:   *b = (PetscInt)(a);
1998:   return 0;
1999: }

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

2005:    Not Collective

2007:    Input Parameter:
2008: .     a - the PetscInt value

2010:    Output Parameter:
2011: .     b - the resulting PetscBLASInt value

2013:    Level: advanced

2015:    Notes:
2016:       Not available from Fortran
2017:       Errors if the integer is negative since PETSc calls to BLAS/LAPACK never need to cast negative integer inputs

2019: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscIntCast(), PetscCountCast()
2020: @*/
2021: static inline PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
2022: {
2023:   *b = 0;
2024:   if (PetscDefined(USE_64BIT_INDICES) && !PetscDefined(HAVE_64BIT_BLAS_INDICES)) {
2026:   }
2028:   *b = (PetscBLASInt)(a);
2029:   return 0;
2030: }

2032: /*@C
2033:     PetscCuBLASIntCast - like PetscBLASIntCast(), but for PetscCuBLASInt.

2035:    Not Collective

2037:    Input Parameter:
2038: .     a - the PetscInt value

2040:    Output Parameter:
2041: .     b - the resulting PetscCuBLASInt value

2043:    Level: advanced

2045:    Notes:
2046:       Errors if the integer is negative since PETSc calls to cuBLAS and friends never need to cast negative integer inputs

2048: .seealso: PetscCuBLASInt, PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscMPIIntCast(), PetscIntCast()
2049: @*/
2050: static inline PetscErrorCode PetscCuBLASIntCast(PetscInt a,PetscCuBLASInt *b)
2051: {
2052:   *b = 0;
2055:   *b = (PetscCuBLASInt)(a);
2056:   return 0;
2057: }

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

2063:    Not Collective

2065:    Input Parameter:
2066: .     a - the PetscInt value

2068:    Output Parameter:
2069: .     b - the resulting PetscMPIInt value

2071:    Level: advanced

2073:    Not available from Fortran

2075: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntCast()
2076: @*/
2077: static inline PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
2078: {
2079:   *b = 0;
2081:   *b = (PetscMPIInt)(a);
2082:   return 0;
2083: }

2085: #define PetscInt64Mult(a,b)   ((PetscInt64)(a))*((PetscInt64)(b))

2087: /*@C

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

2091:    Not Collective

2093:    Input Parameters:
2094: +     a - the PetscReal value
2095: -     b - the second value

2097:    Returns:
2098:       the result as a PetscInt value

2100:    Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2101:    Use PetscIntMultTruncate() to compute the product of two positive PetscInt and truncate to fit a PetscInt
2102:    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

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

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

2109:    Not available from Fortran

2111:    Level: advanced

2113: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError(), PetscIntSumError()
2114: @*/
2115: static inline PetscInt PetscRealIntMultTruncate(PetscReal a,PetscInt b)
2116: {
2117:   PetscInt64 r = (PetscInt64)(a*(PetscReal)b);
2118:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2119:   return (PetscInt)r;
2120: }

2122: /*@C

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

2126:    Not Collective

2128:    Input Parameters:
2129: +     a - the PetscInt value
2130: -     b - the second value

2132:    Returns:
2133:       the result as a PetscInt value

2135:    Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2136:    Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2137:    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

2139:    Not available from Fortran

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

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

2145:    Level: advanced

2147: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError(), PetscIntSumError()
2148: @*/
2149: static inline PetscInt PetscIntMultTruncate(PetscInt a,PetscInt b)
2150: {
2151:   PetscInt64 r = PetscInt64Mult(a,b);
2152:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2153:   return (PetscInt)r;
2154: }

2156: /*@C

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

2160:    Not Collective

2162:    Input Parameters:
2163: +     a - the PetscInt value
2164: -     b - the second value

2166:    Returns:
2167:      the result as a PetscInt value

2169:    Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2170:    Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2171:    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

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

2175:    Not available from Fortran

2177:    Level: advanced

2179: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError()
2180: @*/
2181: static inline PetscInt PetscIntSumTruncate(PetscInt a,PetscInt b)
2182: {
2183:   PetscInt64 r = ((PetscInt64)a) + ((PetscInt64)b);
2184:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2185:   return (PetscInt)r;
2186: }

2188: /*@C

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

2192:    Not Collective

2194:    Input Parameters:
2195: +     a - the PetscInt value
2196: -     b - the second value

2198:    Output Parameter:
2199: .     result - the result as a PetscInt value, or NULL if you do not want the result, you just want to check if it overflows

2201:    Use PetscInt64Mult() to compute the product of two 32 bit PetscInt and store in a PetscInt64
2202:    Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt

2204:    Not available from Fortran

2206:    Developers Note: We currently assume that PetscInt addition does not overflow, this is obviously wrong but requires many more checks.

2208:    Level: advanced

2210: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64(), PetscIntSumError()
2211: @*/
2212: static inline PetscErrorCode PetscIntMultError(PetscInt a,PetscInt b,PetscInt *result)
2213: {
2214:   PetscInt64 r;

2216:   r = PetscInt64Mult(a,b);
2217: #if !defined(PETSC_USE_64BIT_INDICES)
2219: #endif
2220:   if (result) *result = (PetscInt)r;
2221:   return 0;
2222: }

2224: /*@C

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

2228:    Not Collective

2230:    Input Parameters:
2231: +     a - the PetscInt value
2232: -     b - the second value

2234:    Output Parameter:
2235: .     c - the result as a PetscInt value,  or NULL if you do not want the result, you just want to check if it overflows

2237:    Use PetscInt64Mult() to compute the product of two 32 bit PetscInt and store in a PetscInt64
2238:    Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt

2240:    Not available from Fortran

2242:    Level: advanced

2244: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError()
2245: @*/
2246: static inline PetscErrorCode PetscIntSumError(PetscInt a,PetscInt b,PetscInt *result)
2247: {
2248:   PetscInt64 r;

2250:   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2251: #if !defined(PETSC_USE_64BIT_INDICES)
2253: #endif
2254:   if (result) *result = (PetscInt) r;
2255:   return 0;
2256: }

2258: /*
2259:      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2260: */
2261: #if defined(hz)
2262: #  undef hz
2263: #endif

2265: #include <limits.h>

2267: /* The number of bits in a byte */

2269: #define PETSC_BITS_PER_BYTE CHAR_BIT

2271: #if defined(PETSC_HAVE_SYS_TYPES_H)
2272: #  include <sys/types.h>
2273: #endif

2275: /*MC

2277:     PETSC_VERSION - This manual page provides information about how PETSc documents and uses its version information. This information is available to both C/C++
2278:                     and Fortran compilers when petscsys.h is included.

2280:     The current PETSc version and the API for accessing it are defined in petscversion.h

2282:     The complete version number is given as the triple  PETSC_VERSION_MAJOR.PETSC_VERSION_MINOR.PETSC_VERSION_SUBMINOR (in short hand x.y.z)

2284:     A change in the minor version number (y) indicates possible/likely changes in the PETSc API. Note this is different than with the semantic versioning convention
2285:     where only a change in the major version number (x) indicates a change in the API.

2287:     A subminor greater than zero indicates a patch release. Version x.y.z maintains source and binary compatibility with version x.y.w for all z and w

2289:     Use the macros PETSC_VERSION_EQ(x,y,z), PETSC_VERSION_LT(x,y,z), PETSC_VERSION_LE(x,y,z), PETSC_VERSION_GT(x,y,z),
2290:     PETSC_VERSION_GE(x,y,z) to determine if the current version is equal to, less than, less than or equal to, greater than or greater than or equal to a given
2291:     version number (x.y.z).

2293:     PETSC_RELEASE_DATE is the date the x.y version was released (i.e. the version before any patch releases)

2295:     PETSC_VERSION_DATE is the date the x.y.z version was released

2297:     PETSC_VERSION_GIT is the last git commit to the repository given in the form vx.y.z-wwwww

2299:     PETSC_VERSION_DATE_GIT is the date of the last git commit to the repository

2301:     PETSC_VERSION_() is deprecated and will eventually be removed.

2303:     Level: intermediate

2305: M*/

2307: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2308: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2309: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2310: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2311: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2312: PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2313: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2314: PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt*,PetscInt*,PetscInt*,PetscInt*);

2316: PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt,const PetscInt[],PetscBool*);
2317: PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt,const PetscMPIInt[],PetscBool*);
2318: PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt,const PetscReal[],PetscBool*);
2319: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2320: PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt,PetscInt[]);
2321: PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt*,PetscInt[]);
2322: PETSC_EXTERN PetscErrorCode PetscSortedCheckDupsInt(PetscInt,const PetscInt[],PetscBool*);
2323: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2325: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2326: PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt*);
2327: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2328: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2329: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2330: PETSC_EXTERN PetscErrorCode PetscSortIntWithCountArray(PetscCount,PetscInt[],PetscCount[]);
2331: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt[],PetscInt[],PetscInt[]);
2332: PETSC_EXTERN PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount,PetscInt[],PetscInt[],PetscCount[]);
2333: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2334: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2335: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2336: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt,PetscMPIInt[],PetscInt[]);
2337: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2338: PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt,PetscInt[],void*,size_t,void*);
2339: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2340: PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt,PetscReal[],PetscInt[]);
2341: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2342: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2343: PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal,PetscInt,const PetscReal[],PetscReal,PetscInt*);
2344: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2345: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2346: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2347: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt[],const PetscInt[],PetscInt,const PetscInt[],const PetscInt[],PetscInt*,PetscInt**,PetscInt**);
2348: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt*,PetscInt**);
2349: PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt,const PetscMPIInt[],PetscInt,const PetscMPIInt[],PetscInt*,PetscMPIInt**);
2350: PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *);

2352: PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt,void*,size_t,int(*)(const void*,const void*,void*),void*);
2353: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt,PetscInt[]);
2354: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt,PetscMPIInt[]);
2355: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt,PetscReal[]);
2356: PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt,void*,size_t,void*,size_t,int(*)(const void*,const void*,void*),void*);
2357: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt,PetscInt[],PetscInt[]);
2358: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt,PetscMPIInt[],PetscMPIInt[]);
2359: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt,PetscReal[],PetscInt[]);

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

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

2367:    Level: beginner

2369:    Notes:
2370:    To use SPRNG or RANDOM123 you must have ./configure PETSc
2371:    with the option --download-sprng or --download-random123

2373: .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2374: J*/
2375: typedef const char* PetscRandomType;
2376: #define PETSCRAND       "rand"
2377: #define PETSCRAND48     "rand48"
2378: #define PETSCSPRNG      "sprng"
2379: #define PETSCRANDER48   "rander48"
2380: #define PETSCRANDOM123  "random123"
2381: #define PETSCCURAND     "curand"

2383: /* Logging support */
2384: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;

2386: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);

2388: /* Dynamic creation and loading functions */
2389: PETSC_EXTERN PetscFunctionList PetscRandomList;

2391: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2392: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom,PetscRandomType);
2393: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2394: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom,PetscRandomType*);
2395: PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom,PetscObject,const char[]);
2396: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);

2398: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2399: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2400: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2401: PETSC_EXTERN PetscErrorCode PetscRandomGetValues(PetscRandom,PetscInt,PetscScalar*);
2402: PETSC_EXTERN PetscErrorCode PetscRandomGetValuesReal(PetscRandom,PetscInt,PetscReal*);
2403: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2404: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2405: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2406: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2407: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2408: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);

2410: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2411: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2412: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2413: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2414: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2415: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2416: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);
2417: PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2418: PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2419: PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);

2421: static inline PetscBool PetscBinaryBigEndian(void) {long _petsc_v = 1; return ((char*)&_petsc_v)[0] ? PETSC_FALSE : PETSC_TRUE;}

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

2438: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2439: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2440: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);

2442: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2443: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool);
2444: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2445: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2446: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2447: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2448: PETSC_EXTERN PetscErrorCode PetscWaitOnError(void);

2450: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2451: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2452: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2453: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2454: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2455: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt*,const void*,PetscMPIInt*,PetscMPIInt**,void*)
2456:   PetscAttrMPIPointerWithType(6,3);
2457: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2458:                                                     PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2459:                                                     PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2460:   PetscAttrMPIPointerWithType(6,3);
2461: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2462:                                                        MPI_Request**,MPI_Request**,
2463:                                                        PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2464:                                                        PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2465:   PetscAttrMPIPointerWithType(6,3);

2467: PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2468: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2469: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);

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

2473: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);

2475: PETSC_EXTERN const char *const PetscSubcommTypes[];

2477: struct _n_PetscSubcomm {
2478:   MPI_Comm         parent;           /* parent communicator */
2479:   MPI_Comm         dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2480:   MPI_Comm         child;            /* the sub-communicator */
2481:   PetscMPIInt      n;                /* num of subcommunicators under the parent communicator */
2482:   PetscMPIInt      color;            /* color of processors belong to this communicator */
2483:   PetscMPIInt      *subsize;         /* size of subcommunicator[color] */
2484:   PetscSubcommType type;
2485:   char             *subcommprefix;
2486: };

2488: static inline MPI_Comm PetscSubcommParent(PetscSubcomm scomm) {return scomm->parent;}
2489: static inline MPI_Comm PetscSubcommChild(PetscSubcomm scomm) {return scomm->child;}
2490: static inline MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm) {return scomm->dupparent;}
2491: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2492: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2493: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2494: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2495: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2496: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2497: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2498: PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm,const char[]);
2499: PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm,MPI_Comm*);
2500: PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm,MPI_Comm*);
2501: PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm,MPI_Comm*);

2503: PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt,PetscHeap*);
2504: PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap,PetscInt,PetscInt);
2505: PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap,PetscInt*,PetscInt*);
2506: PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap,PetscInt*,PetscInt*);
2507: PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap,PetscInt,PetscInt);
2508: PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2509: PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap*);
2510: PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap,PetscViewer);

2512: PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2513: PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm,PetscShmComm*);
2514: PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2515: PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2516: PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm,MPI_Comm*);

2518: /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2519: PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm,PetscInt,PetscOmpCtrl*);
2520: PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl,MPI_Comm*,MPI_Comm*,PetscBool*);
2521: PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl*);
2522: PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2523: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2524: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);

2526: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2527: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2528: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2529: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2530: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2531: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2532: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2533: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);

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

2540: extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2541: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted*);
2542: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted*);
2543: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted,const char*,const char*,PetscBool*);

2545: #include <stdarg.h>
2546: PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list);
2547: PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE*,const char[],va_list);

2549: PETSC_EXTERN PetscSegBuffer PetscCitationsList;

2551: /*@C
2552:       PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code.

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

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

2560:    Level: intermediate

2562:    Not available from Fortran

2564:      Options Database:
2565: .     -citations [filename]   - print out the bibtex entries for the given computation
2566: @*/
2567: static inline PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2568: {
2569:   size_t  len;
2570:   char   *vstring;

2572:   if (set && *set) return 0;
2573:   PetscStrlen(cit,&len);
2574:   PetscSegBufferGet(PetscCitationsList,len,&vstring);
2575:   PetscArraycpy(vstring,cit,len);
2576:   if (set) *set = PETSC_TRUE;
2577:   return 0;
2578: }

2580: PETSC_EXTERN PetscErrorCode PetscURLShorten(const char[],char[],size_t);
2581: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm,char[],char[],size_t);
2582: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm,const char[],char[],size_t);
2583: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm,const char[],const char []);

2585: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm,char[],char[],size_t);
2586: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm,const char[],char[],char[],size_t);

2588: PETSC_EXTERN PetscErrorCode PetscGlobusGetTransfers(MPI_Comm,const char[],char[],size_t);

2590: PETSC_EXTERN PetscErrorCode PetscTextBelt(MPI_Comm,const char[],const char[],PetscBool*);
2591: PETSC_EXTERN PetscErrorCode PetscTellMyCell(MPI_Comm,const char[],const char[],PetscBool*);

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

2596: #if defined(PETSC_USE_DEBUG)
2597: static inline unsigned int PetscStrHash(const char *str)
2598: {
2599:   unsigned int c,hash = 5381;

2601:   while ((c = (unsigned int)*str++)) hash = ((hash << 5) + hash) + c; /* hash * 33 + c */
2602:   return hash;
2603: }

2605: /*MC
2606:    MPIU_Allreduce - a PETSc replacement for MPI_Allreduce() that tries to determine if the call from all the MPI processes occur from the
2607:                     same place in the PETSc code. This helps to detect bugs where different MPI processes follow different code paths
2608:                     resulting in inconsistent and incorrect calls to MPI_Allreduce().

2610:    Collective

2612:    Synopsis:
2613: #include <petscsys.h>
2614:      PetscErrorCode MPIU_Allreduce(void *indata,void *outdata,PetscMPIInt count,MPI_Datatype datatype, MPI_Op op, MPI_Comm comm);

2616:    Input Parameters:
2617: +  a - pointer to the input data to be reduced
2618: .  c - the number of MPI data items in a and b
2619: .  d - the MPI datatype, for example MPI_INT
2620: .  e - the MPI operation, for example MPI_SUM
2621: -  fcomm - the MPI communicator on which the operation occurs

2623:    Output Parameter:
2624: .  b - the reduced values

2626:    Notes:
2627:      In optimized mode this directly calls MPI_Allreduce()

2629:      This is defined as a macro that can return error codes internally so it cannot be used in a subroutine that returns void.

2631:      The error code this returns should be checked with PetscCall() even though it looks like an MPI function because it always returns PETSc error codes

2633:    Level: developer

2635: .seealso: MPI_Allreduce()
2636: M*/
2637: #define MPIU_Allreduce(a,b,c,d,e,fcomm) PetscMacroReturnStandard(                           \
2638:     PetscMPIInt a_b1[6],a_b2[6];                                                               \
2639:     int         _mpiu_allreduce_c_int = (int)c;                                                \
2640:     a_b1[0] = -(PetscMPIInt)__LINE__;                          a_b1[1] = -a_b1[0];             \
2641:     a_b1[2] = -(PetscMPIInt)PetscStrHash(PETSC_FUNCTION_NAME); a_b1[3] = -a_b1[2];             \
2642:     a_b1[4] = -(PetscMPIInt)(c);                               a_b1[5] = -a_b1[4];             \
2643:     MPI_Allreduce(a_b1,a_b2,6,MPI_INT,MPI_MAX,fcomm);                               \
2647:     MPI_Allreduce((a),(b),(c),d,e,(fcomm));                                         \
2648:   )
2649: #else
2650: #define MPIU_Allreduce(a,b,c,d,e,fcomm) PetscMacroReturnStandard(PetscCallMPI(MPI_Allreduce((a),(b),(c),d,e,(fcomm))))
2651: #endif

2653: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
2654: PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint,PetscMPIInt,MPI_Info,MPI_Comm,void*,MPI_Win*);
2655: PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win,PetscMPIInt,MPI_Aint*,PetscMPIInt*,void*);
2656: #endif

2658: /* this is a vile hack */
2659: #if defined(PETSC_HAVE_NECMPI)
2660: #if !defined(PETSC_NECMPI_VERSION_MAJOR) ||                              \
2661:     !defined(PETSC_NECMPI_VERSION_MINOR) ||                              \
2662:     PETSC_NECMPI_VERSION_MAJOR < 2       ||                              \
2663:     (PETSC_NECMPI_VERSION_MAJOR == 2 && PETSC_NECMPI_VERSION_MINOR < 18)
2664: #define MPI_Type_free(a) (*(a) = MPI_DATATYPE_NULL,0);
2665: #endif
2666: #endif

2668: /*
2669:     List of external packages and queries on it
2670: */
2671: PETSC_EXTERN PetscErrorCode  PetscHasExternalPackage(const char[],PetscBool*);

2673: /*
2674:  OpenMP support
2675: */
2676: #if defined(_OPENMP)
2677: #define PetscPragmaOMP(...) _Pragma(PetscStringize(omp __VA_ARGS__))
2678: #else // no OpenMP so no threads
2679: #define PetscPragmaOMP(...)
2680: #endif

2682: #endif