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: #pragma once

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

 22: /* SUBMANSEC = Sys */

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

 43: #include <petscsystypes.h>

 45: /* ========================================================================== */

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

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

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

112: /*
113:     Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler
114:     see the top of mpicxx.h in the MPICH2 distribution.
115: */
116: #include <stdio.h>

118: /* MSMPI on 32-bit Microsoft Windows requires this yukky hack - that breaks MPI standard compliance */
119: #if !defined(MPIAPI)
120:   #define MPIAPI
121: #endif

123: PETSC_EXTERN MPI_Datatype MPIU_ENUM PETSC_ATTRIBUTE_MPI_TYPE_TAG(PetscEnum);
124: PETSC_EXTERN MPI_Datatype MPIU_BOOL PETSC_ATTRIBUTE_MPI_TYPE_TAG(PetscBool);

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

129:    Level: beginner

131:    Note:
132:    In MPI calls that require an MPI datatype that matches a `PetscInt` or array of `PetscInt` values, pass this value.

134: .seealso: `PetscReal`, `PetscScalar`, `PetscComplex`, `PetscInt`, `MPIU_COUNT`, `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`
135: M*/

137: PETSC_EXTERN MPI_Datatype MPIU_FORTRANADDR;

139: #if defined(PETSC_USE_64BIT_INDICES)
140:   #define MPIU_INT MPIU_INT64
141: #else
142:   #define MPIU_INT MPI_INT
143: #endif

145: /*MC
146:    MPIU_COUNT - Portable MPI datatype corresponding to `PetscCount` independent of the precision of `PetscCount`

148:    Level: beginner

150:    Note:
151:    In MPI calls that require an MPI datatype that matches a `PetscCount` or array of `PetscCount` values, pass this value.

153:   Developer Note:
154:   It seems `MPI_AINT` is unsigned so this may be the wrong choice here since `PetscCount` is signed

156: .seealso: `PetscReal`, `PetscScalar`, `PetscComplex`, `PetscInt`, `MPIU_INT`, `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`
157: M*/
158: #define MPIU_COUNT MPI_AINT

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

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

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

177: /*
178:   Handle inclusion when using clang compiler with CUDA support
179:   __float128 is not available for the device
180: */
181: #if defined(__clang__) && defined(__CUDA_ARCH__)
182:   #define PETSC_SKIP_REAL___FLOAT128
183: #endif

185: /*
186:     Declare extern C stuff after including external header files
187: */

189: PETSC_EXTERN PetscBool PETSC_RUNNING_ON_VALGRIND;
190: /*
191:     Defines elementary mathematics functions and constants.
192: */
193: #include <petscmath.h>

195: /*MC
196:     PETSC_IGNORE - same as `NULL`, means PETSc will ignore this argument

198:    Level: beginner

200:    Note:
201:    Accepted by many PETSc functions to not set a parameter and instead use a default value

203:    Fortran Note:
204:    Use `PETSC_NULL_INTEGER`, `PETSC_NULL_DOUBLE_PRECISION` etc

206: .seealso: `PETSC_DECIDE`, `PETSC_DEFAULT`, `PETSC_DETERMINE`
207: M*/
208: #define PETSC_IGNORE PETSC_NULLPTR
209: #define PETSC_NULL   PETSC_DEPRECATED_MACRO(3, 19, 0, "PETSC_NULLPTR", ) PETSC_NULLPTR

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

215:    Level: beginner

217: .seealso: `PETSC_DEFAULT`, `PETSC_IGNORE`, `PETSC_DETERMINE`
218: M*/

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

223:    Level: beginner

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

229: .seealso: `PETSC_DECIDE`, `PETSC_DEFAULT`, `PETSC_IGNORE`, `VecSetSizes()`
230: M*/

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

235:    Level: beginner

237:    Fortran Note:
238:    You need to use `PETSC_DEFAULT_INTEGER` or `PETSC_DEFAULT_REAL`.

240: .seealso: `PETSC_DECIDE`, `PETSC_IGNORE`, `PETSC_DETERMINE`
241: M*/

243: /* These MUST be preprocessor defines! see https://gitlab.com/petsc/petsc/-/issues/1370 */
244: #define PETSC_DECIDE    (-1)
245: #define PETSC_DETERMINE PETSC_DECIDE
246: #define PETSC_DEFAULT   (-2)

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

251:    Level: beginner

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

259:    The value of `PETSC_COMM_WORLD` should never be used or accessed before `PetscInitialize()`
260:    is called because it may not have a valid value yet.

262: .seealso: `PETSC_COMM_SELF`
263: M*/
264: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;

266: /*MC
267:    PETSC_COMM_SELF - This is always `MPI_COMM_SELF`

269:    Level: beginner

271:    Note:
272:    Do not USE/access or set this variable before `PetscInitialize()` has been called.

274: .seealso: `PETSC_COMM_WORLD`
275: M*/
276: #define PETSC_COMM_SELF MPI_COMM_SELF

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

281:    No Fortran Support

283:    Level: beginner

285:    Note:
286:    By default `PETSC_MPI_THREAD_REQUIRED` equals `MPI_THREAD_FUNNELED` when the MPI implementation provides `MPI_Init_thread()`, otherwise it equals `MPI_THREAD_SINGLE`

288: .seealso: `PetscInitialize()`
289: M*/
290: PETSC_EXTERN PetscMPIInt PETSC_MPI_THREAD_REQUIRED;

292: /*MC
293:    PetscBeganMPI - indicates if PETSc initialized MPI during `PetscInitialize()` or if MPI was already initialized.

295:    Synopsis:
296: #include <petscsys.h>
297:    PetscBool PetscBeganMPI;

299:    No Fortran Support

301:    Level: developer

303: .seealso: `PetscInitialize()`, `PetscInitializeCalled`
304: M*/
305: PETSC_EXTERN PetscBool PetscBeganMPI;

307: PETSC_EXTERN PetscBool PetscErrorHandlingInitialized;
308: PETSC_EXTERN PetscBool PetscInitializeCalled;
309: PETSC_EXTERN PetscBool PetscFinalizeCalled;
310: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;

312: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm), PetscErrorCode (*)(MPI_Comm));
313: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm, MPI_Comm *, int *);
314: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm *);
315: PETSC_EXTERN PetscErrorCode PetscCommGetComm(MPI_Comm, MPI_Comm *);
316: PETSC_EXTERN PetscErrorCode PetscCommRestoreComm(MPI_Comm, MPI_Comm *);

318: #if defined(PETSC_HAVE_KOKKOS)
319: PETSC_EXTERN PetscErrorCode PetscKokkosInitializeCheck(void); /* Initialize Kokkos if not yet. */
320: #endif

322: #if defined(PETSC_HAVE_NVSHMEM)
323: PETSC_EXTERN PetscBool      PetscBeganNvshmem;
324: PETSC_EXTERN PetscBool      PetscNvshmemInitialized;
325: PETSC_EXTERN PetscErrorCode PetscNvshmemFinalize(void);
326: #endif

328: #if defined(PETSC_HAVE_ELEMENTAL)
329: PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
330: PETSC_EXTERN PetscErrorCode PetscElementalInitialized(PetscBool *);
331: PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
332: #endif

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

337:    Synopsis:
338: #include <petscsys.h>
339:    PetscErrorCode PetscMalloc(size_t m,void **result)

341:    Not Collective

343:    Input Parameter:
344: .  m - number of bytes to allocate

346:    Output Parameter:
347: .  result - memory allocated

349:    Level: beginner

351:    Notes:
352:    Memory is always allocated at least double aligned

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

356: .seealso: `PetscFree()`, `PetscNew()`
357: M*/
358: #define PetscMalloc(a, b) ((*PetscTrMalloc)((a), PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (void **)(b)))

360: /*MC
361:    PetscRealloc - Reallocates memory

363:    Synopsis:
364: #include <petscsys.h>
365:    PetscErrorCode PetscRealloc(size_t m,void **result)

367:    Not Collective

369:    Input Parameters:
370: +  m - number of bytes to allocate
371: -  result - previous memory

373:    Output Parameter:
374: .  result - new memory allocated

376:    Level: developer

378:    Note:
379:    Memory is always allocated at least double aligned

381: .seealso: `PetscMalloc()`, `PetscFree()`, `PetscNew()`
382: M*/
383: #define PetscRealloc(a, b) ((*PetscTrRealloc)((a), __LINE__, PETSC_FUNCTION_NAME, __FILE__, (void **)(b)))

