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:    Portions of this code are under:
  5:    Copyright (c) 2022 Advanced Micro Devices, Inc. All rights reserved.
  6: */
  7: #ifndef PETSCSYS_H
  8: #define PETSCSYS_H

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

 23: /* SUBMANSEC = Sys */

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

 44: #include <petscsystypes.h>

 46: /* ========================================================================== */

 48: /*
 49:     Defines the interface to MPI allowing the use of all MPI functions.

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

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

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

120: #if !defined(PETSC_MPI_WIN_FMT)
121:   #define PETSC_MPI_WIN_FMT "d"
122: #endif

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

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

135: PETSC_EXTERN MPI_Datatype MPIU_ENUM PETSC_ATTRIBUTE_MPI_TYPE_TAG(PetscEnum);
136: PETSC_EXTERN MPI_Datatype MPIU_BOOL PETSC_ATTRIBUTE_MPI_TYPE_TAG(PetscBool);

138: /*MC
139:    MPIU_INT - Portable MPI datatype corresponding to `PetscInt` independent of the precision of `PetscInt`

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

144:    Level: beginner

146: .seealso: `PetscReal`, `PetscScalar`, `PetscComplex`, `PetscInt`, `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`
147: M*/

149: PETSC_EXTERN MPI_Datatype MPIU_FORTRANADDR;

151: #if defined(PETSC_USE_64BIT_INDICES)
152:   #define MPIU_INT     MPIU_INT64
153:   #define PetscInt_FMT PetscInt64_FMT
154: #else
155:   #define MPIU_INT     MPI_INT
156:   #define PetscInt_FMT "d"
157: #endif

159: /*
160:     For the rare cases when one needs to send a size_t object with MPI
161: */
162: PETSC_EXTERN MPI_Datatype MPIU_SIZE_T PETSC_ATTRIBUTE_MPI_TYPE_TAG(size_t);

164: /*
165:       You can use PETSC_STDOUT as a replacement of stdout. You can also change
166:     the value of PETSC_STDOUT to redirect all standard output elsewhere
167: */
168: PETSC_EXTERN FILE *PETSC_STDOUT;

170: /*
171:       You can use PETSC_STDERR as a replacement of stderr. You can also change
172:     the value of PETSC_STDERR to redirect all standard error elsewhere
173: */
174: PETSC_EXTERN FILE *PETSC_STDERR;

176: /* PetscPragmaSIMD - from CeedPragmaSIMD */

178: #if defined(__NEC__)
179:   #define PetscPragmaSIMD _Pragma("_NEC ivdep")
180: #elif defined(__INTEL_COMPILER) && !defined(_WIN32)
181:   #define PetscPragmaSIMD _Pragma("vector")
182: #elif defined(__GNUC__) && __GNUC__ >= 5 && !defined(__PGI)
183:   #define PetscPragmaSIMD _Pragma("GCC ivdep")
184: #elif defined(_OPENMP) && _OPENMP >= 201307
185:   #if defined(_MSC_VER)
186:     #define PetscPragmaSIMD __pragma(omp simd)
187:   #else
188:     #define PetscPragmaSIMD _Pragma("omp simd")
189:   #endif
190: #elif defined(PETSC_HAVE_CRAY_VECTOR)
191:   #define PetscPragmaSIMD _Pragma("_CRI ivdep")
192: #else
193:   #define PetscPragmaSIMD
194: #endif

196: /*
197:     Declare extern C stuff after including external header files
198: */

200: PETSC_EXTERN const char *const PetscBools[];

202: PETSC_EXTERN PetscBool PETSC_RUNNING_ON_VALGRIND;
203: /*
204:     Defines elementary mathematics functions and constants.
205: */
206: #include <petscmath.h>

208: PETSC_EXTERN const char *const PetscCopyModes[];

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

213:    Level: beginner

215:    Note:
216:    Accepted by many PETSc functions to not set a parameter and instead use some default

218:    Fortran Note:
219:    This macro does not exist in Fortran; you must use `PETSC_NULL_INTEGER`, `PETSC_NULL_DOUBLE_PRECISION` etc

221: .seealso: `PETSC_DECIDE`, `PETSC_DEFAULT`, `PETSC_DETERMINE`

223: M*/
224: #define PETSC_IGNORE NULL

226: /* This is deprecated */
227: #define PETSC_NULL NULL

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

233:    Level: beginner

235: .seealso: `PETSC_DEFAULT`, `PETSC_IGNORE`, `PETSC_DETERMINE`

237: M*/
238: #define PETSC_DECIDE -1

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

244:    Level: beginner

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

250: .seealso: `PETSC_DECIDE`, `PETSC_DEFAULT`, `PETSC_IGNORE`, `VecSetSizes()`

252: M*/
253: #define PETSC_DETERMINE PETSC_DECIDE

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

259:    Level: beginner

261:    Fortran Note:
262:    You need to use `PETSC_DEFAULT_INTEGER` or `PETSC_DEFAULT_REAL`.

264: .seealso: `PETSC_DECIDE`, `PETSC_IGNORE`, `PETSC_DETERMINE`

266: M*/
267: #define PETSC_DEFAULT -2

269: /*MC
270:     PETSC_COMM_WORLD - the equivalent of the `MPI_COMM_WORLD` communicator which represents
271:            all the processes that PETSc knows about.

273:    Level: beginner

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

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

284: .seealso: `PETSC_COMM_SELF`

286: M*/
287: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;

289: /*MC
290:     PETSC_COMM_SELF - This is always `MPI_COMM_SELF`

292:    Level: beginner

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

297: .seealso: `PETSC_COMM_WORLD`

299: M*/
300: #define PETSC_COMM_SELF MPI_COMM_SELF

302: /*MC
303:     PETSC_MPI_THREAD_REQUIRED - the required threading support used if PETSc initializes
304:            MPI with `MPI_Init_thread()`.

306:    Level: beginner

308:    Notes:
309:    By default `PETSC_MPI_THREAD_REQUIRED` equals `MPI_THREAD_FUNNELED`.

311: .seealso: `PetscInitialize()`

313: M*/
314: PETSC_EXTERN PetscMPIInt PETSC_MPI_THREAD_REQUIRED;

316: PETSC_EXTERN PetscBool PetscBeganMPI;
317: PETSC_EXTERN PetscBool PetscErrorHandlingInitialized;
318: PETSC_EXTERN PetscBool PetscInitializeCalled;
319: PETSC_EXTERN PetscBool PetscFinalizeCalled;
320: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;

322: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm), PetscErrorCode (*)(MPI_Comm));
323: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm, MPI_Comm *, int *);
324: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm *);
325: PETSC_EXTERN PetscErrorCode PetscCommGetComm(MPI_Comm, MPI_Comm *);
326: PETSC_EXTERN PetscErrorCode PetscCommRestoreComm(MPI_Comm, MPI_Comm *);

328: #if defined(PETSC_HAVE_KOKKOS)
329: PETSC_EXTERN PetscErrorCode PetscKokkosInitializeCheck(void); /* Initialize Kokkos if not yet. */
330: #endif

332: #if defined(PETSC_HAVE_NVSHMEM)
333: PETSC_EXTERN PetscBool      PetscBeganNvshmem;
334: PETSC_EXTERN PetscBool      PetscNvshmemInitialized;
335: PETSC_EXTERN PetscErrorCode PetscNvshmemFinalize(void);
336: #endif

338: #if defined(PETSC_HAVE_ELEMENTAL)
339: PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
340: PETSC_EXTERN PetscErrorCode PetscElementalInitialized(PetscBool *);
341: PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
342: #endif

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

347:    Synopsis:
348: #include <petscsys.h>
349:    PetscErrorCode PetscMalloc(size_t m,void **result)

351:    Not Collective

353:    Input Parameter:
354: .  m - number of bytes to allocate

356:    Output Parameter:
357: .  result - memory allocated

359:    Level: beginner

361:    Notes:
362:    Memory is always allocated at least double aligned

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

366: .seealso: `PetscFree()`, `PetscNew()`

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

371: /*MC
372:    PetscRealloc - Reallocates memory

374:    Synopsis:
375: #include <petscsys.h>
376:    PetscErrorCode PetscRealloc(size_t m,void **result)

378:    Not Collective

380:    Input Parameters:
381: +  m - number of bytes to allocate
382: -  result - previous memory

384:    Output Parameter:
385: .  result - new memory allocated

387:    Level: developer

389:    Notes:
390:    Memory is always allocated at least double aligned

392: .seealso: `PetscMalloc()`, `PetscFree()`, `PetscNew()`

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

397: /*MC
398:    PetscAddrAlign - Rounds up an address to `PETSC_MEMALIGN` alignment

400:    Synopsis:
401: #include <petscsys.h>
402:    void *PetscAddrAlign(void *addr)

404:    Not Collective

406:    Input Parameters:
407: .  addr - address to align (any pointer type)

409:    Level: developer

411: .seealso: `PetscMallocAlign()`

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

416: /*MC
417:    PetscCalloc - Allocates a cleared (zeroed) memory region aligned to `PETSC_MEMALIGN`

419:    Synopsis:
420: #include <petscsys.h>
421:    PetscErrorCode PetscCalloc(size_t m,void **result)

423:    Not Collective

425:    Input Parameter:
426: .  m - number of bytes to allocate

428:    Output Parameter:
429: .  result - memory allocated

431:    Level: beginner

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

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

438: .seealso: `PetscFree()`, `PetscNew()`

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

443: /*MC
444:    PetscMalloc1 - Allocates an array of memory aligned to `PETSC_MEMALIGN`

446:    Synopsis:
447: #include <petscsys.h>
448:    PetscErrorCode PetscMalloc1(size_t m1,type **r1)

450:    Not Collective

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

455:    Output Parameter:
456: .  r1 - memory allocated

458:    Note:
459:    This uses the sizeof() of the memory type requested to determine the total memory to be allocated, therefore you should not
460:          multiply the number of elements requested by the `sizeof()` the type. For example use
461: $  PetscInt *id;
462: $  PetscMalloc1(10,&id);
463:         not
464: $  PetscInt *id;
465: $  PetscMalloc1(10*sizeof(PetscInt),&id);

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

469:    Level: beginner

471: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc2()`

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

476: /*MC
477:    PetscCalloc1 - Allocates a cleared (zeroed) array of memory aligned to `PETSC_MEMALIGN`

479:    Synopsis:
480: #include <petscsys.h>
481:    PetscErrorCode PetscCalloc1(size_t m1,type **r1)

483:    Not Collective

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

488:    Output Parameter:
489: .  r1 - memory allocated

491:    Notes:
492:    See `PetsMalloc1()` for more details on usage.

494:    Level: beginner

496: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc2()`

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

501: /*MC
502:    PetscMalloc2 - Allocates 2 arrays of memory both aligned to `PETSC_MEMALIGN`

504:    Synopsis:
505: #include <petscsys.h>
506:    PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)

508:    Not Collective

510:    Input Parameters:
511: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
512: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

514:    Output Parameters:
515: +  r1 - memory allocated in first chunk
516: -  r2 - memory allocated in second chunk

518:    Level: developer

520: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc2()`

522: M*/
523: #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))

525: /*MC
526:    PetscCalloc2 - Allocates 2 cleared (zeroed) arrays of memory both aligned to `PETSC_MEMALIGN`

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

532:    Not Collective

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

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

542:    Level: developer

544: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc2()`

546: M*/
547: #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))

549: /*MC
550:    PetscMalloc3 - Allocates 3 arrays of memory, all aligned to `PETSC_MEMALIGN`

552:    Synopsis:
553: #include <petscsys.h>
554:    PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

556:    Not Collective

558:    Input Parameters:
559: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
560: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
561: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

563:    Output Parameters:
564: +  r1 - memory allocated in first chunk
565: .  r2 - memory allocated in second chunk
566: -  r3 - memory allocated in third chunk

568:    Level: developer

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

572: M*/
573: #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))

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

578:    Synopsis:
579: #include <petscsys.h>
580:    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

582:    Not Collective

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

589:    Output Parameters:
590: +  r1 - memory allocated in first chunk
591: .  r2 - memory allocated in second chunk
592: -  r3 - memory allocated in third chunk

594:    Level: developer

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

598: M*/
599: #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))

601: /*MC
602:    PetscMalloc4 - Allocates 4 arrays of memory, all aligned to `PETSC_MEMALIGN`

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

608:    Not Collective

610:    Input Parameters:
611: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
612: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
613: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
614: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

616:    Output Parameters:
617: +  r1 - memory allocated in first chunk
618: .  r2 - memory allocated in second chunk
619: .  r3 - memory allocated in third chunk
620: -  r4 - memory allocated in fourth chunk

622:    Level: developer

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

626: M*/
627: #define PetscMalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \
628:   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))

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

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

637:    Not Collective

639:    Input Parameters:
640: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
641: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
642: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
643: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

645:    Output Parameters:
646: +  r1 - memory allocated in first chunk
647: .  r2 - memory allocated in second chunk
648: .  r3 - memory allocated in third chunk
649: -  r4 - memory allocated in fourth chunk

651:    Level: developer

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

655: M*/
656: #define PetscCalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \
657:   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))

659: /*MC
660:    PetscMalloc5 - Allocates 5 arrays of memory, all aligned to `PETSC_MEMALIGN`

662:    Synopsis:
663: #include <petscsys.h>
664:    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)

666:    Not Collective

668:    Input Parameters:
669: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
670: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
671: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
672: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
673: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

675:    Output Parameters:
676: +  r1 - memory allocated in first chunk
677: .  r2 - memory allocated in second chunk
678: .  r3 - memory allocated in third chunk
679: .  r4 - memory allocated in fourth chunk
680: -  r5 - memory allocated in fifth chunk

682:    Level: developer

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

686: M*/
687: #define PetscMalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \
688:   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))

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

693:    Synopsis:
694: #include <petscsys.h>
695:    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)

697:    Not Collective

699:    Input Parameters:
700: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
701: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
702: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
703: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
704: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

706:    Output Parameters:
707: +  r1 - memory allocated in first chunk
708: .  r2 - memory allocated in second chunk
709: .  r3 - memory allocated in third chunk
710: .  r4 - memory allocated in fourth chunk
711: -  r5 - memory allocated in fifth chunk

713:    Level: developer

715: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc5()`, `PetscFree5()`

717: M*/
718: #define PetscCalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \
719:   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))

721: /*MC
722:    PetscMalloc6 - Allocates 6 arrays of memory, all aligned to `PETSC_MEMALIGN`

724:    Synopsis:
725: #include <petscsys.h>
726:    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)

728:    Not Collective

730:    Input Parameters:
731: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
732: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
733: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
734: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
735: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
736: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

738:    Output Parameteasr:
739: +  r1 - memory allocated in first chunk
740: .  r2 - memory allocated in second chunk
741: .  r3 - memory allocated in third chunk
742: .  r4 - memory allocated in fourth chunk
743: .  r5 - memory allocated in fifth chunk
744: -  r6 - memory allocated in sixth chunk

746:    Level: developer

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

750: M*/
751: #define PetscMalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \
752:   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))

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

757:    Synopsis:
758: #include <petscsys.h>
759:    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)

761:    Not Collective

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

771:    Output Parameters:
772: +  r1 - memory allocated in first chunk
773: .  r2 - memory allocated in second chunk
774: .  r3 - memory allocated in third chunk
775: .  r4 - memory allocated in fourth chunk
776: .  r5 - memory allocated in fifth chunk
777: -  r6 - memory allocated in sixth chunk

779:    Level: developer

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

783: M*/
784: #define PetscCalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \
785:   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))

787: /*MC
788:    PetscMalloc7 - Allocates 7 arrays of memory, all aligned to `PETSC_MEMALIGN`

790:    Synopsis:
791: #include <petscsys.h>
792:    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)

794:    Not Collective

796:    Input Parameters:
797: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
798: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
799: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
800: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
801: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
802: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
803: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

805:    Output Parameters:
806: +  r1 - memory allocated in first chunk
807: .  r2 - memory allocated in second chunk
808: .  r3 - memory allocated in third chunk
809: .  r4 - memory allocated in fourth chunk
810: .  r5 - memory allocated in fifth chunk
811: .  r6 - memory allocated in sixth chunk
812: -  r7 - memory allocated in seventh chunk

814:    Level: developer

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

818: M*/
819: #define PetscMalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \
820:   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))

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

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

829:    Not Collective

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

840:    Output Parameters:
841: +  r1 - memory allocated in first chunk
842: .  r2 - memory allocated in second chunk
843: .  r3 - memory allocated in third chunk
844: .  r4 - memory allocated in fourth chunk
845: .  r5 - memory allocated in fifth chunk
846: .  r6 - memory allocated in sixth chunk
847: -  r7 - memory allocated in seventh chunk

849:    Level: developer

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

853: M*/
854: #define PetscCalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \
855:   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))

857: /*MC
858:    PetscNew - Allocates memory of a particular type, zeros the memory! Aligned to `PETSC_MEMALIGN`

860:    Synopsis:
861: #include <petscsys.h>
862:    PetscErrorCode PetscNew(type **result)

864:    Not Collective

866:    Output Parameter:
867: .  result - memory allocated, sized to match pointer type

869:    Level: beginner

871: .seealso: `PetscFree()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc1()`

873: M*/
874: #define PetscNew(b) PetscCalloc1(1, (b))

876: #define PetscNewLog(o, b) PETSC_DEPRECATED_MACRO("GCC warning \"PetscNewLog is deprecated, use PetscNew() instead (since version 3.18)\"") PetscNew((b))

878: /*MC
879:    PetscFree - Frees memory

881:    Synopsis:
882: #include <petscsys.h>
883:    PetscErrorCode PetscFree(void *memory)

885:    Not Collective

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

890:    Level: beginner

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

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

897: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc1()`

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