385: /*MC
386:    PetscAddrAlign - Rounds up an address to `PETSC_MEMALIGN` alignment

388:    Synopsis:
389: #include <petscsys.h>
390:    void *PetscAddrAlign(void *addr)

392:    Not Collective

394:    Input Parameter:
395: .  addr - address to align (any pointer type)

397:    Level: developer

399: .seealso: `PetscMallocAlign()`
400: M*/
401: #define PetscAddrAlign(a) ((void *)((((PETSC_UINTPTR_T)(a)) + (PETSC_MEMALIGN - 1)) & ~(PETSC_MEMALIGN - 1)))

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

406:    Synopsis:
407: #include <petscsys.h>
408:    PetscErrorCode PetscCalloc(size_t m,void **result)

410:    Not Collective

412:    Input Parameter:
413: .  m - number of bytes to allocate

415:    Output Parameter:
416: .  result - memory allocated

418:    Level: beginner

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

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

425: .seealso: `PetscFree()`, `PetscNew()`
426: M*/
427: #define PetscCalloc(m, result) PetscMallocA(1, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)m), (result))

429: /*MC
430:    PetscMalloc1 - Allocates an array of memory aligned to `PETSC_MEMALIGN`

432:    Synopsis:
433: #include <petscsys.h>
434:    PetscErrorCode PetscMalloc1(size_t m1,type **r1)

436:    Not Collective

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

441:    Output Parameter:
442: .  r1 - memory allocated

444:    Level: beginner

446:    Note:
447:    This uses the sizeof() of the memory type requested to determine the total memory to be allocated, therefore you should not
448:          multiply the number of elements requested by the `sizeof()` the type. For example use
449: .vb
450:   PetscInt *id;
451:   PetscMalloc1(10,&id);
452: .ve
453:        not
454: .vb
455:   PetscInt *id;
456:   PetscMalloc1(10*sizeof(PetscInt),&id);
457: .ve

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

461: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc2()`
462: M*/
463: #define PetscMalloc1(m1, r1) PetscMallocA(1, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1))

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

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

472:    Not Collective

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

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

480:    Level: beginner

482:    Note:
483:    See `PetsMalloc1()` for more details on usage.

485: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc2()`
486: M*/
487: #define PetscCalloc1(m1, r1) PetscMallocA(1, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1))

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

492:    Synopsis:
493: #include <petscsys.h>
494:    PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)

496:    Not Collective

498:    Input Parameters:
499: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
500: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

502:    Output Parameters:
503: +  r1 - memory allocated in first chunk
504: -  r2 - memory allocated in second chunk

506:    Level: developer

508: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc2()`
509: M*/
510: #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))

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

515:    Synopsis:
516: #include <petscsys.h>
517:    PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)

519:    Not Collective

521:    Input Parameters:
522: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
523: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

525:    Output Parameters:
526: +  r1 - memory allocated in first chunk
527: -  r2 - memory allocated in second chunk

529:    Level: developer

531: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc2()`
532: M*/
533: #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))

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

538:    Synopsis:
539: #include <petscsys.h>
540:    PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

542:    Not Collective

544:    Input Parameters:
545: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
546: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
547: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

549:    Output Parameters:
550: +  r1 - memory allocated in first chunk
551: .  r2 - memory allocated in second chunk
552: -  r3 - memory allocated in third chunk

554:    Level: developer

556: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc3()`, `PetscFree3()`
557: M*/
558: #define PetscMalloc3(m1, r1, m2, r2, m3, r3) \
559:   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))

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

564:    Synopsis:
565: #include <petscsys.h>
566:    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

568:    Not Collective

570:    Input Parameters:
571: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
572: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
573: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

575:    Output Parameters:
576: +  r1 - memory allocated in first chunk
577: .  r2 - memory allocated in second chunk
578: -  r3 - memory allocated in third chunk

580:    Level: developer

582: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc2()`, `PetscMalloc3()`, `PetscFree3()`
583: M*/
584: #define PetscCalloc3(m1, r1, m2, r2, m3, r3) \
585:   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))

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

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

594:    Not Collective

596:    Input Parameters:
597: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
598: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
599: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
600: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

602:    Output Parameters:
603: +  r1 - memory allocated in first chunk
604: .  r2 - memory allocated in second chunk
605: .  r3 - memory allocated in third chunk
606: -  r4 - memory allocated in fourth chunk

608:    Level: developer

610: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc4()`, `PetscFree4()`
611: M*/
612: #define PetscMalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \
613:   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))

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

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

622:    Not Collective

624:    Input Parameters:
625: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
626: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
627: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
628: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

630:    Output Parameters:
631: +  r1 - memory allocated in first chunk
632: .  r2 - memory allocated in second chunk
633: .  r3 - memory allocated in third chunk
634: -  r4 - memory allocated in fourth chunk

636:    Level: developer

638: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc4()`, `PetscFree4()`
639: M*/
640: #define PetscCalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \
641:   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))

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

646:    Synopsis:
647: #include <petscsys.h>
648:    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)

650:    Not Collective

652:    Input Parameters:
653: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
654: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
655: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
656: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
657: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

659:    Output Parameters:
660: +  r1 - memory allocated in first chunk
661: .  r2 - memory allocated in second chunk
662: .  r3 - memory allocated in third chunk
663: .  r4 - memory allocated in fourth chunk
664: -  r5 - memory allocated in fifth chunk

666:    Level: developer

668: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc5()`, `PetscFree5()`
669: M*/
670: #define PetscMalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \
671:   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))

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

676:    Synopsis:
677: #include <petscsys.h>
678:    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)

680:    Not Collective

682:    Input Parameters:
683: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
684: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
685: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
686: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
687: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

689:    Output Parameters:
690: +  r1 - memory allocated in first chunk
691: .  r2 - memory allocated in second chunk
692: .  r3 - memory allocated in third chunk
693: .  r4 - memory allocated in fourth chunk
694: -  r5 - memory allocated in fifth chunk

696:    Level: developer

698: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc5()`, `PetscFree5()`
699: M*/
700: #define PetscCalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \
701:   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))

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

706:    Synopsis:
707: #include <petscsys.h>
708:    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)

710:    Not Collective

712:    Input Parameters:
713: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
714: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
715: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
716: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
717: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
718: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

720:    Output Parameteasr:
721: +  r1 - memory allocated in first chunk
722: .  r2 - memory allocated in second chunk
723: .  r3 - memory allocated in third chunk
724: .  r4 - memory allocated in fourth chunk
725: .  r5 - memory allocated in fifth chunk
726: -  r6 - memory allocated in sixth chunk

728:    Level: developer

730: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc6()`, `PetscFree3()`, `PetscFree4()`, `PetscFree5()`, `PetscFree6()`
731: M*/
732: #define PetscMalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \
733:   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))

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

738:    Synopsis:
739: #include <petscsys.h>
740:    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)

742:    Not Collective

744:    Input Parameters:
745: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
746: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
747: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
748: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
749: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
750: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

752:    Output Parameters:
753: +  r1 - memory allocated in first chunk
754: .  r2 - memory allocated in second chunk
755: .  r3 - memory allocated in third chunk
756: .  r4 - memory allocated in fourth chunk
757: .  r5 - memory allocated in fifth chunk
758: -  r6 - memory allocated in sixth chunk

760:    Level: developer

762: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscMalloc6()`, `PetscFree6()`
763: M*/
764: #define PetscCalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \
765:   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))

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

770:    Synopsis:
771: #include <petscsys.h>
772:    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)

774:    Not Collective

776:    Input Parameters:
777: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
778: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
779: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
780: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
781: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
782: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
783: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

785:    Output Parameters:
786: +  r1 - memory allocated in first chunk
787: .  r2 - memory allocated in second chunk
788: .  r3 - memory allocated in third chunk
789: .  r4 - memory allocated in fourth chunk
790: .  r5 - memory allocated in fifth chunk
791: .  r6 - memory allocated in sixth chunk
792: -  r7 - memory allocated in seventh chunk

794:    Level: developer

796: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc7()`, `PetscFree7()`
797: M*/
798: #define PetscMalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \
799:   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))

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

804:    Synopsis:
805: #include <petscsys.h>
806:    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)

808:    Not Collective

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

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

828:    Level: developer