902: /*MC
903:    PetscFree2 - Frees 2 chunks of memory obtained with `PetscMalloc2()`

905:    Synopsis:
906: #include <petscsys.h>
907:    PetscErrorCode PetscFree2(void *memory1,void *memory2)

909:    Not Collective

911:    Input Parameters:
912: +   memory1 - memory to free
913: -   memory2 - 2nd memory to free

915:    Level: developer

917:    Notes:
918:     Memory must have been obtained with `PetscMalloc2()`

920: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`

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

925: /*MC
926:    PetscFree3 - Frees 3 chunks of memory obtained with `PetscMalloc3()`

928:    Synopsis:
929: #include <petscsys.h>
930:    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)

932:    Not Collective

934:    Input Parameters:
935: +   memory1 - memory to free
936: .   memory2 - 2nd memory to free
937: -   memory3 - 3rd memory to free

939:    Level: developer

941:    Notes:
942:     Memory must have been obtained with `PetscMalloc3()`

944: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`

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

949: /*MC
950:    PetscFree4 - Frees 4 chunks of memory obtained with `PetscMalloc4()`

952:    Synopsis:
953: #include <petscsys.h>
954:    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)

956:    Not Collective

958:    Input Parameters:
959: +   m1 - memory to free
960: .   m2 - 2nd memory to free
961: .   m3 - 3rd memory to free
962: -   m4 - 4th memory to free

964:    Level: developer

966:    Notes:
967:     Memory must have been obtained with `PetscMalloc4()`

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

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

974: /*MC
975:    PetscFree5 - Frees 5 chunks of memory obtained with `PetscMalloc5()`

977:    Synopsis:
978: #include <petscsys.h>
979:    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)

981:    Not Collective

983:    Input Parameters:
984: +   m1 - memory to free
985: .   m2 - 2nd memory to free
986: .   m3 - 3rd memory to free
987: .   m4 - 4th memory to free
988: -   m5 - 5th memory to free

990:    Level: developer

992:    Notes:
993:     Memory must have been obtained with `PetscMalloc5()`

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

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

1000: /*MC
1001:    PetscFree6 - Frees 6 chunks of memory obtained with `PetscMalloc6()`

1003:    Synopsis:
1004: #include <petscsys.h>
1005:    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)

1007:    Not Collective

1009:    Input Parameters:
1010: +   m1 - memory to free
1011: .   m2 - 2nd memory to free
1012: .   m3 - 3rd memory to free
1013: .   m4 - 4th memory to free
1014: .   m5 - 5th memory to free
1015: -   m6 - 6th memory to free

1017:    Level: developer

1019:    Notes:
1020:     Memory must have been obtained with `PetscMalloc6()`

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

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

1027: /*MC
1028:    PetscFree7 - Frees 7 chunks of memory obtained with `PetscMalloc7()`

1030:    Synopsis:
1031: #include <petscsys.h>
1032:    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)

1034:    Not Collective

1036:    Input Parameters:
1037: +   m1 - memory to free
1038: .   m2 - 2nd memory to free
1039: .   m3 - 3rd memory to free
1040: .   m4 - 4th memory to free
1041: .   m5 - 5th memory to free
1042: .   m6 - 6th memory to free
1043: -   m7 - 7th memory to free

1045:    Level: developer

1047:    Notes:
1048:     Memory must have been obtained with `PetscMalloc7()`

1050: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()`,
1051:           `PetscMalloc7()`

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

1056: PETSC_EXTERN PetscErrorCode PetscMallocA(int, PetscBool, int, const char *, const char *, size_t, void *, ...);
1057: PETSC_EXTERN PetscErrorCode PetscFreeA(int, int, const char *, const char *, void *, ...);
1058: PETSC_EXTERN                PetscErrorCode (*PetscTrMalloc)(size_t, PetscBool, int, const char[], const char[], void **);
1059: PETSC_EXTERN                PetscErrorCode (*PetscTrFree)(void *, int, const char[], const char[]);
1060: PETSC_EXTERN                PetscErrorCode (*PetscTrRealloc)(size_t, int, const char[], const char[], void **);
1061: PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1062: 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 **));
1063: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);

1065: /*
1066:   Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1067: */
1068: PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1069: PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);
1070: #if defined(PETSC_HAVE_CUDA)
1071: PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void);
1072: PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void);
1073: #endif
1074: #if defined(PETSC_HAVE_HIP)
1075: PETSC_EXTERN PetscErrorCode PetscMallocSetHIPHost(void);
1076: PETSC_EXTERN PetscErrorCode PetscMallocResetHIPHost(void);
1077: #endif

1079: #define MPIU_PETSCLOGDOUBLE  MPI_DOUBLE
1080: #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION

1082: /*
1083:    Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1084: */
1085: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1086: PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1087: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1088: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1089: PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1090: PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int, PetscLogDouble *);
1091: PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool, PetscBool);
1092: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool *, PetscBool *, PetscBool *);
1093: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int, const char[], const char[]);
1094: PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1095: PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool *);
1096: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool);
1097: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool *);

1099: PETSC_EXTERN const char *const PetscDataTypes[];
1100: PETSC_EXTERN PetscErrorCode    PetscDataTypeToMPIDataType(PetscDataType, MPI_Datatype *);
1101: PETSC_EXTERN PetscErrorCode    PetscMPIDataTypeToPetscDataType(MPI_Datatype, PetscDataType *);
1102: PETSC_EXTERN PetscErrorCode    PetscDataTypeGetSize(PetscDataType, size_t *);
1103: PETSC_EXTERN PetscErrorCode    PetscDataTypeFromString(const char *, PetscDataType *, PetscBool *);

1105: /*
1106:     Basic memory and string operations. These are usually simple wrappers
1107:    around the basic Unix system calls, but a few of them have additional
1108:    functionality and/or error checking.
1109: */
1110: PETSC_EXTERN PetscErrorCode PetscMemcmp(const void *, const void *, size_t, PetscBool *);
1111: PETSC_EXTERN PetscErrorCode PetscStrlen(const char[], size_t *);
1112: PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[], char, int *, char ***);
1113: PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int, char **);
1114: PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[], const char[], PetscBool *);
1115: PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[], const char[], PetscBool *);
1116: PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[], const char[], PetscBool *);
1117: PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[], const char[], size_t, PetscBool *);
1118: PETSC_EXTERN PetscErrorCode PetscStrcpy(char[], const char[]);
1119: PETSC_EXTERN PetscErrorCode PetscStrcat(char[], const char[]);
1120: PETSC_EXTERN PetscErrorCode PetscStrlcat(char[], const char[], size_t);
1121: PETSC_EXTERN PetscErrorCode PetscStrncpy(char[], const char[], size_t);
1122: PETSC_EXTERN PetscErrorCode PetscStrchr(const char[], char, char *[]);
1123: PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1124: PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1125: PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[], char, char *[]);
1126: PETSC_EXTERN PetscErrorCode PetscStrstr(const char[], const char[], char *[]);
1127: PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[], const char[], char *[]);
1128: PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[], const char[], PetscBool *);
1129: PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[], const char[], PetscBool *);
1130: PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[], const char *const *, PetscInt *);
1131: PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[], char *[]);
1132: PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const *, char ***);
1133: PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char ***);
1134: PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt, const char *const *, char ***);
1135: PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt, char ***);
1136: PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm, const char[], char[], size_t);

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