830: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscMalloc7()`, `PetscFree7()`
831: M*/
832: #define PetscCalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \
833:   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))

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

838:    Synopsis:
839: #include <petscsys.h>
840:    PetscErrorCode PetscNew(type **result)

842:    Not Collective

844:    Output Parameter:
845: .  result - memory allocated, sized to match pointer type

847:    Level: beginner

849: .seealso: `PetscFree()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc1()`
850: M*/
851: #define PetscNew(b) PetscCalloc1(1, (b))

853: #define PetscNewLog(o, b) PETSC_DEPRECATED_MACRO(3, 18, 0, "PetscNew()", ) PetscNew(b)

855: /*MC
856:    PetscFree - Frees memory

858:    Synopsis:
859: #include <petscsys.h>
860:    PetscErrorCode PetscFree(void *memory)

862:    Not Collective

864:    Input Parameter:
865: .   memory - memory to free (the pointer is ALWAYS set to `NULL` upon success)

867:    Level: beginner

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

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

874: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc1()`
875: M*/
876: #define PetscFree(a) ((PetscErrorCode)((*PetscTrFree)((void *)(a), __LINE__, PETSC_FUNCTION_NAME, __FILE__) || ((a) = PETSC_NULLPTR, PETSC_SUCCESS)))

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

881:    Synopsis:
882: #include <petscsys.h>
883:    PetscErrorCode PetscFree2(void *memory1,void *memory2)

885:    Not Collective

887:    Input Parameters:
888: +   memory1 - memory to free
889: -   memory2 - 2nd memory to free

891:    Level: developer

893:    Notes:
894:     Memory must have been obtained with `PetscMalloc2()`

896:     The arguments need to be in the same order as they were in the call to `PetscMalloc2()`

898: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`
899: M*/
900: #define PetscFree2(m1, m2) PetscFreeA(2, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2))

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

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

909:    Not Collective

911:    Input Parameters:
912: +   memory1 - memory to free
913: .   memory2 - 2nd memory to free
914: -   memory3 - 3rd memory to free

916:    Level: developer

918:    Notes:
919:     Memory must have been obtained with `PetscMalloc3()`

921:     The arguments need to be in the same order as they were in the call to `PetscMalloc3()`

923: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`
924: M*/
925: #define PetscFree3(m1, m2, m3) PetscFreeA(3, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3))

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

930:    Synopsis:
931: #include <petscsys.h>
932:    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)

934:    Not Collective

936:    Input Parameters:
937: +   m1 - memory to free
938: .   m2 - 2nd memory to free
939: .   m3 - 3rd memory to free
940: -   m4 - 4th memory to free

942:    Level: developer

944:    Notes:
945:     Memory must have been obtained with `PetscMalloc4()`

947:     The arguments need to be in the same order as they were in the call to `PetscMalloc4()`

949: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`
950: M*/
951: #define PetscFree4(m1, m2, m3, m4) PetscFreeA(4, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4))

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

956:    Synopsis:
957: #include <petscsys.h>
958:    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)

960:    Not Collective

962:    Input Parameters:
963: +   m1 - memory to free
964: .   m2 - 2nd memory to free
965: .   m3 - 3rd memory to free
966: .   m4 - 4th memory to free
967: -   m5 - 5th memory to free

969:    Level: developer

971:    Notes:
972:     Memory must have been obtained with `PetscMalloc5()`

974:     The arguments need to be in the same order as they were in the call to `PetscMalloc5()`

976: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`
977: M*/
978: #define PetscFree5(m1, m2, m3, m4, m5) PetscFreeA(5, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5))

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

983:    Synopsis:
984: #include <petscsys.h>
985:    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)

987:    Not Collective

989:    Input Parameters:
990: +   m1 - memory to free
991: .   m2 - 2nd memory to free
992: .   m3 - 3rd memory to free
993: .   m4 - 4th memory to free
994: .   m5 - 5th memory to free
995: -   m6 - 6th memory to free

997:    Level: developer

999:    Notes:
1000:     Memory must have been obtained with `PetscMalloc6()`

1002:     The arguments need to be in the same order as they were in the call to `PetscMalloc6()`

1004: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()`
1005: M*/
1006: #define PetscFree6(m1, m2, m3, m4, m5, m6) PetscFreeA(6, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6))

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

1011:    Synopsis:
1012: #include <petscsys.h>
1013:    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)

1015:    Not Collective

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

1026:    Level: developer

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

1031:     The arguments need to be in the same order as they were in the call to `PetscMalloc7()`

1033: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()`,
1034:           `PetscMalloc7()`
1035: M*/
1036: #define PetscFree7(m1, m2, m3, m4, m5, m6, m7) PetscFreeA(7, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6), &(m7))

1038: PETSC_EXTERN PetscErrorCode PetscMallocA(int, PetscBool, int, const char *, const char *, size_t, void *, ...);
1039: PETSC_EXTERN PetscErrorCode PetscFreeA(int, int, const char *, const char *, void *, ...);
1040: PETSC_EXTERN                PetscErrorCode (*PetscTrMalloc)(size_t, PetscBool, int, const char[], const char[], void **);
1041: PETSC_EXTERN                PetscErrorCode (*PetscTrFree)(void *, int, const char[], const char[]);
1042: PETSC_EXTERN                PetscErrorCode (*PetscTrRealloc)(size_t, int, const char[], const char[], void **);
1043: PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1044: 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 **));
1045: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);

1047: /*
1048:   Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1049: */
1050: PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1051: PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);
1052: #if defined(PETSC_HAVE_CUDA)
1053: PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void);
1054: PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void);
1055: #endif
1056: #if defined(PETSC_HAVE_HIP)
1057: PETSC_EXTERN PetscErrorCode PetscMallocSetHIPHost(void);
1058: PETSC_EXTERN PetscErrorCode PetscMallocResetHIPHost(void);
1059: #endif

1061: #define MPIU_PETSCLOGDOUBLE  MPI_DOUBLE
1062: #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION

1064: /*
1065:    Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1066: */
1067: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1068: PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1069: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1070: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1071: PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1072: PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int, PetscLogDouble *);
1073: PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool, PetscBool);
1074: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool *, PetscBool *, PetscBool *);
1075: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int, const char[], const char[]);
1076: PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1077: PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool *);
1078: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool);
1079: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool *);

1081: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType, MPI_Datatype *);
1082: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype, PetscDataType *);
1083: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType, size_t *);
1084: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char *, PetscDataType *, PetscBool *);

1086: /*
1087:    These are MPI operations for MPI_Allreduce() etc
1088: */
1089: PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1090: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1091: PETSC_EXTERN MPI_Op MPIU_SUM;
1092: PETSC_EXTERN MPI_Op MPIU_MAX;
1093: PETSC_EXTERN MPI_Op MPIU_MIN;
1094: #else
1095:   #define MPIU_SUM MPI_SUM
1096:   #define MPIU_MAX MPI_MAX
1097:   #define MPIU_MIN MPI_MIN
1098: #endif
1099: PETSC_EXTERN MPI_Op         Petsc_Garbage_SetIntersectOp;
1100: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm, const PetscInt[], PetscInt *, PetscInt *);

1102: #if (defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128)) || (defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16))
1103: /*MC
1104:    MPIU_SUM___FP16___FLOAT128 - MPI_Op that acts as a replacement for `MPI_SUM` with
1105:    custom `MPI_Datatype` `MPIU___FLOAT128`, `MPIU___COMPLEX128`, and `MPIU___FP16`.

1107:    Level: advanced

1109:    Developer Note:
1110:    This should be unified with `MPIU_SUM`

1112: .seealso: `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`
1113: M*/
1114: PETSC_EXTERN MPI_Op MPIU_SUM___FP16___FLOAT128;
1115: #endif
1116: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm, const PetscInt[], PetscInt *, PetscInt *);

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

1121: /*
1122:     Defines PETSc error handling.
1123: */
1124: #include <petscerror.h>

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

1131: #define PETSC_SMALLEST_CLASSID ((PetscClassId)1211211)
1132: PETSC_EXTERN PetscClassId   PETSC_LARGEST_CLASSID;
1133: PETSC_EXTERN PetscClassId   PETSC_OBJECT_CLASSID;
1134: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[], PetscClassId *);
1135: PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject, PetscObjectId *);
1136: PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject, PetscObjectId, PetscBool *);

1138: /*
1139:    Routines that get memory usage information from the OS
1140: */
1141: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1142: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1143: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1144: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);

1146: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);

1148: /*
1149:    Initialization of PETSc
1150: */
1151: PETSC_EXTERN PetscErrorCode PetscInitialize(int *, char ***, const char[], const char[]);
1152: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int, char **, const char[], const char[]);
1153: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1154: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1155: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1156: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1157: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1158: PETSC_EXTERN PetscErrorCode PetscGetArgs(int *, char ***);
1159: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1160: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);

1162: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1163: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);
1164: PETSC_EXTERN PetscErrorCode PetscSysFinalizePackage(void);

1166: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[], const char[]);
1167: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1168: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1169: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject, const char[]);

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

1173: /*
1174:      These are so that in extern C code we can caste function pointers to non-extern C
1175:    function pointers. Since the regular C++ code expects its function pointers to be C++
1176: */

1178: /*S
1179:   PetscVoidFn - A prototype of a void (fn)(void) function

1181:   Level: developer

1183:   Notes:
1184:   The deprecated `PetscVoidFunction` works as a replacement for `PetscVoidFn` *.

1186:   The deprecated `PetscVoidStarFunction` works as a replacement for `PetscVoidFn` **.

1188: .seealso: `PetscObject`, `PetscObjectDestroy()`
1189: S*/
1190: PETSC_EXTERN_TYPEDEF typedef void(PetscVoidFn)(void);

1192: PETSC_EXTERN_TYPEDEF typedef PetscVoidFn  *PetscVoidFunction;
1193: PETSC_EXTERN_TYPEDEF typedef PetscVoidFn **PetscVoidStarFunction;

1195: /*S
1196:   PetscErrorCodeFn - A prototype of a PetscErrorCode (fn)(void) function

1198:   Level: developer

1200:   Notes:
1201:   The deprecated `PetscErrorCodeFunction` works as a replacement for `PetscErrorCodeFn` *.

1203: .seealso: `PetscObject`, `PetscObjectDestroy()`
1204: S*/
1205: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(PetscErrorCodeFn)(void);

1207: PETSC_EXTERN_TYPEDEF typedef PetscErrorCodeFn *PetscErrorCodeFunction;

1209: /*
1210:     Functions that can act on any PETSc object.
1211: */
1212: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject *);
1213: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject, MPI_Comm *);
1214: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject, PetscClassId *);
1215: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject, const char *[]);
1216: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject, const char *[]);
1217: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject, const char[]);
1218: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject, const char *[]);
1219: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject, PetscInt);
1220: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject, PetscInt *);
1221: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject, PetscObject, PetscInt);
1222: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1223: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject, PetscInt *);
1224: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1225: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject, PetscMPIInt *);
1226: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject, const char[], PetscObject);
1227: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject, const char[]);
1228: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject, const char[], PetscObject *);
1229: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject, const char[], void (*)(void));
1230: #define PetscObjectComposeFunction(a, b, ...) PetscObjectComposeFunction_Private((a), (b), (PetscVoidFn *)(__VA_ARGS__))
1231: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1232: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1233: PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1234: PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject, PetscObject);
1235: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm, PetscMPIInt *);

1237: #include <petscviewertypes.h>
1238: #include <petscoptions.h>

1240: PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer, PetscBool, PetscLogDouble);
1241: PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool *);

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

1245: PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer, const char[]);
1246: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject, PetscViewer);
1247: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject, PetscViewer);
1248: #define PetscObjectQueryFunction(obj, name, fptr) PetscObjectQueryFunction_Private((obj), (name), (PetscVoidFn **)(fptr))
1249: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject, const char[], void (**)(void));
1250: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject, const char[]);
1251: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject, const char[]);
1252: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject, const char[]);
1253: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject, const char *[]);
1254: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject, const char[]);
1255: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1256: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1257: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject, PetscObject, const char[]);
1258: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1259: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject, const char[], PetscBool *);
1260: PETSC_EXTERN PetscErrorCode PetscObjectObjectTypeCompare(PetscObject, PetscObject, PetscBool *);
1261: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject, const char[], PetscBool *);
1262: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject, PetscBool *, const char[], ...);
1263: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject, PetscBool *, const char[], ...);
1264: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1265: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);

1267: #if defined(PETSC_HAVE_SAWS)
1268: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1269: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1270: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject, PetscBool);
1271: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1272: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1273: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1274: PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1275: PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1276: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1277: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);

1279: #else
1280:   #define PetscSAWsBlock()                  PETSC_SUCCESS
1281:   #define PetscObjectSAWsViewOff(obj)       PETSC_SUCCESS
1282:   #define PetscObjectSAWsSetBlock(obj, flg) PETSC_SUCCESS
1283:   #define PetscObjectSAWsBlock(obj)         PETSC_SUCCESS
1284:   #define PetscObjectSAWsGrantAccess(obj)   PETSC_SUCCESS
1285:   #define PetscObjectSAWsTakeAccess(obj)    PETSC_SUCCESS
1286:   #define PetscStackViewSAWs()              PETSC_SUCCESS
1287:   #define PetscStackSAWsViewOff()           PETSC_SUCCESS
1288:   #define PetscStackSAWsTakeAccess()
1289:   #define PetscStackSAWsGrantAccess()

1291: #endif

1293: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[], PetscDLMode, PetscDLHandle *);
1294: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1295: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle, const char[], void **);
1296: PETSC_EXTERN PetscErrorCode PetscDLAddr(void (*)(void), char **);
1297: #ifdef PETSC_HAVE_CXX
1298: PETSC_EXTERN PetscErrorCode PetscDemangleSymbol(const char *, char **);
1299: #endif

1301: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void *, PetscStack **);

1303: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE *, PetscBool);
1304: PETSC_EXTERN PetscErrorCode PetscObjectsView(PetscViewer);
1305: PETSC_EXTERN PetscErrorCode PetscObjectsGetObject(const char *, PetscObject *, char **);
1306: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList *);
1307: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList, const char[], PetscObject *);
1308: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList, PetscObject, char **, PetscBool *);
1309: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *, const char[], PetscObject);
1310: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *, const char[]);
1311: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList, PetscObjectList *);

1313: /*
1314:     Dynamic library lists. Lists of names of routines in objects or in dynamic
1315:   link libraries that will be loaded as needed.
1316: */

1318: #define PetscFunctionListAdd(list, name, fptr) PetscFunctionListAdd_Private((list), (name), (PetscVoidFn *)(fptr))
1319: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList *, const char[], PetscVoidFn *);
1320: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList *);
1321: PETSC_EXTERN PetscErrorCode PetscFunctionListClear(PetscFunctionList);
1322: #define PetscFunctionListFind(list, name, fptr) PetscFunctionListFind_Private((list), (name), (PetscVoidFn **)(fptr))
1323: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList, const char[], PetscVoidFn **);
1324: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm, FILE *, const char[], const char[], const char[], const char[], PetscFunctionList, const char[], const char[]);
1325: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList, PetscFunctionList *);
1326: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList, PetscViewer);
1327: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList, const char ***, int *);
1328: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintNonEmpty(PetscFunctionList);
1329: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintAll(void);

1331: PETSC_EXTERN PetscDLLibrary PetscDLLibrariesLoaded;
1332: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm, PetscDLLibrary *, const char[]);
1333: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm, PetscDLLibrary *, const char[]);
1334: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm, PetscDLLibrary *, const char[], const char[], void **);
1335: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1336: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm, const char[], char *, size_t, PetscBool *);
1337: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm, const char[], PetscDLLibrary *);
1338: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);

1340: /*
1341:      Useful utility routines
1342: */
1343: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm, PetscInt *, PetscInt *);
1344: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm, PetscInt, PetscInt *, PetscInt *);
1345: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm, PetscInt *, PetscInt *);
1346: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm, PetscMPIInt);
1347: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm, PetscMPIInt);
1348: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1349: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE *);
1350: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm, const PetscInt[2], PetscInt[2]);
1351: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm, const PetscReal[2], PetscReal[2]);

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

1356:     Level: beginner

1358:     Note:
1359:     This is useful in cases like
1360: .vb
1361:      int        *a;
1362:      PetscBool  flag = PetscNot(a)
1363: .ve
1364:      where !a would not return a `PetscBool` because we cannot provide a cast from int to `PetscBool` in C.

1366: .seealso: `PetscBool`, `PETSC_TRUE`, `PETSC_FALSE`
1367: M*/
1368: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)

1370: /*MC
1371:    PetscHelpPrintf - Prints help messages.

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

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

1379:    Input Parameters:
1380: +  comm - the MPI communicator over which the help message is printed
1381: .  format - the usual printf() format string
1382: -  args - arguments to be printed

1384:    Level: developer

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

1389:    To use, write your own function, for example,
1390: .vb
1391:    PetscErrorCode mypetschelpprintf(MPI_Comm comm,const char format[],....)
1392:    {
1393:      PetscFunctionReturn(PETSC_SUCCESS);
1394:    }
1395: .ve
1396: then do the assignment
1397: .vb
1398:   PetscHelpPrintf = mypetschelpprintf;
1399: .ve

1401:   You can do the assignment before `PetscInitialize()`.

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

1405: .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscErrorPrintf()`, `PetscHelpPrintfDefault()`
1406: M*/
1407: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);