1140: PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[], const char, PetscToken *);
1141: PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken, char *[]);
1142: PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken *);

1144: PETSC_EXTERN PetscErrorCode PetscStrInList(const char[], const char[], char, PetscBool *);
1145: PETSC_EXTERN const char    *PetscBasename(const char[]);
1146: PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt, const char *const *, const char *, PetscInt *, PetscBool *);
1147: PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const *, const char *, PetscEnum *, PetscBool *);

1149: /*
1150:    These are MPI operations for MPI_Allreduce() etc
1151: */
1152: PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1153: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1154: PETSC_EXTERN MPI_Op MPIU_SUM;
1155: PETSC_EXTERN MPI_Op MPIU_MAX;
1156: PETSC_EXTERN MPI_Op MPIU_MIN;
1157: #else
1158:   #define MPIU_SUM MPI_SUM
1159:   #define MPIU_MAX MPI_MAX
1160:   #define MPIU_MIN MPI_MIN
1161: #endif
1162: PETSC_EXTERN MPI_Op         Petsc_Garbage_SetIntersectOp;
1163: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm, const PetscInt[], PetscInt *, PetscInt *);

1165: #if defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16)
1166: /*MC
1167:     MPIU_SUM___FP16___FLOAT128 - MPI_Op that acts as a replacement for MPI_SUM with
1168:     custom MPI_Datatype `MPIU___FLOAT128`, `MPIU___COMPLEX128`, and `MPIU___FP16`.

1170:    Level: advanced

1172:    Developer Note:
1173:    This should be unified with `MPIU_SUM`

1175: .seealso: `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`
1176: M*/
1177: PETSC_EXTERN MPI_Op MPIU_SUM___FP16___FLOAT128;
1178: #endif
1179: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm, const PetscInt[], PetscInt *, PetscInt *);

1181: PETSC_EXTERN PetscErrorCode MPIULong_Send(void *, PetscInt, MPI_Datatype, PetscMPIInt, PetscMPIInt, MPI_Comm) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(1, 3);
1182: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void *, PetscInt, MPI_Datatype, PetscMPIInt, PetscMPIInt, MPI_Comm) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(1, 3);

1184: PETSC_EXTERN const char *const PetscFileModes[];

1186: /*
1187:     Defines PETSc error handling.
1188: */
1189: #include <petscerror.h>

1191: PETSC_EXTERN PetscBool   PetscCIEnabled;                    /* code is running in the PETSc test harness CI */
1192: PETSC_EXTERN PetscBool   PetscCIEnabledPortableErrorOutput; /* error output is stripped to ensure portability of error messages across systems */
1193: PETSC_EXTERN const char *PetscCIFilename(const char *);
1194: PETSC_EXTERN int         PetscCILinenumber(int);

1196: #define PETSC_SMALLEST_CLASSID 1211211
1197: PETSC_EXTERN PetscClassId   PETSC_LARGEST_CLASSID;
1198: PETSC_EXTERN PetscClassId   PETSC_OBJECT_CLASSID;
1199: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[], PetscClassId *);
1200: PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject, PetscObjectId *);
1201: PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject, PetscObjectId, PetscBool *);

1203: /*
1204:    Routines that get memory usage information from the OS
1205: */
1206: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1207: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1208: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1209: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);

1211: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);

1213: /*
1214:    Initialization of PETSc
1215: */
1216: PETSC_EXTERN PetscErrorCode PetscInitialize(int *, char ***, const char[], const char[]);
1217: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int, char **, const char[], const char[]);
1218: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1219: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1220: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1221: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1222: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1223: PETSC_EXTERN PetscErrorCode PetscGetArgs(int *, char ***);
1224: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1225: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);

1227: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1228: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);

1230: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[], const char[]);
1231: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1232: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1233: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject, const char[]);

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

1237: /*
1238:      These are so that in extern C code we can caste function pointers to non-extern C
1239:    function pointers. Since the regular C++ code expects its function pointers to be C++
1240: */
1241: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1242: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1243: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);

1245: /*
1246:     Functions that can act on any PETSc object.
1247: */
1248: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject *);
1249: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject, MPI_Comm *);
1250: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject, PetscClassId *);
1251: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject, const char *[]);
1252: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject, const char *[]);
1253: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject, const char[]);
1254: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject, const char *[]);
1255: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject, PetscInt);
1256: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject, PetscInt *);
1257: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject, PetscObject, PetscInt);
1258: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1259: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject, PetscInt *);
1260: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1261: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject, PetscMPIInt *);
1262: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject, const char[], PetscObject);
1263: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject, const char[]);
1264: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject, const char[], PetscObject *);
1265: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject, const char[], void (*)(void));
1266: #define PetscObjectComposeFunction(a, b, d) PetscObjectComposeFunction_Private(a, b, (PetscVoidFunction)(d))
1267: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1268: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1269: PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1270: PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject, PetscObject);
1271: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm, PetscMPIInt *);

1273: #include <petscviewertypes.h>
1274: #include <petscoptions.h>

1276: PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer, PetscBool, PetscLogDouble);
1277: PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool *);

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

1281: PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer, const char[]);
1282: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject, PetscViewer);
1283: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject, PetscViewer);
1284: #define PetscObjectQueryFunction(obj, name, fptr) PetscObjectQueryFunction_Private((obj), (name), (PetscVoidFunction *)(fptr))
1285: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject, const char[], void (**)(void));
1286: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject, const char[]);
1287: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject, const char[]);
1288: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject, const char[]);
1289: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject, const char *[]);
1290: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject, const char[]);
1291: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1292: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1293: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject, PetscObject, const char[]);
1294: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1295: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject, const char[], PetscBool *);
1296: PETSC_EXTERN PetscErrorCode PetscObjectObjectTypeCompare(PetscObject, PetscObject, PetscBool *);
1297: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject, const char[], PetscBool *);
1298: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject, PetscBool *, const char[], ...);
1299: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject, PetscBool *, const char[], ...);
1300: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1301: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);

1303: #if defined(PETSC_HAVE_SAWS)
1304: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1305: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1306: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject, PetscBool);
1307: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1308: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1309: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1310: PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1311: PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1312: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1313: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);

1315: #else
1316:   #define PetscSAWsBlock()                  0
1317:   #define PetscObjectSAWsViewOff(obj)       0
1318:   #define PetscObjectSAWsSetBlock(obj, flg) 0
1319:   #define PetscObjectSAWsBlock(obj)         0
1320:   #define PetscObjectSAWsGrantAccess(obj)   0
1321:   #define PetscObjectSAWsTakeAccess(obj)    0
1322:   #define PetscStackViewSAWs()              0
1323:   #define PetscStackSAWsViewOff()           0
1324:   #define PetscStackSAWsTakeAccess()
1325:   #define PetscStackSAWsGrantAccess()

1327: #endif

1329: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[], PetscDLMode, PetscDLHandle *);
1330: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1331: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle, const char[], void **);
1332: PETSC_EXTERN PetscErrorCode PetscDLAddr(void (*)(void), char **);
1333: #ifdef PETSC_HAVE_CXX
1334: PETSC_EXTERN PetscErrorCode PetscDemangleSymbol(const char *, char **);
1335: #endif

1337: #if defined(PETSC_USE_DEBUG)
1338: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void *, PetscStack **);
1339: #endif
1340: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE *, PetscBool);

1342: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList *);
1343: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList, const char[], PetscObject *);
1344: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList, PetscObject, char **, PetscBool *);
1345: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *, const char[], PetscObject);
1346: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *, const char[]);
1347: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList, PetscObjectList *);

1349: /*
1350:     Dynamic library lists. Lists of names of routines in objects or in dynamic
1351:   link libraries that will be loaded as needed.
1352: */

1354: #define PetscFunctionListAdd(list, name, fptr) PetscFunctionListAdd_Private((list), (name), (PetscVoidFunction)(fptr))
1355: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList *, const char[], void (*)(void));
1356: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList *);
1357: PETSC_EXTERN PetscErrorCode PetscFunctionListClear(PetscFunctionList);
1358: #define PetscFunctionListFind(list, name, fptr) PetscFunctionListFind_Private((list), (name), (PetscVoidFunction *)(fptr))
1359: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList, const char[], void (**)(void));
1360: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm, FILE *, const char[], const char[], const char[], const char[], PetscFunctionList, const char[], const char[]);
1361: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList, PetscFunctionList *);
1362: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList, PetscViewer);
1363: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList, const char ***, int *);
1364: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintNonEmpty(PetscFunctionList);
1365: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintAll(void);