1409: /*
1410:      Defines PETSc profiling.
1411: */
1412: #include <petsclog.h>

1414: /*
1415:       Simple PETSc parallel IO for ASCII printing
1416: */
1417: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[], char[]);
1418: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm, const char[], const char[], FILE **);
1419: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm, FILE *);
1420: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1421: PETSC_EXTERN PetscErrorCode PetscFFlush(FILE *);
1422: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1423: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char *, size_t, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1424: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char *, size_t, const char[], size_t *, ...) PETSC_ATTRIBUTE_FORMAT(3, 5);
1425: PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[], size_t, const char *, PetscInt, const PetscReal[]);

1427: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1428: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1429: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);

1431: PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char *, size_t *);
1432: PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char *, char *);

1434: #if defined(PETSC_HAVE_POPEN)
1435: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm, const char[], const char[], const char[], FILE **);
1436: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm, FILE *);
1437: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1438: #endif

1440: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1441: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1442: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm, FILE *);
1443: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm, FILE *, size_t, char[]);
1444: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm, const char[], const char[], FILE **);
1445: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char *[]);

1447: PETSC_EXTERN PetscClassId   PETSC_CONTAINER_CLASSID;
1448: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer, void **);
1449: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer, void *);
1450: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer *);
1451: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm, PetscContainer *);
1452: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void *));
1453: PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void *);

1455: /*
1456:    For use in debuggers
1457: */
1458: PETSC_EXTERN PetscMPIInt    PetscGlobalRank;
1459: PETSC_EXTERN PetscMPIInt    PetscGlobalSize;
1460: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt, const PetscInt[], PetscViewer);
1461: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt, const PetscReal[], PetscViewer);
1462: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt, const PetscScalar[], PetscViewer);

1464: /*
1465:     Basic memory and string operations. These are usually simple wrappers
1466:    around the basic Unix system calls, but a few of them have additional
1467:    functionality and/or error checking.
1468: */
1469: #include <petscstring.h>

1471: #include <stddef.h>
1472: #include <stdlib.h>

1474: #if defined(PETSC_CLANG_STATIC_ANALYZER)
1475:   #define PetscPrefetchBlock(a, b, c, d)
1476: #else
1477:   /*MC
1478:    PetscPrefetchBlock - Prefetches a block of memory

1480:    Synopsis:
1481: #include <petscsys.h>
1482:     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)

1484:    Not Collective

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

1492:    Level: developer

1494:    Notes:
1495:    The last two arguments (`rw` and `t`) must be compile-time constants.

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

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

1507: M*/
1508:   #define PetscPrefetchBlock(a, n, rw, t) \
1509:     do { \
1510:       const char *_p = (const char *)(a), *_end = (const char *)((a) + (n)); \
1511:       for (; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p, (rw), (t)); \
1512:     } while (0)
1513: #endif
1514: /*
1515:       Determine if some of the kernel computation routines use
1516:    Fortran (rather than C) for the numerical calculations. On some machines
1517:    and compilers (like complex numbers) the Fortran version of the routines
1518:    is faster than the C/C++ versions. The flag --with-fortran-kernels
1519:    should be used with ./configure to turn these on.
1520: */
1521: #if defined(PETSC_USE_FORTRAN_KERNELS)

1523:   #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1524:     #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1525:   #endif

1527:   #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1528:     #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
1529:   #endif

1531:   #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1532:     #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
1533:   #endif

1535:   #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
1536:     #define PETSC_USE_FORTRAN_KERNEL_MAXPY
1537:   #endif

1539:   #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1540:     #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
1541:   #endif

1543:   #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
1544:     #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
1545:   #endif

1547:   #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1548:     #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
1549:   #endif

1551:   #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
1552:     #define PETSC_USE_FORTRAN_KERNEL_MDOT
1553:   #endif

1555:   #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
1556:     #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
1557:   #endif

1559:   #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
1560:     #define PETSC_USE_FORTRAN_KERNEL_AYPX
1561:   #endif

1563:   #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
1564:     #define PETSC_USE_FORTRAN_KERNEL_WAXPY
1565:   #endif

1567: #endif

1569: /*
1570:     Macros for indicating code that should be compiled with a C interface,
1571:    rather than a C++ interface. Any routines that are dynamically loaded
1572:    (such as the PCCreate_XXX() routines) must be wrapped so that the name
1573:    mangler does not change the functions symbol name. This just hides the
1574:    ugly extern "C" {} wrappers.
1575: */
1576: #if defined(__cplusplus)
1577:   #define EXTERN_C_BEGIN extern "C" {
1578:   #define EXTERN_C_END   }
1579: #else
1580:   #define EXTERN_C_BEGIN
1581:   #define EXTERN_C_END
1582: #endif

1584: /*MC
1585:    MPI_Comm - the basic object used by MPI to determine which processes are involved in a
1586:    communication

1588:    Level: beginner

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

1593: .seealso: `PETSC_COMM_WORLD`, `PETSC_COMM_SELF`
1594: M*/

1596: #if defined(PETSC_HAVE_MPIIO)
1597: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4);
1598: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4);
1599: 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);
1600: 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);
1601: 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);
1602: 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);
1603: #endif

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

1609:    Not Collective; No Fortran Support

1611:    Input Parameter:
1612: .  a - the `PetscInt64` value

1614:    Output Parameter:
1615: .  b - the resulting `PetscInt` value

1617:    Level: advanced

1619:    Note:
1620:    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 integers

1622: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscIntMultError()`, `PetscIntSumError()`
1623: @*/
1624: static inline PetscErrorCode PetscIntCast(PetscInt64 a, PetscInt *b)
1625: {
1626:   PetscFunctionBegin;
1627:   *b = 0;
1628:   // if using 64-bit indices already then this comparison is tautologically true
1629:   PetscCheck(a < PETSC_MAX_INT, 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);
1630:   *b = (PetscInt)a;
1631:   PetscFunctionReturn(PETSC_SUCCESS);
1632: }

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

1638:    Not Collective; No Fortran Support

1640:    Input Parameter:
1641: .  a - the `PetscCount` value

1643:    Output Parameter:
1644: .  b - the resulting `PetscInt` value

1646:    Level: advanced

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

1651: .seealso: `PetscCount`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscIntMultError()`, `PetscIntSumError()`, `PetscIntCast()`
1652: @*/
1653: static inline PetscErrorCode PetscCountCast(PetscCount a, PetscInt *b)
1654: {
1655:   PetscFunctionBegin;
1656:   *b = 0;
1657:   PetscCheck(sizeof(PetscCount) <= sizeof(PetscInt) || a <= PETSC_MAX_INT, 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);
1658:   *b = (PetscInt)a;
1659:   PetscFunctionReturn(PETSC_SUCCESS);
1660: }

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

1666:    Not Collective; No Fortran Support

1668:    Input Parameter:
1669: .  a - the `PetscInt` value

1671:    Output Parameter:
1672: .  b - the resulting `PetscBLASInt` value

1674:    Level: advanced

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

1679: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscIntCast()`, `PetscCountCast()`
1680: @*/
1681: static inline PetscErrorCode PetscBLASIntCast(PetscInt a, PetscBLASInt *b)
1682: {
1683:   PetscFunctionBegin;
1684:   *b = 0;
1685:   if (PetscDefined(USE_64BIT_INDICES) && !PetscDefined(HAVE_64BIT_BLAS_INDICES)) {
1686:     PetscCheck(a <= PETSC_BLAS_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " is too big for BLAS/LAPACK, which is restricted to 32-bit integers. Either you have an invalidly large integer error in your code or you must ./configure PETSc with --with-64-bit-blas-indices for the case you are running", a);
1687:   }
1688:   PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer to BLAS/LAPACK routine");
1689:   *b = (PetscBLASInt)a;
1690:   PetscFunctionReturn(PETSC_SUCCESS);
1691: }

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

1696:    Not Collective; No Fortran Support

1698:    Input Parameter:
1699: .  a - the `PetscInt` value

1701:    Output Parameter:
1702: .  b - the resulting `PetscCuBLASInt` value

1704:    Level: advanced

1706:    Note:
1707:    Errors if the integer is negative since PETSc calls to cuBLAS and friends never need to cast negative integer inputs

1709: .seealso: `PetscCuBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscIntCast()`
1710: @*/
1711: static inline PetscErrorCode PetscCuBLASIntCast(PetscInt a, PetscCuBLASInt *b)
1712: {
1713:   PetscFunctionBegin;
1714:   *b = 0;
1715:   PetscCheck(a <= PETSC_CUBLAS_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " is too big for cuBLAS, which is restricted to 32-bit integers.", a);
1716:   PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer %" PetscInt_FMT "to cuBLAS routine", a);
1717:   *b = (PetscCuBLASInt)a;
1718:   PetscFunctionReturn(PETSC_SUCCESS);
1719: }

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

1724:    Not Collective; No Fortran Support

1726:    Input Parameter:
1727: .  a - the `PetscInt` value

1729:    Output Parameter:
1730: .  b - the resulting `PetscHipBLASInt` value

1732:    Level: advanced

1734:    Note:
1735:    Errors if the integer is negative since PETSc calls to hipBLAS and friends never need to cast negative integer inputs

1737: .seealso: `PetscHipBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscIntCast()`
1738: @*/
1739: static inline PetscErrorCode PetscHipBLASIntCast(PetscInt a, PetscHipBLASInt *b)
1740: {
1741:   PetscFunctionBegin;
1742:   *b = 0;
1743:   PetscCheck(a <= PETSC_HIPBLAS_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " is too big for hipBLAS, which is restricted to 32-bit integers.", a);
1744:   PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer %" PetscInt_FMT "to hipBLAS routine", a);
1745:   *b = (PetscHipBLASInt)a;
1746:   PetscFunctionReturn(PETSC_SUCCESS);
1747: }

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

1753:    Not Collective; No Fortran Support

1755:    Input Parameter:
1756: .  a - the `PetscInt` value

1758:    Output Parameter:
1759: .  b - the resulting `PetscMPIInt` value

1761:    Level: advanced

1763: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscIntCast()`
1764: @*/
1765: static inline PetscErrorCode PetscMPIIntCast(PetscInt a, PetscMPIInt *b)
1766: {
1767:   PetscFunctionBegin;
1768:   *b = 0;
1769:   PetscCheck(a <= PETSC_MPI_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " is too big for MPI buffer length. Maximum supported value is %d", a, PETSC_MPI_INT_MAX);
1770:   *b = (PetscMPIInt)a;
1771:   PetscFunctionReturn(PETSC_SUCCESS);
1772: }

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

1776: /*@C
1777:   PetscRealIntMultTruncate - Computes the product of a positive `PetscReal` and a positive
1778:   `PetscInt` and truncates the value to slightly less than the maximal possible value.

1780:   Not Collective; No Fortran Support

1782:   Input Parameters:
1783: + a - The `PetscReal` value
1784: - b - The `PetscInt` value

1786:   Level: advanced

1788:   Notes:
1789:   Returns the result as a `PetscInt` value.

1791:   Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`.

1793:   Use `PetscIntMultTruncate()` to compute the product of two positive `PetscInt` and truncate
1794:   to fit a `PetscInt`.

1796:   Use `PetscIntMultError()` to compute the product of two `PetscInt` if you wish to generate an
1797:   error if the result will not fit in a `PetscInt`.

1799:   Developer Notes:
1800:   We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but
1801:   requires many more checks.

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

1806: .seealso: `PetscReal`, `PetscInt`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()`
1807: @*/
1808: static inline PetscInt PetscRealIntMultTruncate(PetscReal a, PetscInt b)
1809: {
1810:   PetscInt64 r = (PetscInt64)(a * (PetscReal)b);
1811:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
1812:   return (PetscInt)r;
1813: }

1815: /*@C
1816:    PetscIntMultTruncate - Computes the product of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value

1818:    Not Collective; No Fortran Support

1820:    Input Parameters:
1821: +  a - the `PetscInt` value
1822: -  b - the second value

1824:    Returns:
1825:    The result as a `PetscInt` value

1827:    Level: advanced

1829:    Notes:
1830:    Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`

1832:    Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt`

1834:    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`

1836:    Developer Notes:
1837:    We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but requires many more checks.

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

1841: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()`,
1842:           `PetscIntSumTruncate()`
1843: @*/
1844: static inline PetscInt PetscIntMultTruncate(PetscInt a, PetscInt b)
1845: {
1846:   PetscInt64 r = PetscInt64Mult(a, b);
1847:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
1848:   return (PetscInt)r;
1849: }

1851: /*@C
1852:    PetscIntSumTruncate - Computes the sum of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value

1854:    Not Collective; No Fortran Support

1856:    Input Parameters:
1857: +  a - the `PetscInt` value
1858: -  b - the second value

1860:    Returns:
1861:    The result as a `PetscInt` value

1863:    Level: advanced

1865:    Notes:
1866:    Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`

1868:    Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt`

1870:    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`

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

1875: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`
1876: @*/
1877: static inline PetscInt PetscIntSumTruncate(PetscInt a, PetscInt b)
1878: {
1879:   PetscInt64 r = ((PetscInt64)a) + ((PetscInt64)b);
1880:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
1881:   return (PetscInt)r;
1882: }

1884: /*@C
1885:    PetscIntMultError - Computes the product of two positive `PetscInt` and generates an error with overflow.

1887:    Not Collective; No Fortran Support

1889:    Input Parameters:
1890: +  a - the `PetscInt` value
1891: -  b - the second value

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

1896:    Level: advanced

1898:    Notes:
1899:    Use `PetscInt64Mult()` to compute the product of two `PetscInt` and store in a `PetscInt64`

1901:    Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt`

1903:    Developer Note:
1904:    In most places in the source code we currently assume that `PetscInt` addition does not overflow, this is obviously wrong but requires many more checks.
1905:    `PetscIntSumError()` can be used to check for this situation.

1907: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscIntMult64()`, `PetscIntSumError()`
1908: @*/
1909: static inline PetscErrorCode PetscIntMultError(PetscInt a, PetscInt b, PetscInt *result)
1910: {
1911:   PetscInt64 r = PetscInt64Mult(a, b);

1913:   PetscFunctionBegin;
1914:   if (result) *result = (PetscInt)r;
1915:   if (!PetscDefined(USE_64BIT_INDICES)) {
1916:     PetscCheck(r <= PETSC_MAX_INT, PETSC_COMM_SELF, PETSC_ERR_SUP, "Product of two integers %" PetscInt_FMT " %" PetscInt_FMT " overflow, either you have an invalidly large integer error in your code or you must ./configure PETSc with --with-64-bit-indices for the case you are running", a, b);
1917:   }
1918:   PetscFunctionReturn(PETSC_SUCCESS);
1919: }

1921: /*@C

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

1925:    Not Collective; No Fortran Support

1927:    Input Parameters:
1928: +  a - the `PetscInt` value
1929: -  b - the second value

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

1934:    Level: advanced

1936:    Notes:
1937:    Use `PetscInt64Mult()` to compute the product of two 32-bit `PetscInt` and store in a `PetscInt64`

1939:    Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt`

1941: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`
1942: @*/
1943: static inline PetscErrorCode PetscIntSumError(PetscInt a, PetscInt b, PetscInt *result)
1944: {
1945:   PetscInt64 r = ((PetscInt64)a) + ((PetscInt64)b);

1947:   PetscFunctionBegin;
1948:   if (result) *result = (PetscInt)r;
1949:   if (!PetscDefined(USE_64BIT_INDICES)) {
1950:     PetscCheck(r <= PETSC_MAX_INT, PETSC_COMM_SELF, PETSC_ERR_SUP, "Sum of two integers %" PetscInt_FMT " %" PetscInt_FMT " overflow, either you have an invalidly large integer error in your code or you must ./configure PETSc with --with-64-bit-indices for the case you are running", a, b);
1951:   }
1952:   PetscFunctionReturn(PETSC_SUCCESS);
1953: }

1955: /*
1956:      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
1957: */
1958: #if defined(hz)
1959:   #undef hz
1960: #endif

1962: #if defined(PETSC_HAVE_SYS_TYPES_H)
1963:   #include <sys/types.h>
1964: #endif

1966: /*MC

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

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

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

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

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

1980:     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),
1981:     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
1982:     version number (x.y.z).

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

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

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

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

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

1994:     Level: intermediate
1995: M*/

1997: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[], size_t);
1998: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[], size_t);
1999: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[], size_t);
2000: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[], size_t);
2001: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2002: PETSC_EXTERN PetscErrorCode PetscGetDate(char[], size_t);
2003: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2004: PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt *, PetscInt *, PetscInt *, PetscInt *);