1367: PETSC_EXTERN PetscDLLibrary PetscDLLibrariesLoaded;
1368: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm, PetscDLLibrary *, const char[]);
1369: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm, PetscDLLibrary *, const char[]);
1370: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm, PetscDLLibrary *, const char[], const char[], void **);
1371: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1372: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm, const char[], char *, size_t, PetscBool *);
1373: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm, const char[], PetscDLLibrary *);
1374: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);

1376: /*
1377:      Useful utility routines
1378: */
1379: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm, PetscInt *, PetscInt *);
1380: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm, PetscInt, PetscInt *, PetscInt *);
1381: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm, PetscInt *, PetscInt *);
1382: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm, PetscMPIInt);
1383: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm, PetscMPIInt);
1384: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1385: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE *);
1386: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm, const PetscInt[2], PetscInt[2]);
1387: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm, const PetscReal[2], PetscReal[2]);

1389: /*MC
1390:     PetscNot - negates a logical type value and returns result as a `PetscBool`

1392:     Level: beginner

1394:     Note:
1395:     This is useful in cases like
1396: .vb
1397:      int        *a;
1398:      PetscBool  flag = PetscNot(a)
1399: .ve
1400:      where !a would not return a `PetscBool` because we cannot provide a cast from int to `PetscBool` in C.

1402: .seealso: `PetscBool`, `PETSC_TRUE`, `PETSC_FALSE`
1403: M*/
1404: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)

1406: /*MC
1407:     PetscHelpPrintf - Prints help messages.

1409:     Not Collective, only applies on rank 0; No Fortran Support

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

1415:     Input Parameters:
1416: +  comm - the MPI communicator over which the help message is printed
1417: .  format - the usual printf() format string
1418: -  args - arguments to be printed

1420:    Level: developer

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

1425:       To use, write your own function, for example,
1426: $PetscErrorCode mypetschelpprintf(MPI_Comm comm,const char format[],....)
1427: ${
1428: $ return 0;
1429: $}
1430: then do the assignment
1431: $    PetscHelpPrintf = mypetschelpprintf;
1432:    You can do the assignment before `PetscInitialize()`.

1434:   The default routine used is called `PetscHelpPrintfDefault()`.

1436: .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscErrorPrintf()`
1437: M*/
1438: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);

1440: /*
1441:      Defines PETSc profiling.
1442: */
1443: #include <petsclog.h>

1445: /*
1446:       Simple PETSc parallel IO for ASCII printing
1447: */
1448: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[], char[]);
1449: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm, const char[], const char[], FILE **);
1450: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm, FILE *);
1451: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1452: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1453: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char *, size_t, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1454: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char *, size_t, const char[], size_t *, ...) PETSC_ATTRIBUTE_FORMAT(3, 5);
1455: PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[], size_t, const char *, PetscInt, const PetscReal[]);

1457: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1458: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1459: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);

1461: PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char *, size_t *);
1462: PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char *, char *);

1464: #if defined(PETSC_HAVE_POPEN)
1465: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm, const char[], const char[], const char[], FILE **);
1466: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm, FILE *);
1467: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1468: #endif

1470: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1471: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1472: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm, FILE *);
1473: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm, FILE *, size_t, char[]);
1474: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm, const char[], const char[], FILE **);
1475: PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm, const char[], const char[], FILE **);
1476: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char *[]);

1478: PETSC_EXTERN PetscClassId   PETSC_CONTAINER_CLASSID;
1479: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer, void **);
1480: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer, void *);
1481: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer *);
1482: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm, PetscContainer *);
1483: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void *));
1484: PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void *);

1486: /*
1487:    For use in debuggers
1488: */
1489: PETSC_EXTERN PetscMPIInt    PetscGlobalRank;
1490: PETSC_EXTERN PetscMPIInt    PetscGlobalSize;
1491: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt, const PetscInt[], PetscViewer);
1492: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt, const PetscReal[], PetscViewer);
1493: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt, const PetscScalar[], PetscViewer);

1495: #include <stddef.h>
1496: #include <string.h> /* for memcpy, memset */
1497: #include <stdlib.h>

1499: #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1500:   #include <xmmintrin.h>
1501: #endif

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

1508:    Not Collective

1510:    Input Parameters:
1511: +  b - pointer to initial memory space
1512: .  a - pointer to copy space
1513: -  n - length (in bytes) of space to copy

1515:    Level: intermediate

1517:    Notes:
1518:    `PetscArraymove()` is preferred

1520:    This routine is analogous to `memmove()`.

1522:    Developers Note:
1523:    This is inlined for performance

1525: .seealso: `PetscMemcpy()`, `PetscMemcmp()`, `PetscArrayzero()`, `PetscMemzero()`, `PetscArraycmp()`, `PetscArraycpy()`, `PetscStrallocpy()`,
1526:           `PetscArraymove()`
1527: @*/
1528: static inline PetscErrorCode PetscMemmove(void *a, const void *b, size_t n)
1529: {
1530:   if (n > 0) {
1533:   }
1534: #if !defined(PETSC_HAVE_MEMMOVE)
1535:   if (a < b) {
1536:     if (a <= b - n) memcpy(a, b, n);
1537:     else {
1538:       memcpy(a, b, (int)(b - a));
1539:       PetscMemmove(b, b + (int)(b - a), n - (int)(b - a));
1540:     }
1541:   } else {
1542:     if (b <= a - n) memcpy(a, b, n);
1543:     else {
1544:       memcpy(b + n, b + (n - (int)(a - b)), (int)(a - b));
1545:       PetscMemmove(a, b, n - (int)(a - b));
1546:     }
1547:   }
1548: #else
1549:   memmove((char *)(a), (char *)(b), n);
1550: #endif
1551:   return 0;
1552: }

1554: /*@C
1555:    PetscMemcpy - Copies n bytes, beginning at location b, to the space
1556:    beginning at location a. The two memory regions CANNOT overlap, use
1557:    `PetscMemmove()` in that case.

1559:    Not Collective

1561:    Input Parameters:
1562: +  b - pointer to initial memory space
1563: -  n - length (in bytes) of space to copy

1565:    Output Parameter:
1566: .  a - pointer to copy space

1568:    Level: intermediate

1570:    Compile Option:
1571:     `PETSC_PREFER_DCOPY_FOR_MEMCPY` will cause the BLAS dcopy() routine to be used
1572:                                   for memory copies on double precision values.
1573:     `PETSC_PREFER_COPY_FOR_MEMCPY` will cause C code to be used
1574:                                   for memory copies on double precision values.
1575:     `PETSC_PREFER_FORTRAN_FORMEMCPY` will cause Fortran code to be used
1576:                                   for memory copies on double precision values.

1578:    Notes:
1579:    Prefer `PetscArraycpy()`

1581:    This routine is analogous to `memcpy()`.

1583:    Not available from Fortran

1585:    Developer Note:
1586:    This is inlined for fastest performance

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

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

1632: /*@C
1633:    PetscMemzero - Zeros the specified memory.

1635:    Not Collective

1637:    Input Parameters:
1638: +  a - pointer to beginning memory location
1639: -  n - length (in bytes) of memory to initialize

1641:    Level: intermediate

1643:    Compile Option:
1644:    `PETSC_PREFER_BZERO` - on certain machines (the IBM RS6000) the bzero() routine happens
1645:   to be faster than the memset() routine. This flag causes the bzero() routine to be used.

1647:    Notes:
1648:    Not available from Fortran

1650:    Prefer `PetscArrayzero()`

1652:    Developer Note:
1653:    This is inlined for fastest performance

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

1687: /*MC
1688:    PetscArraycmp - Compares two arrays in memory.

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

1694:    Not Collective

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

1701:    Output Parameters:
1702: .   e - `PETSC_TRUE` if equal else `PETSC_FALSE`.

1704:    Level: intermediate

1706:    Notes:
1707:    This routine is a preferred replacement to `PetscMemcmp()`

1709:    The arrays must be of the same type

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

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

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

1724:    Not Collective

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

1731:    Level: intermediate

1733:    Notes:
1734:    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:    Notes:
1759:    This routine is a preferred replacement to `PetscMemcpy()`

1761:    The arrays must be of the same type

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

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

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

1774:    Not Collective

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

1780:    Level: intermediate

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

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

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

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

1796:    Not Collective

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

1804:    Level: developer

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

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

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

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

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

1835:   #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1836:     #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1837:   #endif

1839:   #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
1840:     #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
1841:   #endif

1843:   #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1844:     #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
1845:   #endif

1847:   #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1848:     #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
1849:   #endif

1851:   #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
1852:     #define PETSC_USE_FORTRAN_KERNEL_NORM
1853:   #endif

1855:   #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
1856:     #define PETSC_USE_FORTRAN_KERNEL_MAXPY
1857:   #endif

1859:   #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1860:     #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
1861:   #endif

1863:   #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1864:     #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
1865:   #endif

1867:   #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
1868:     #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
1869:   #endif

1871:   #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1872:     #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
1873:   #endif

1875:   #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
1876:     #define PETSC_USE_FORTRAN_KERNEL_MDOT
1877:   #endif

1879:   #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
1880:     #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
1881:   #endif

1883:   #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
1884:     #define PETSC_USE_FORTRAN_KERNEL_AYPX
1885:   #endif

1887:   #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
1888:     #define PETSC_USE_FORTRAN_KERNEL_WAXPY
1889:   #endif

1891: #endif

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

1908: /* --------------------------------------------------------------------*/

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