2006: PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt, const PetscInt[], PetscBool *);
2007: PETSC_EXTERN PetscErrorCode PetscSortedInt64(PetscInt, const PetscInt64[], PetscBool *);
2008: PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt, const PetscMPIInt[], PetscBool *);
2009: PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt, const PetscReal[], PetscBool *);
2010: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt, PetscInt[]);
2011: PETSC_EXTERN PetscErrorCode PetscSortInt64(PetscInt, PetscInt64[]);
2012: PETSC_EXTERN PetscErrorCode PetscSortCount(PetscInt, PetscCount[]);
2013: PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt, PetscInt[]);
2014: PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *, PetscInt[]);
2015: PETSC_EXTERN PetscErrorCode PetscSortedCheckDupsInt(PetscInt, const PetscInt[], PetscBool *);
2016: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt *, PetscInt[]);
2017: PETSC_EXTERN PetscErrorCode PetscCheckDupsInt(PetscInt, const PetscInt[], PetscBool *);
2018: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt *);
2019: PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt *);
2020: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt, const PetscInt[], PetscInt[]);
2021: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt, const char *[], PetscInt[]);
2022: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt, PetscInt[], PetscInt[]);
2023: PETSC_EXTERN PetscErrorCode PetscSortIntWithCountArray(PetscCount, PetscInt[], PetscCount[]);
2024: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt, PetscInt[], PetscInt[], PetscInt[]);
2025: PETSC_EXTERN PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount, PetscInt[], PetscInt[], PetscCount[]);
2026: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt, PetscMPIInt[]);
2027: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *, PetscMPIInt[]);
2028: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt, PetscMPIInt[], PetscMPIInt[]);
2029: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt, PetscMPIInt[], PetscInt[]);
2030: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt, PetscInt[], PetscScalar[]);
2031: PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt, PetscInt[], void *, size_t, void *);
2032: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt, PetscReal[]);
2033: PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt, PetscReal[], PetscInt[]);
2034: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt, const PetscReal[], PetscInt[]);
2035: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt *, PetscReal[]);
2036: PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal, PetscInt, const PetscReal[], PetscReal, PetscInt *);
2037: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt, PetscInt, PetscScalar[], PetscInt[]);
2038: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt, PetscInt, PetscReal[], PetscInt[]);
2039: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt, const PetscBool[], const PetscInt[], PetscInt *, PetscInt **, PetscInt **, PetscInt **, PetscInt **);
2040: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt, const PetscInt[], const PetscInt[], PetscInt, const PetscInt[], const PetscInt[], PetscInt *, PetscInt **, PetscInt **);
2041: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscInt *, PetscInt **);
2042: PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt, const PetscMPIInt[], PetscInt, const PetscMPIInt[], PetscInt *, PetscMPIInt **);
2043: PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *);

2045: PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt, void *, size_t, int (*)(const void *, const void *, void *), void *);
2046: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt, PetscInt[]);
2047: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt, PetscMPIInt[]);
2048: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt, PetscReal[]);
2049: PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt, void *, size_t, void *, size_t, int (*)(const void *, const void *, void *), void *);
2050: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt, PetscInt[], PetscInt[]);
2051: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt, PetscMPIInt[], PetscMPIInt[]);
2052: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt, PetscReal[], PetscInt[]);

2054: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2055: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[], size_t);

2057: /*J
2058:     PetscRandomType - String with the name of a PETSc randomizer

2060:    Level: beginner

2062:    Note:
2063:    To use `PETSCSPRNG` or `PETSCRANDOM123` you must have ./configure PETSc
2064:    with the option `--download-sprng` or `--download-random123`. We recommend the default provided with PETSc.

2066: .seealso: `PetscRandomSetType()`, `PetscRandom`, `PetscRandomCreate()`
2067: J*/
2068: typedef const char *PetscRandomType;
2069: #define PETSCRAND      "rand"
2070: #define PETSCRAND48    "rand48"
2071: #define PETSCSPRNG     "sprng"
2072: #define PETSCRANDER48  "rander48"
2073: #define PETSCRANDOM123 "random123"
2074: #define PETSCCURAND    "curand"

2076: /* Logging support */
2077: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;

2079: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);
2080: PETSC_EXTERN PetscErrorCode PetscRandomFinalizePackage(void);

2082: /* Dynamic creation and loading functions */
2083: PETSC_EXTERN PetscFunctionList PetscRandomList;

2085: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[], PetscErrorCode (*)(PetscRandom));
2086: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2087: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2088: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType *);
2089: PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom, PetscObject, const char[]);
2090: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom, PetscViewer);

2092: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm, PetscRandom *);
2093: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom, PetscScalar *);
2094: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom, PetscReal *);
2095: PETSC_EXTERN PetscErrorCode PetscRandomGetValues(PetscRandom, PetscInt, PetscScalar *);
2096: PETSC_EXTERN PetscErrorCode PetscRandomGetValuesReal(PetscRandom, PetscInt, PetscReal *);
2097: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom, PetscScalar *, PetscScalar *);
2098: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom, PetscScalar, PetscScalar);
2099: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom, unsigned long);
2100: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom, unsigned long *);
2101: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2102: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom *);

2104: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[], char[], size_t);
2105: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[], char[], size_t);
2106: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[], size_t);
2107: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[], char[]);
2108: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[], size_t);
2109: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[], char, PetscBool *);
2110: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[], char, PetscBool *);
2111: PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2112: PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2113: PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);

2115: /*MC
2116:    PetscBinaryBigEndian - indicates if values in memory are stored with big endian format

2118:    Synopsis:
2119: #include <petscsys.h>
2120:    PetscBool PetscBinaryBigEndian(void);

2122:    No Fortran Support

2124:    Level: developer

2126: .seealso: `PetscInitialize()`, `PetscFinalize()`, `PetscInitializeCalled`
2127: M*/
2128: static inline PetscBool PetscBinaryBigEndian(void)
2129: {
2130:   long _petsc_v = 1;
2131:   return ((char *)&_petsc_v)[0] ? PETSC_FALSE : PETSC_TRUE;
2132: }

2134: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int, void *, PetscInt, PetscInt *, PetscDataType);
2135: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm, int, void *, PetscInt, PetscInt *, PetscDataType);
2136: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int, const void *, PetscInt, PetscDataType);
2137: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm, int, const void *, PetscInt, PetscDataType);
2138: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[], PetscFileMode, int *);
2139: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2140: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm, PetscBool *);
2141: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm, PetscBool *);
2142: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm, char[], size_t);
2143: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm, const char[], char[], size_t, PetscBool *);
2144: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm, const char[], char[], size_t, PetscBool *);
2145: #if defined(PETSC_USE_SOCKET_VIEWER)
2146: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[], int, int *);
2147: #endif

2149: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int, off_t, PetscBinarySeekType, off_t *);
2150: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm, int, off_t, PetscBinarySeekType, off_t *);
2151: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *, PetscDataType, PetscInt);

2153: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2154: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[], PetscBool);
2155: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2156: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char *);
2157: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2158: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2159: PETSC_EXTERN PetscErrorCode PetscWaitOnError(void);

2161: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt *);
2162: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], PetscMPIInt **, PetscMPIInt **);
2163: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt **, PetscMPIInt **, PetscMPIInt **);
2164: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscInt ***, MPI_Request **);
2165: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscScalar ***, MPI_Request **);
2166: 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);
2167: 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);
2168: 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);

2170: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm, PetscBuildTwoSidedType);
2171: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm, PetscBuildTwoSidedType *);

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

2175: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);

2177: struct _n_PetscSubcomm {
2178:   MPI_Comm         parent;    /* parent communicator */
2179:   MPI_Comm         dupparent; /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2180:   MPI_Comm         child;     /* the sub-communicator */
2181:   PetscMPIInt      n;         /* num of subcommunicators under the parent communicator */
2182:   PetscMPIInt      color;     /* color of processors belong to this communicator */
2183:   PetscMPIInt     *subsize;   /* size of subcommunicator[color] */
2184:   PetscSubcommType type;
2185:   char            *subcommprefix;
2186: };

2188: static inline MPI_Comm PetscSubcommParent(PetscSubcomm scomm)
2189: {
2190:   return scomm->parent;
2191: }
2192: static inline MPI_Comm PetscSubcommChild(PetscSubcomm scomm)
2193: {
2194:   return scomm->child;
2195: }
2196: static inline MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm)
2197: {
2198:   return scomm->dupparent;
2199: }
2200: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm, PetscSubcomm *);
2201: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm *);
2202: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm, PetscInt);
2203: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm, PetscSubcommType);
2204: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm, PetscMPIInt, PetscMPIInt);
2205: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm, PetscViewer);
2206: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2207: PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm, const char[]);
2208: PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm, MPI_Comm *);
2209: PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm, MPI_Comm *);
2210: PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm, MPI_Comm *);

2212: PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt, PetscHeap *);
2213: PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap, PetscInt, PetscInt);
2214: PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap, PetscInt *, PetscInt *);
2215: PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap, PetscInt *, PetscInt *);
2216: PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap, PetscInt, PetscInt);
2217: PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2218: PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap *);
2219: PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap, PetscViewer);

2221: PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2222: PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm, PetscShmComm *);
2223: PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm, PetscMPIInt, PetscMPIInt *);
2224: PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm, PetscMPIInt, PetscMPIInt *);
2225: PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm, MPI_Comm *);

2227: /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2228: PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm, PetscInt, PetscOmpCtrl *);
2229: PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl, MPI_Comm *, MPI_Comm *, PetscBool *);
2230: PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl *);
2231: PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2232: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2233: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);

2235: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t, size_t, PetscSegBuffer *);
2236: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer *);
2237: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer, size_t, void *);
2238: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer, void *);
2239: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer, void *);
2240: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer, void *);
2241: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer, size_t *);
2242: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer, size_t);

2244: /*MC
2245:   PetscSegBufferGetInts - access an array of `PetscInt` from a `PetscSegBuffer`

2247:   Synopsis:
2248: #include <petscsys.h>
2249:   PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg, size_t count, PetscInt *PETSC_RESTRICT *slot);

2251:   No Fortran Support

2253:   Input Parameters:
2254: + seg   - `PetscSegBuffer` buffer
2255: - count - number of entries needed

2257:   Output Parameter:
2258: . buf - address of new buffer for contiguous data

2260:   Level: intermediate

2262:   Developer Note:
2263:   Type-safe wrapper to encourage use of PETSC_RESTRICT. Does not use PetscFunctionBegin because the error handling
2264:   prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as
2265:   possible.

2267: .seealso: `PetscSegBuffer`, `PetscSegBufferGet()`, `PetscInitialize()`, `PetscFinalize()`, `PetscInitializeCalled`
2268: M*/
2269: /* */
2270: static inline PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg, size_t count, PetscInt *PETSC_RESTRICT *slot)
2271: {
2272:   return PetscSegBufferGet(seg, count, (void **)slot);
2273: }

2275: extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2276: PETSC_EXTERN PetscErrorCode    PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted *);
2277: PETSC_EXTERN PetscErrorCode    PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted *);
2278: PETSC_EXTERN PetscErrorCode    PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted, const char *, const char *, PetscBool *);

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

2284: PETSC_EXTERN PetscSegBuffer PetscCitationsList;

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

2289:      Not Collective; No Fortran Support

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

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

2298:      Level: intermediate
2299: @*/
2300: static inline PetscErrorCode PetscCitationsRegister(const char cit[], PetscBool *set)
2301: {
2302:   size_t len;
2303:   char  *vstring;

2305:   PetscFunctionBegin;
2306:   if (set && *set) PetscFunctionReturn(PETSC_SUCCESS);
2307:   PetscCall(PetscStrlen(cit, &len));
2308:   PetscCall(PetscSegBufferGet(PetscCitationsList, len, &vstring));
2309:   PetscCall(PetscArraycpy(vstring, cit, len));
2310:   if (set) *set = PETSC_TRUE;
2311:   PetscFunctionReturn(PETSC_SUCCESS);
2312: }

2314: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm, char[], char[], size_t);
2315: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm, const char[], char[], size_t);
2316: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm, const char[], const char[]);

2318: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm, char[], char[], size_t);
2319: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm, const char[], char[], char[], size_t);
2320: PETSC_EXTERN PetscErrorCode PetscBoxUpload(MPI_Comm, const char[], const char[]);

2322: PETSC_EXTERN PetscErrorCode PetscGlobusGetTransfers(MPI_Comm, const char[], char[], size_t);
2323: PETSC_EXTERN PetscErrorCode PetscGlobusAuthorize(MPI_Comm, char[], size_t);
2324: PETSC_EXTERN PetscErrorCode PetscGlobusUpload(MPI_Comm, const char[], const char[]);

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

2329: #if defined(PETSC_USE_DEBUG)
2330: static inline unsigned int PetscStrHash(const char *str)
2331: {
2332:   unsigned int c, hash = 5381;

2334:   while ((c = (unsigned int)*str++)) hash = ((hash << 5) + hash) + c; /* hash * 33 + c */
2335:   return hash;
2336: }

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

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

2347:    Collective

2349:    Input Parameters:
2350: +  a - pointer to the input data to be reduced
2351: .  c - the number of MPI data items in a and b
2352: .  d - the MPI datatype, for example `MPI_INT`
2353: .  e - the MPI operation, for example `MPI_SUM`
2354: -  fcomm - the MPI communicator on which the operation occurs

2356:    Output Parameter:
2357: .  b - the reduced values

2359:    Level: developer

2361:    Notes:
2362:    In optimized mode this directly calls `MPI_Allreduce()`

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

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

2368: .seealso: `MPI_Allreduce()`
2369: M*/
2370:   #define MPIU_Allreduce(a, b, c, d, e, fcomm) \
2371:     PetscMacroReturnStandard( \
2372:     PetscMPIInt a_b1[6], a_b2[6]; \
2373:     int _mpiu_allreduce_c_int = (int)(c); \
2374:     a_b1[0] = -(PetscMPIInt)__LINE__; \
2375:     a_b1[1] = -a_b1[0]; \
2376:     a_b1[2] = -(PetscMPIInt)PetscStrHash(PETSC_FUNCTION_NAME); \
2377:     a_b1[3] = -a_b1[2]; \
2378:     a_b1[4] = -(PetscMPIInt)(c); \
2379:     a_b1[5] = -a_b1[4]; \
2380:     \
2381:     PetscCallMPI(MPI_Allreduce(a_b1, a_b2, 6, MPI_INT, MPI_MAX, fcomm)); \
2382:     PetscCheck(-a_b2[0] == a_b2[1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "MPI_Allreduce() called in different locations (code lines) on different processors"); \
2383:     PetscCheck(-a_b2[2] == a_b2[3], PETSC_COMM_SELF, PETSC_ERR_PLIB, "MPI_Allreduce() called in different locations (functions) on different processors"); \
2384:     PetscCheck(-a_b2[4] == a_b2[5], PETSC_COMM_SELF, PETSC_ERR_PLIB, "MPI_Allreduce() called with different counts %d on different processors", _mpiu_allreduce_c_int); \
2385:     PetscCallMPI(MPI_Allreduce((a), (b), (c), (d), (e), (fcomm)));)
2386: #else
2387:   #define MPIU_Allreduce(a, b, c, d, e, fcomm) PetscMacroReturnStandard(PetscCallMPI(MPI_Allreduce((a), (b), (c), (d), (e), (fcomm))))
2388: #endif

2390: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
2391: PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint, PetscMPIInt, MPI_Info, MPI_Comm, void *, MPI_Win *);
2392: PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win, PetscMPIInt, MPI_Aint *, PetscMPIInt *, void *);
2393: #endif

2395: /* this is a vile hack */
2396: #if defined(PETSC_HAVE_NECMPI)
2397:   #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)
2398:     #define MPI_Type_free(a) (*(a) = MPI_DATATYPE_NULL, 0);
2399:   #endif
2400: #endif

2402: /*
2403:     List of external packages and queries on it
2404: */
2405: PETSC_EXTERN PetscErrorCode PetscHasExternalPackage(const char[], PetscBool *);

2407: /* this cannot go here because it may be in a different shared library */
2408: PETSC_EXTERN PetscErrorCode PCMPIServerBegin(void);
2409: PETSC_EXTERN PetscErrorCode PCMPIServerEnd(void);
2410: PETSC_EXTERN PetscErrorCode PCMPICommsDestroy(void);
2411: PETSC_EXTERN PetscBool      PCMPIServerActive;

2413: #define PETSC_HAVE_FORTRAN PETSC_DEPRECATED_MACRO(3, 20, 0, "PETSC_USE_FORTRAN_BINDINGS", ) PETSC_USE_FORTRAN_BINDINGS

2415: PETSC_EXTERN PetscErrorCode PetscBLASSetNumThreads(PetscInt);
2416: PETSC_EXTERN PetscErrorCode PetscBLASGetNumThreads(PetscInt *);

2418: /*MC
2419:    PetscSafePointerPlusOffset - Checks that a pointer is not `NULL` before applying an offset

2421:    Level: beginner

2423:    Note:
2424:    This is needed to avoid errors with undefined-behavior sanitizers such as
2425:    UBSan, assuming PETSc has been configured with `-fsanitize=undefined` as part of the compiler flags
2426: M*/
2427: #define PetscSafePointerPlusOffset(ptr, offset) ((ptr) ? (ptr) + (offset) : NULL)