1914:    Level: beginner

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

1918: .seealso: `PETSC_COMM_WORLD`, `PETSC_COMM_SELF`
1919: M*/

1921: #if defined(PETSC_HAVE_MPIIO)
1922: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4);
1923: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4);
1924: PETSC_EXTERN PetscErrorCode MPIU_File_write_at(MPI_File, MPI_Offset, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(3, 5);
1925: PETSC_EXTERN PetscErrorCode MPIU_File_read_at(MPI_File, MPI_Offset, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(3, 5);
1926: PETSC_EXTERN PetscErrorCode MPIU_File_write_at_all(MPI_File, MPI_Offset, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(3, 5);
1927: PETSC_EXTERN PetscErrorCode MPIU_File_read_at_all(MPI_File, MPI_Offset, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(3, 5);
1928: #endif

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

1932: /* Limit MPI to 32-bits */
1933: #define PETSC_MPI_INT_MAX 2147483647
1934: #define PETSC_MPI_INT_MIN -2147483647
1935: /* Limit BLAS to 32-bits */
1936: #define PETSC_BLAS_INT_MAX    2147483647
1937: #define PETSC_BLAS_INT_MIN    -2147483647
1938: #define PETSC_CUBLAS_INT_MAX  2147483647
1939: #define PETSC_HIPBLAS_INT_MAX 2147483647

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

1945:    Not Collective

1947:    Input Parameter:
1948: .     a - the `PetscInt64` value

1950:    Output Parameter:
1951: .     b - the resulting `PetscInt` value

1953:    Level: advanced

1955:    Notes:
1956:    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

1958:    Not available from Fortran

1960: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscIntMultError()`, `PetscIntSumError()`
1961: @*/
1962: static inline PetscErrorCode PetscIntCast(PetscInt64 a, PetscInt *b)
1963: {
1964: #if !defined(PETSC_USE_64BIT_INDICES)
1965:   if (a > PETSC_MAX_INT) {
1966:     *b = 0;
1967:     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);
1968:   }
1969: #endif
1970:   *b = (PetscInt)(a);
1971:   return 0;
1972: }

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

1978:    Not Collective

1980:    Input Parameter:
1981: .     a - the `PetscCount` value

1983:    Output Parameter:
1984: .     b - the resulting `PetscInt` value

1986:    Level: advanced

1988:    Notes:
1989:      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

1991:    Not available from Fortran

1993: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscIntMultError()`, `PetscIntSumError()`, `PetscIntCast()`
1994: @*/
1995: static inline PetscErrorCode PetscCountCast(PetscCount a, PetscInt *b)
1996: {
1997:   if (sizeof(PetscCount) > sizeof(PetscInt) && a > PETSC_MAX_INT) {
1998:     *b = 0;
1999:     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);
2000:   }
2001:   *b = (PetscInt)(a);
2002:   return 0;
2003: }

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

2009:    Not Collective

2011:    Input Parameter:
2012: .     a - the `PetscInt` value

2014:    Output Parameter:
2015: .     b - the resulting `PetscBLASInt` value

2017:    Level: advanced

2019:    Notes:
2020:    Not available from Fortran

2022:    Errors if the integer is negative since PETSc calls to BLAS/LAPACK never need to cast negative integer inputs

2024: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscIntCast()`, `PetscCountCast()`
2025: @*/
2026: static inline PetscErrorCode PetscBLASIntCast(PetscInt a, PetscBLASInt *b)
2027: {
2028:   *b = 0;
2029:   if (PetscDefined(USE_64BIT_INDICES) && !PetscDefined(HAVE_64BIT_BLAS_INDICES)) {
2031:   }
2033:   *b = (PetscBLASInt)(a);
2034:   return 0;
2035: }

2037: /*@C
2038:     PetscCuBLASIntCast - like `PetscBLASIntCast()`, but for `PetscCuBLASInt`.

2040:    Not Collective

2042:    Input Parameter:
2043: .     a - the `PetscInt` value

2045:    Output Parameter:
2046: .     b - the resulting `PetscCuBLASInt` value

2048:    Level: advanced

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

2053: .seealso: `PetscCuBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscIntCast()`
2054: @*/
2055: static inline PetscErrorCode PetscCuBLASIntCast(PetscInt a, PetscCuBLASInt *b)
2056: {
2057:   *b = 0;
2060:   *b = (PetscCuBLASInt)(a);
2061:   return 0;
2062: }

2064: /*@C
2065:     PetscHipBLASIntCast - like `PetscBLASIntCast()`, but for `PetscHipBLASInt`.

2067:    Not Collective

2069:    Input Parameter:
2070: .     a - the `PetscInt` value

2072:    Output Parameter:
2073: .     b - the resulting `PetscHipBLASInt` value

2075:    Level: advanced

2077:    Notes:
2078:       Errors if the integer is negative since PETSc calls to hipBLAS and friends never need to cast negative integer inputs

2080: .seealso: `PetscHipBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscIntCast()`
2081: @*/
2082: static inline PetscErrorCode PetscHipBLASIntCast(PetscInt a, PetscHipBLASInt *b)
2083: {
2084:   *b = 0;
2087:   *b = (PetscHipBLASInt)(a);
2088:   return 0;
2089: }

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

2095:    Not Collective

2097:    Input Parameter:
2098: .     a - the `PetscInt` value

2100:    Output Parameter:
2101: .     b - the resulting `PetscMPIInt` value

2103:    Level: advanced

2105:    Not available from Fortran

2107: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscIntCast()`
2108: @*/
2109: static inline PetscErrorCode PetscMPIIntCast(PetscInt a, PetscMPIInt *b)
2110: {
2111:   *b = 0;
2113:   *b = (PetscMPIInt)(a);
2114:   return 0;
2115: }

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

2119: /*@C

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

2123:    Not Collective

2125:    Input Parameters:
2126: +     a - the `PetscReal` value
2127: -     b - the second value

2129:    Returns:
2130:       the result as a `PetscInt` value

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

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

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

2141:    Not available from Fortran

2143:    Level: advanced

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

2154: /*@C

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

2158:    Not Collective

2160:    Input Parameters:
2161: +     a - the PetscInt value
2162: -     b - the second value

2164:    Returns:
2165:       the result as a `PetscInt` value

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

2171:    Not available from Fortran

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

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

2178:    Level: advanced

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

2189: /*@C

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

2193:    Not Collective

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

2199:    Returns:
2200:      the result as a `PetscInt` value

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

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

2208:    Not available from Fortran

2210:    Level: advanced

2212: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`
2213: @*/
2214: static inline PetscInt PetscIntSumTruncate(PetscInt a, PetscInt b)
2215: {
2216:   PetscInt64 r = ((PetscInt64)a) + ((PetscInt64)b);
2217:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2218:   return (PetscInt)r;
2219: }

2221: /*@C

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

2225:    Not Collective

2227:    Input Parameters:
2228: +     a - the `PetscInt` value
2229: -     b - the second value

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

2234:    Use `PetscInt64Mult()` to compute the product of two `PetscInt` and store in a `PetscInt64`
2235:    Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt`

2237:    Not available from Fortran

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

2242:    Level: advanced

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

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

2258: /*@C

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

2262:    Not Collective

2264:    Input Parameters:
2265: +     a - the `PetscInt` value
2266: -     b - the second value

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

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

2274:    Not available from Fortran

2276:    Level: advanced

2278: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`
2279: @*/
2280: static inline PetscErrorCode PetscIntSumError(PetscInt a, PetscInt b, PetscInt *result)
2281: {
2282:   PetscInt64 r;

2284:   r = ((PetscInt64)a) + ((PetscInt64)b);
2285: #if !defined(PETSC_USE_64BIT_INDICES)
2287: #endif
2288:   if (result) *result = (PetscInt)r;
2289:   return 0;
2290: }

2292: /*
2293:      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2294: */
2295: #if defined(hz)
2296:   #undef hz
2297: #endif

2299: #include <limits.h>

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

2303: #define PETSC_BITS_PER_BYTE CHAR_BIT

2305: #if defined(PETSC_HAVE_SYS_TYPES_H)
2306:   #include <sys/types.h>
2307: #endif

2309: /*MC

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

2314:     The current PETSc version and the API for accessing it are defined in <A HREF="PETSC_DOC_OUT_.._PLACEHOLDER/include/petscversion.h.html">include/petscverson.html</A>

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

2318:     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
2319:     where only a change in the major version number (x) indicates a change in the API.

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

2323:     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),
2324:     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
2325:     version number (x.y.z).

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

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

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

2333:     `PETSC_VERSION_DATE_GIT` is the date of the last git commit to the repository

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

2337:     Level: intermediate

2339: M*/

2341: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[], size_t);
2342: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[], size_t);
2343: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[], size_t);
2344: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[], size_t);
2345: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2346: PETSC_EXTERN PetscErrorCode PetscGetDate(char[], size_t);
2347: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2348: PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt *, PetscInt *, PetscInt *, PetscInt *);

2350: PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt, const PetscInt[], PetscBool *);
2351: PETSC_EXTERN PetscErrorCode PetscSortedInt64(PetscInt, const PetscInt64[], PetscBool *);
2352: PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt, const PetscMPIInt[], PetscBool *);
2353: PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt, const PetscReal[], PetscBool *);
2354: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt, PetscInt[]);
2355: PETSC_EXTERN PetscErrorCode PetscSortInt64(PetscInt, PetscInt64[]);
2356: PETSC_EXTERN PetscErrorCode PetscSortCount(PetscInt, PetscCount[]);
2357: PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt, PetscInt[]);
2358: PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *, PetscInt[]);
2359: PETSC_EXTERN PetscErrorCode PetscSortedCheckDupsInt(PetscInt, const PetscInt[], PetscBool *);
2360: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt *, PetscInt[]);
2362: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt *);
2363: PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt *);
2364: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt, const PetscInt[], PetscInt[]);
2365: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt, const char *[], PetscInt[]);
2366: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt, PetscInt[], PetscInt[]);
2367: PETSC_EXTERN PetscErrorCode PetscSortIntWithCountArray(PetscCount, PetscInt[], PetscCount[]);
2368: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt, PetscInt[], PetscInt[], PetscInt[]);
2369: PETSC_EXTERN PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount, PetscInt[], PetscInt[], PetscCount[]);
2370: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt, PetscMPIInt[]);
2371: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *, PetscMPIInt[]);
2372: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt, PetscMPIInt[], PetscMPIInt[]);
2373: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt, PetscMPIInt[], PetscInt[]);
2374: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt, PetscInt[], PetscScalar[]);
2375: PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt, PetscInt[], void *, size_t, void *);
2376: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt, PetscReal[]);
2377: PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt, PetscReal[], PetscInt[]);
2378: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt, const PetscReal[], PetscInt[]);
2379: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt *, PetscReal[]);
2380: PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal, PetscInt, const PetscReal[], PetscReal, PetscInt *);
2381: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt, PetscInt, PetscScalar[], PetscInt[]);
2382: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt, PetscInt, PetscReal[], PetscInt[]);
2383: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt, const PetscBool[], const PetscInt[], PetscInt *, PetscInt **, PetscInt **, PetscInt **, PetscInt **);
2384: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt, const PetscInt[], const PetscInt[], PetscInt, const PetscInt[], const PetscInt[], PetscInt *, PetscInt **, PetscInt **);
2385: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscInt *, PetscInt **);
2386: PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt, const PetscMPIInt[], PetscInt, const PetscMPIInt[], PetscInt *, PetscMPIInt **);
2387: PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *);

2389: PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt, void *, size_t, int (*)(const void *, const void *, void *), void *);
2390: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt, PetscInt[]);
2391: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt, PetscMPIInt[]);
2392: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt, PetscReal[]);
2393: PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt, void *, size_t, void *, size_t, int (*)(const void *, const void *, void *), void *);
2394: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt, PetscInt[], PetscInt[]);
2395: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt, PetscMPIInt[], PetscMPIInt[]);
2396: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt, PetscReal[], PetscInt[]);

2398: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2399: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[], size_t);

2401: /*J
2402:     PetscRandomType - String with the name of a PETSc randomizer

2404:    Level: beginner

2406:    Notes:
2407:    To use `PETSCSPRNG` or `PETSCRANDOM123` you must have ./configure PETSc
2408:    with the option --download-sprng or --download-random123

2410: .seealso: `PetscRandomSetType()`, `PetscRandom`, `PetscRandomCreate()`
2411: J*/
2412: typedef const char *PetscRandomType;
2413: #define PETSCRAND      "rand"
2414: #define PETSCRAND48    "rand48"
2415: #define PETSCSPRNG     "sprng"
2416: #define PETSCRANDER48  "rander48"
2417: #define PETSCRANDOM123 "random123"
2418: #define PETSCCURAND    "curand"

2420: /* Logging support */
2421: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;

2423: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);

2425: /* Dynamic creation and loading functions */
2426: PETSC_EXTERN PetscFunctionList PetscRandomList;

2428: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[], PetscErrorCode (*)(PetscRandom));
2429: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2430: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2431: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType *);
2432: PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom, PetscObject, const char[]);
2433: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom, PetscViewer);

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

2447: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[], char[], size_t);
2448: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[], char[], size_t);
2449: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[], size_t);
2450: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[], char[]);
2451: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[], size_t);
2452: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[], char, PetscBool *);
2453: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[], char, PetscBool *);
2454: PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2455: PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2456: PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);

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

2464: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int, void *, PetscInt, PetscInt *, PetscDataType);
2465: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm, int, void *, PetscInt, PetscInt *, PetscDataType);
2466: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int, const void *, PetscInt, PetscDataType);
2467: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm, int, const void *, PetscInt, PetscDataType);
2468: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[], PetscFileMode, int *);
2469: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2470: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm, PetscBool *);
2471: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm, PetscBool *);
2472: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm, char[], size_t);
2473: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm, const char[], char[], size_t, PetscBool *);
2474: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm, const char[], char[], size_t, PetscBool *);
2475: #if defined(PETSC_USE_SOCKET_VIEWER)
2476: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[], int, int *);
2477: #endif

2479: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int, off_t, PetscBinarySeekType, off_t *);
2480: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm, int, off_t, PetscBinarySeekType, off_t *);
2481: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *, PetscDataType, PetscInt);

2483: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2484: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[], PetscBool);
2485: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2486: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char *);
2487: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2488: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2489: PETSC_EXTERN PetscErrorCode PetscWaitOnError(void);

2491: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt *);
2492: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], PetscMPIInt **, PetscMPIInt **);
2493: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt **, PetscMPIInt **, PetscMPIInt **);
2494: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscInt ***, MPI_Request **);
2495: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscScalar ***, MPI_Request **);
2496: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm, PetscMPIInt, MPI_Datatype, PetscMPIInt, const PetscMPIInt *, const void *, PetscMPIInt *, PetscMPIInt **, void *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(6, 3);
2497: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm, PetscMPIInt, MPI_Datatype, PetscMPIInt, const PetscMPIInt[], const void *, PetscMPIInt *, PetscMPIInt **, void *, PetscMPIInt, PetscErrorCode (*send)(MPI_Comm, const PetscMPIInt[], PetscMPIInt, PetscMPIInt, void *, MPI_Request[], void *), PetscErrorCode (*recv)(MPI_Comm, const PetscMPIInt[], PetscMPIInt, void *, MPI_Request[], void *), void *ctx) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(6, 3);
2498: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm, PetscMPIInt, MPI_Datatype, PetscMPIInt, const PetscMPIInt[], const void *, PetscMPIInt *, PetscMPIInt **, void *, PetscMPIInt, MPI_Request **, MPI_Request **, PetscErrorCode (*send)(MPI_Comm, const PetscMPIInt[], PetscMPIInt, PetscMPIInt, void *, MPI_Request[], void *), PetscErrorCode (*recv)(MPI_Comm, const PetscMPIInt[], PetscMPIInt, void *, MPI_Request[], void *), void *ctx) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(6, 3);

2500: PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2501: PETSC_EXTERN PetscErrorCode    PetscCommBuildTwoSidedSetType(MPI_Comm, PetscBuildTwoSidedType);
2502: PETSC_EXTERN PetscErrorCode    PetscCommBuildTwoSidedGetType(MPI_Comm, PetscBuildTwoSidedType *);

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

2506: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);

2508: PETSC_EXTERN const char *const PetscSubcommTypes[];

2510: struct _n_PetscSubcomm {
2511:   MPI_Comm         parent;    /* parent communicator */
2512:   MPI_Comm         dupparent; /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2513:   MPI_Comm         child;     /* the sub-communicator */
2514:   PetscMPIInt      n;         /* num of subcommunicators under the parent communicator */
2515:   PetscMPIInt      color;     /* color of processors belong to this communicator */
2516:   PetscMPIInt     *subsize;   /* size of subcommunicator[color] */
2517:   PetscSubcommType type;
2518:   char            *subcommprefix;
2519: };

2521: static inline MPI_Comm PetscSubcommParent(PetscSubcomm scomm)
2522: {
2523:   return scomm->parent;
2524: }
2525: static inline MPI_Comm PetscSubcommChild(PetscSubcomm scomm)
2526: {
2527:   return scomm->child;
2528: }
2529: static inline MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm)
2530: {
2531:   return scomm->dupparent;
2532: }
2533: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm, PetscSubcomm *);
2534: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm *);
2535: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm, PetscInt);
2536: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm, PetscSubcommType);
2537: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm, PetscMPIInt, PetscMPIInt);
2538: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm, PetscViewer);
2539: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2540: PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm, const char[]);
2541: PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm, MPI_Comm *);
2542: PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm, MPI_Comm *);
2543: PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm, MPI_Comm *);

2545: PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt, PetscHeap *);
2546: PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap, PetscInt, PetscInt);
2547: PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap, PetscInt *, PetscInt *);
2548: PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap, PetscInt *, PetscInt *);
2549: PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap, PetscInt, PetscInt);
2550: PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2551: PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap *);
2552: PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap, PetscViewer);

2554: PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2555: PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm, PetscShmComm *);
2556: PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm, PetscMPIInt, PetscMPIInt *);
2557: PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm, PetscMPIInt, PetscMPIInt *);
2558: PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm, MPI_Comm *);

2560: /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2561: PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm, PetscInt, PetscOmpCtrl *);
2562: PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl, MPI_Comm *, MPI_Comm *, PetscBool *);
2563: PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl *);
2564: PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2565: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2566: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);

2568: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t, size_t, PetscSegBuffer *);
2569: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer *);
2570: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer, size_t, void *);
2571: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer, void *);
2572: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer, void *);
2573: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer, void *);
2574: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer, size_t *);
2575: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer, size_t);

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

2585: extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2586: PETSC_EXTERN PetscErrorCode    PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted *);
2587: PETSC_EXTERN PetscErrorCode    PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted *);
2588: PETSC_EXTERN PetscErrorCode    PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted, const char *, const char *, PetscBool *);

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

2594: PETSC_EXTERN PetscSegBuffer PetscCitationsList;

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

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

2601:      Input Parameters:
2602: +      cite - the bibtex item, formatted to displayed on multiple lines nicely
2603: -      set - a boolean variable initially set to `PETSC_FALSE`; this is used to insure only a single registration of the citation

2605:      Options Database: Key
2606: .     -citations [filename]   - print out the bibtex entries for the given computation

2608:      Level: intermediate

2610:      Fortran Note:
2611:      Not available from Fortran

2613: @*/
2614: static inline PetscErrorCode PetscCitationsRegister(const char cit[], PetscBool *set)
2615: {
2616:   size_t len;
2617:   char  *vstring;

2619:   if (set && *set) return 0;
2620:   PetscStrlen(cit, &len);
2621:   PetscSegBufferGet(PetscCitationsList, len, &vstring);
2622:   PetscArraycpy(vstring, cit, len);
2623:   if (set) *set = PETSC_TRUE;
2624:   return 0;
2625: }

2627: PETSC_EXTERN                PETSC_DEPRECATED_FUNCTION("Google has discontinued its URL shortener service") PetscErrorCode PetscURLShorten(const char[], char[], size_t c);
2628: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm, char[], char[], size_t);
2629: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm, const char[], char[], size_t);
2630: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm, const char[], const char[]);

2632: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm, char[], char[], size_t);
2633: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm, const char[], char[], char[], size_t);

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

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

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

2643: #if defined(PETSC_USE_DEBUG)
2644: static inline unsigned int PetscStrHash(const char *str)
2645: {
2646:   unsigned int c, hash = 5381;

2648:   while ((c = (unsigned int)*str++)) hash = ((hash << 5) + hash) + c; /* hash * 33 + c */
2649:   return hash;
2650: }

2652:   /*MC
2653:    MPIU_Allreduce - a PETSc replacement for `MPI_Allreduce()` that tries to determine if the call from all the MPI ranks occur from the
2654:                     same place in the PETSc code. This helps to detect bugs where different MPI ranks follow different code paths
2655:                     resulting in inconsistent and incorrect calls to `MPI_Allreduce()`.

2657:    Collective

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

2663:    Input Parameters:
2664: +  a - pointer to the input data to be reduced
2665: .  c - the number of MPI data items in a and b
2666: .  d - the MPI datatype, for example `MPI_INT`
2667: .  e - the MPI operation, for example `MPI_SUM`
2668: -  fcomm - the MPI communicator on which the operation occurs

2670:    Output Parameter:
2671: .  b - the reduced values

2673:    Notes:
2674:      In optimized mode this directly calls `MPI_Allreduce()`

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

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

2680:    Level: developer

2682: .seealso: `MPI_Allreduce()`
2683: M*/
2684:   #define MPIU_Allreduce(a, b, c, d, e, fcomm) \
2686: #else
2687:   #define MPIU_Allreduce(a, b, c, d, e, fcomm) PetscMacroReturnStandard(PetscCallMPI(MPI_Allreduce((a), (b), (c), d, e, (fcomm))))
2688: #endif

2690: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
2691: PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint, PetscMPIInt, MPI_Info, MPI_Comm, void *, MPI_Win *);
2692: PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win, PetscMPIInt, MPI_Aint *, PetscMPIInt *, void *);
2693: #endif

2695: /* this is a vile hack */
2696: #if defined(PETSC_HAVE_NECMPI)
2697:   #if !defined(PETSC_NECMPI_VERSION_MAJOR) || !defined(PETSC_NECMPI_VERSION_MINOR) || PETSC_NECMPI_VERSION_MAJOR < 2 || (PETSC_NECMPI_VERSION_MAJOR == 2 && PETSC_NECMPI_VERSION_MINOR < 18)
2698:     #define MPI_Type_free(a) (*(a) = MPI_DATATYPE_NULL, 0);
2699:   #endif
2700: #endif

2702: /*
2703:     List of external packages and queries on it
2704: */
2705: PETSC_EXTERN PetscErrorCode PetscHasExternalPackage(const char[], PetscBool *);

2707: /*
2708:  OpenMP support
2709: */
2710: #if defined(_OPENMP)
2711:   #define PetscPragmaOMP(...) _Pragma(PetscStringize(omp __VA_ARGS__))
2712: #else // no OpenMP so no threads
2713:   #define PetscPragmaOMP(...)
2714: #endif

2716: /* this cannot go here because it may be in a different shared library */
2717: PETSC_EXTERN PetscErrorCode PCMPIServerBegin(void);
2718: PETSC_EXTERN PetscErrorCode PCMPIServerEnd(void);
2719: PETSC_EXTERN PetscErrorCode PCMPICommsDestroy(void);
2720: #endif