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: Level: beginner
283: Note:
284: By default `PETSC_MPI_THREAD_REQUIRED` equals `MPI_THREAD_FUNNELED` when the MPI implementation provides MPI_Init_thread(), otherwise it equals `MPI_THREAD_SINGLE`
286: .seealso: `PetscInitialize()`
287: M*/
288: PETSC_EXTERN PetscMPIInt PETSC_MPI_THREAD_REQUIRED;
290: PETSC_EXTERN PetscBool PetscBeganMPI;
291: PETSC_EXTERN PetscBool PetscErrorHandlingInitialized;
292: PETSC_EXTERN PetscBool PetscInitializeCalled;
293: PETSC_EXTERN PetscBool PetscFinalizeCalled;
294: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;
296: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm), PetscErrorCode (*)(MPI_Comm));
297: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm, MPI_Comm *, int *);
298: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm *);
299: PETSC_EXTERN PetscErrorCode PetscCommGetComm(MPI_Comm, MPI_Comm *);
300: PETSC_EXTERN PetscErrorCode PetscCommRestoreComm(MPI_Comm, MPI_Comm *);
302: #if defined(PETSC_HAVE_KOKKOS)
303: PETSC_EXTERN PetscErrorCode PetscKokkosInitializeCheck(void); /* Initialize Kokkos if not yet. */
304: #endif
306: #if defined(PETSC_HAVE_NVSHMEM)
307: PETSC_EXTERN PetscBool PetscBeganNvshmem;
308: PETSC_EXTERN PetscBool PetscNvshmemInitialized;
309: PETSC_EXTERN PetscErrorCode PetscNvshmemFinalize(void);
310: #endif
312: #if defined(PETSC_HAVE_ELEMENTAL)
313: PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
314: PETSC_EXTERN PetscErrorCode PetscElementalInitialized(PetscBool *);
315: PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
316: #endif
318: /*MC
319: PetscMalloc - Allocates memory, One should use `PetscNew()`, `PetscMalloc1()` or `PetscCalloc1()` usually instead of this
321: Synopsis:
322: #include <petscsys.h>
323: PetscErrorCode PetscMalloc(size_t m,void **result)
325: Not Collective
327: Input Parameter:
328: . m - number of bytes to allocate
330: Output Parameter:
331: . result - memory allocated
333: Level: beginner
335: Notes:
336: Memory is always allocated at least double aligned
338: It is safe to allocate size 0 and pass the resulting pointer (which may or may not be `NULL`) to `PetscFree()`.
340: .seealso: `PetscFree()`, `PetscNew()`
341: M*/
342: #define PetscMalloc(a, b) ((*PetscTrMalloc)((a), PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (void **)(b)))
344: /*MC
345: PetscRealloc - Reallocates memory
347: Synopsis:
348: #include <petscsys.h>
349: PetscErrorCode PetscRealloc(size_t m,void **result)
351: Not Collective
353: Input Parameters:
354: + m - number of bytes to allocate
355: - result - previous memory
357: Output Parameter:
358: . result - new memory allocated
360: Level: developer
362: Note:
363: Memory is always allocated at least double aligned
365: .seealso: `PetscMalloc()`, `PetscFree()`, `PetscNew()`
366: M*/
367: #define PetscRealloc(a, b) ((*PetscTrRealloc)((a), __LINE__, PETSC_FUNCTION_NAME, __FILE__, (void **)(b)))
369: /*MC
370: PetscAddrAlign - Rounds up an address to `PETSC_MEMALIGN` alignment
372: Synopsis:
373: #include <petscsys.h>
374: void *PetscAddrAlign(void *addr)
376: Not Collective
378: Input Parameter:
379: . addr - address to align (any pointer type)
381: Level: developer
383: .seealso: `PetscMallocAlign()`
384: M*/
385: #define PetscAddrAlign(a) ((void *)((((PETSC_UINTPTR_T)(a)) + (PETSC_MEMALIGN - 1)) & ~(PETSC_MEMALIGN - 1)))
387: /*MC
388: PetscCalloc - Allocates a cleared (zeroed) memory region aligned to `PETSC_MEMALIGN`
390: Synopsis:
391: #include <petscsys.h>
392: PetscErrorCode PetscCalloc(size_t m,void **result)
394: Not Collective
396: Input Parameter:
397: . m - number of bytes to allocate
399: Output Parameter:
400: . result - memory allocated
402: Level: beginner
404: Notes:
405: Memory is always allocated at least double aligned. This macro is useful in allocating memory pointed by void pointers
407: It is safe to allocate size 0 and pass the resulting pointer (which may or may not be `NULL`) to `PetscFree()`.
409: .seealso: `PetscFree()`, `PetscNew()`
410: M*/
411: #define PetscCalloc(m, result) PetscMallocA(1, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)m), (result))
413: /*MC
414: PetscMalloc1 - Allocates an array of memory aligned to `PETSC_MEMALIGN`
416: Synopsis:
417: #include <petscsys.h>
418: PetscErrorCode PetscMalloc1(size_t m1,type **r1)
420: Not Collective
422: Input Parameter:
423: . m1 - number of elements to allocate (may be zero)
425: Output Parameter:
426: . r1 - memory allocated
428: Level: beginner
430: Note:
431: This uses the sizeof() of the memory type requested to determine the total memory to be allocated, therefore you should not
432: multiply the number of elements requested by the `sizeof()` the type. For example use
433: .vb
434: PetscInt *id;
435: PetscMalloc1(10,&id);
436: .ve
437: not
438: .vb
439: PetscInt *id;
440: PetscMalloc1(10*sizeof(PetscInt),&id);
441: .ve
443: Does not zero the memory allocated, use `PetscCalloc1()` to obtain memory that has been zeroed.
445: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc2()`
446: M*/
447: #define PetscMalloc1(m1, r1) PetscMallocA(1, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1))
449: /*MC
450: PetscCalloc1 - Allocates a cleared (zeroed) array of memory aligned to `PETSC_MEMALIGN`
452: Synopsis:
453: #include <petscsys.h>
454: PetscErrorCode PetscCalloc1(size_t m1,type **r1)
456: Not Collective
458: Input Parameter:
459: . m1 - number of elements to allocate in 1st chunk (may be zero)
461: Output Parameter:
462: . r1 - memory allocated
464: Level: beginner
466: Note:
467: See `PetsMalloc1()` for more details on usage.
469: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc2()`
470: M*/
471: #define PetscCalloc1(m1, r1) PetscMallocA(1, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1))
473: /*MC
474: PetscMalloc2 - Allocates 2 arrays of memory both aligned to `PETSC_MEMALIGN`
476: Synopsis:
477: #include <petscsys.h>
478: PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)
480: Not Collective
482: Input Parameters:
483: + m1 - number of elements to allocate in 1st chunk (may be zero)
484: - m2 - number of elements to allocate in 2nd chunk (may be zero)
486: Output Parameters:
487: + r1 - memory allocated in first chunk
488: - r2 - memory allocated in second chunk
490: Level: developer
492: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc2()`
493: M*/
494: #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))
496: /*MC
497: PetscCalloc2 - Allocates 2 cleared (zeroed) arrays of memory both aligned to `PETSC_MEMALIGN`
499: Synopsis:
500: #include <petscsys.h>
501: PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)
503: Not Collective
505: Input Parameters:
506: + m1 - number of elements to allocate in 1st chunk (may be zero)
507: - m2 - number of elements to allocate in 2nd chunk (may be zero)
509: Output Parameters:
510: + r1 - memory allocated in first chunk
511: - r2 - memory allocated in second chunk
513: Level: developer
515: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc2()`
516: M*/
517: #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))
519: /*MC
520: PetscMalloc3 - Allocates 3 arrays of memory, all aligned to `PETSC_MEMALIGN`
522: Synopsis:
523: #include <petscsys.h>
524: PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
526: Not Collective
528: Input Parameters:
529: + m1 - number of elements to allocate in 1st chunk (may be zero)
530: . m2 - number of elements to allocate in 2nd chunk (may be zero)
531: - m3 - number of elements to allocate in 3rd chunk (may be zero)
533: Output Parameters:
534: + r1 - memory allocated in first chunk
535: . r2 - memory allocated in second chunk
536: - r3 - memory allocated in third chunk
538: Level: developer
540: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc3()`, `PetscFree3()`
541: M*/
542: #define PetscMalloc3(m1, r1, m2, r2, m3, r3) \
543: 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))
545: /*MC
546: PetscCalloc3 - Allocates 3 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
548: Synopsis:
549: #include <petscsys.h>
550: PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
552: Not Collective
554: Input Parameters:
555: + m1 - number of elements to allocate in 1st chunk (may be zero)
556: . m2 - number of elements to allocate in 2nd chunk (may be zero)
557: - m3 - number of elements to allocate in 3rd chunk (may be zero)
559: Output Parameters:
560: + r1 - memory allocated in first chunk
561: . r2 - memory allocated in second chunk
562: - r3 - memory allocated in third chunk
564: Level: developer
566: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc2()`, `PetscMalloc3()`, `PetscFree3()`
567: M*/
568: #define PetscCalloc3(m1, r1, m2, r2, m3, r3) \
569: 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))
571: /*MC
572: PetscMalloc4 - Allocates 4 arrays of memory, all aligned to `PETSC_MEMALIGN`
574: Synopsis:
575: #include <petscsys.h>
576: PetscErrorCode PetscMalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
578: Not Collective
580: Input Parameters:
581: + m1 - number of elements to allocate in 1st chunk (may be zero)
582: . m2 - number of elements to allocate in 2nd chunk (may be zero)
583: . m3 - number of elements to allocate in 3rd chunk (may be zero)
584: - m4 - number of elements to allocate in 4th chunk (may be zero)
586: Output Parameters:
587: + r1 - memory allocated in first chunk
588: . r2 - memory allocated in second chunk
589: . r3 - memory allocated in third chunk
590: - r4 - memory allocated in fourth chunk
592: Level: developer
594: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc4()`, `PetscFree4()`
595: M*/
596: #define PetscMalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \
597: 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))
599: /*MC
600: PetscCalloc4 - Allocates 4 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
602: Synopsis:
603: #include <petscsys.h>
604: PetscErrorCode PetscCalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
606: Not Collective
608: Input Parameters:
609: + m1 - number of elements to allocate in 1st chunk (may be zero)
610: . m2 - number of elements to allocate in 2nd chunk (may be zero)
611: . m3 - number of elements to allocate in 3rd chunk (may be zero)
612: - m4 - number of elements to allocate in 4th chunk (may be zero)
614: Output Parameters:
615: + r1 - memory allocated in first chunk
616: . r2 - memory allocated in second chunk
617: . r3 - memory allocated in third chunk
618: - r4 - memory allocated in fourth chunk
620: Level: developer
622: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc4()`, `PetscFree4()`
623: M*/
624: #define PetscCalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \
625: 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))
627: /*MC
628: PetscMalloc5 - Allocates 5 arrays of memory, all aligned to `PETSC_MEMALIGN`
630: Synopsis:
631: #include <petscsys.h>
632: 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)
634: Not Collective
636: Input Parameters:
637: + m1 - number of elements to allocate in 1st chunk (may be zero)
638: . m2 - number of elements to allocate in 2nd chunk (may be zero)
639: . m3 - number of elements to allocate in 3rd chunk (may be zero)
640: . m4 - number of elements to allocate in 4th chunk (may be zero)
641: - m5 - number of elements to allocate in 5th chunk (may be zero)
643: Output Parameters:
644: + r1 - memory allocated in first chunk
645: . r2 - memory allocated in second chunk
646: . r3 - memory allocated in third chunk
647: . r4 - memory allocated in fourth chunk
648: - r5 - memory allocated in fifth chunk
650: Level: developer
652: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc5()`, `PetscFree5()`
653: M*/
654: #define PetscMalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \
655: 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))
657: /*MC
658: PetscCalloc5 - Allocates 5 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
660: Synopsis:
661: #include <petscsys.h>
662: 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)
664: Not Collective
666: Input Parameters:
667: + m1 - number of elements to allocate in 1st chunk (may be zero)
668: . m2 - number of elements to allocate in 2nd chunk (may be zero)
669: . m3 - number of elements to allocate in 3rd chunk (may be zero)
670: . m4 - number of elements to allocate in 4th chunk (may be zero)
671: - m5 - number of elements to allocate in 5th chunk (may be zero)
673: Output Parameters:
674: + r1 - memory allocated in first chunk
675: . r2 - memory allocated in second chunk
676: . r3 - memory allocated in third chunk
677: . r4 - memory allocated in fourth chunk
678: - r5 - memory allocated in fifth chunk
680: Level: developer
682: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc5()`, `PetscFree5()`
683: M*/
684: #define PetscCalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \
685: 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))
687: /*MC
688: PetscMalloc6 - Allocates 6 arrays of memory, all aligned to `PETSC_MEMALIGN`
690: Synopsis:
691: #include <petscsys.h>
692: 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)
694: Not Collective
696: Input Parameters:
697: + m1 - number of elements to allocate in 1st chunk (may be zero)
698: . m2 - number of elements to allocate in 2nd chunk (may be zero)
699: . m3 - number of elements to allocate in 3rd chunk (may be zero)
700: . m4 - number of elements to allocate in 4th chunk (may be zero)
701: . m5 - number of elements to allocate in 5th chunk (may be zero)
702: - m6 - number of elements to allocate in 6th chunk (may be zero)
704: Output Parameteasr:
705: + r1 - memory allocated in first chunk
706: . r2 - memory allocated in second chunk
707: . r3 - memory allocated in third chunk
708: . r4 - memory allocated in fourth chunk
709: . r5 - memory allocated in fifth chunk
710: - r6 - memory allocated in sixth chunk
712: Level: developer
714: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc6()`, `PetscFree3()`, `PetscFree4()`, `PetscFree5()`, `PetscFree6()`
715: M*/
716: #define PetscMalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \
717: 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))
719: /*MC
720: PetscCalloc6 - Allocates 6 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
722: Synopsis:
723: #include <petscsys.h>
724: 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)
726: Not Collective
728: Input Parameters:
729: + m1 - number of elements to allocate in 1st chunk (may be zero)
730: . m2 - number of elements to allocate in 2nd chunk (may be zero)
731: . m3 - number of elements to allocate in 3rd chunk (may be zero)
732: . m4 - number of elements to allocate in 4th chunk (may be zero)
733: . m5 - number of elements to allocate in 5th chunk (may be zero)
734: - m6 - number of elements to allocate in 6th chunk (may be zero)
736: Output Parameters:
737: + r1 - memory allocated in first chunk
738: . r2 - memory allocated in second chunk
739: . r3 - memory allocated in third chunk
740: . r4 - memory allocated in fourth chunk
741: . r5 - memory allocated in fifth chunk
742: - r6 - memory allocated in sixth chunk
744: Level: developer
746: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscMalloc6()`, `PetscFree6()`
747: M*/
748: #define PetscCalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \
749: 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))
751: /*MC
752: PetscMalloc7 - Allocates 7 arrays of memory, all aligned to `PETSC_MEMALIGN`
754: Synopsis:
755: #include <petscsys.h>
756: 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)
758: Not Collective
760: Input Parameters:
761: + m1 - number of elements to allocate in 1st chunk (may be zero)
762: . m2 - number of elements to allocate in 2nd chunk (may be zero)
763: . m3 - number of elements to allocate in 3rd chunk (may be zero)
764: . m4 - number of elements to allocate in 4th chunk (may be zero)
765: . m5 - number of elements to allocate in 5th chunk (may be zero)
766: . m6 - number of elements to allocate in 6th chunk (may be zero)
767: - m7 - number of elements to allocate in 7th chunk (may be zero)
769: Output Parameters:
770: + r1 - memory allocated in first chunk
771: . r2 - memory allocated in second chunk
772: . r3 - memory allocated in third chunk
773: . r4 - memory allocated in fourth chunk
774: . r5 - memory allocated in fifth chunk
775: . r6 - memory allocated in sixth chunk
776: - r7 - memory allocated in seventh chunk
778: Level: developer
780: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc7()`, `PetscFree7()`
781: M*/
782: #define PetscMalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \
783: 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))
785: /*MC
786: PetscCalloc7 - Allocates 7 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
788: Synopsis:
789: #include <petscsys.h>
790: 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)
792: Not Collective
794: Input Parameters:
795: + m1 - number of elements to allocate in 1st chunk (may be zero)
796: . m2 - number of elements to allocate in 2nd chunk (may be zero)
797: . m3 - number of elements to allocate in 3rd chunk (may be zero)
798: . m4 - number of elements to allocate in 4th chunk (may be zero)
799: . m5 - number of elements to allocate in 5th chunk (may be zero)
800: . m6 - number of elements to allocate in 6th chunk (may be zero)
801: - m7 - number of elements to allocate in 7th chunk (may be zero)
803: Output Parameters:
804: + r1 - memory allocated in first chunk
805: . r2 - memory allocated in second chunk
806: . r3 - memory allocated in third chunk
807: . r4 - memory allocated in fourth chunk
808: . r5 - memory allocated in fifth chunk
809: . r6 - memory allocated in sixth chunk
810: - r7 - memory allocated in seventh chunk
812: Level: developer
814: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscMalloc7()`, `PetscFree7()`
815: M*/
816: #define PetscCalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \
817: 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))
819: /*MC
820: PetscNew - Allocates memory of a particular type, zeros the memory! Aligned to `PETSC_MEMALIGN`
822: Synopsis:
823: #include <petscsys.h>
824: PetscErrorCode PetscNew(type **result)
826: Not Collective
828: Output Parameter:
829: . result - memory allocated, sized to match pointer type
831: Level: beginner
833: .seealso: `PetscFree()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc1()`
834: M*/
835: #define PetscNew(b) PetscCalloc1(1, (b))
837: #define PetscNewLog(o, b) PETSC_DEPRECATED_MACRO(3, 18, 0, "PetscNew()", ) PetscNew((b))
839: /*MC
840: PetscFree - Frees memory
842: Synopsis:
843: #include <petscsys.h>
844: PetscErrorCode PetscFree(void *memory)
846: Not Collective
848: Input Parameter:
849: . memory - memory to free (the pointer is ALWAYS set to `NULL` upon success)
851: Level: beginner
853: Notes:
854: Do not free memory obtained with `PetscMalloc2()`, `PetscCalloc2()` etc, they must be freed with `PetscFree2()` etc.
856: It is safe to call `PetscFree()` on a `NULL` pointer.
858: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc1()`
859: M*/
860: #define PetscFree(a) ((PetscErrorCode)((*PetscTrFree)((void *)(a), __LINE__, PETSC_FUNCTION_NAME, __FILE__) || ((a) = PETSC_NULLPTR, PETSC_SUCCESS)))
862: /*MC
863: PetscFree2 - Frees 2 chunks of memory obtained with `PetscMalloc2()`
865: Synopsis:
866: #include <petscsys.h>
867: PetscErrorCode PetscFree2(void *memory1,void *memory2)
869: Not Collective
871: Input Parameters:
872: + memory1 - memory to free
873: - memory2 - 2nd memory to free
875: Level: developer
877: Notes:
878: Memory must have been obtained with `PetscMalloc2()`
880: The arguments need to be in the same order as they were in the call to `PetscMalloc2()`
882: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`
883: M*/
884: #define PetscFree2(m1, m2) PetscFreeA(2, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2))
886: /*MC
887: PetscFree3 - Frees 3 chunks of memory obtained with `PetscMalloc3()`
889: Synopsis:
890: #include <petscsys.h>
891: PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)
893: Not Collective
895: Input Parameters:
896: + memory1 - memory to free
897: . memory2 - 2nd memory to free
898: - memory3 - 3rd memory to free
900: Level: developer
902: Notes:
903: Memory must have been obtained with `PetscMalloc3()`
905: The arguments need to be in the same order as they were in the call to `PetscMalloc3()`
907: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`
908: M*/
909: #define PetscFree3(m1, m2, m3) PetscFreeA(3, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3))
911: /*MC
912: PetscFree4 - Frees 4 chunks of memory obtained with `PetscMalloc4()`
914: Synopsis:
915: #include <petscsys.h>
916: PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)
918: Not Collective
920: Input Parameters:
921: + m1 - memory to free
922: . m2 - 2nd memory to free
923: . m3 - 3rd memory to free
924: - m4 - 4th memory to free
926: Level: developer
928: Notes:
929: Memory must have been obtained with `PetscMalloc4()`
931: The arguments need to be in the same order as they were in the call to `PetscMalloc4()`
933: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`
934: M*/
935: #define PetscFree4(m1, m2, m3, m4) PetscFreeA(4, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4))
937: /*MC
938: PetscFree5 - Frees 5 chunks of memory obtained with `PetscMalloc5()`
940: Synopsis:
941: #include <petscsys.h>
942: PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)
944: Not Collective
946: Input Parameters:
947: + m1 - memory to free
948: . m2 - 2nd memory to free
949: . m3 - 3rd memory to free
950: . m4 - 4th memory to free
951: - m5 - 5th memory to free
953: Level: developer
955: Notes:
956: Memory must have been obtained with `PetscMalloc5()`
958: The arguments need to be in the same order as they were in the call to `PetscMalloc5()`
960: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`
961: M*/
962: #define PetscFree5(m1, m2, m3, m4, m5) PetscFreeA(5, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5))
964: /*MC
965: PetscFree6 - Frees 6 chunks of memory obtained with `PetscMalloc6()`
967: Synopsis:
968: #include <petscsys.h>
969: PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)
971: Not Collective
973: Input Parameters:
974: + m1 - memory to free
975: . m2 - 2nd memory to free
976: . m3 - 3rd memory to free
977: . m4 - 4th memory to free
978: . m5 - 5th memory to free
979: - m6 - 6th memory to free
981: Level: developer
983: Notes:
984: Memory must have been obtained with `PetscMalloc6()`
986: The arguments need to be in the same order as they were in the call to `PetscMalloc6()`
988: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()`
989: M*/
990: #define PetscFree6(m1, m2, m3, m4, m5, m6) PetscFreeA(6, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6))
992: /*MC
993: PetscFree7 - Frees 7 chunks of memory obtained with `PetscMalloc7()`
995: Synopsis:
996: #include <petscsys.h>
997: PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)
999: Not Collective
1001: Input Parameters:
1002: + m1 - memory to free
1003: . m2 - 2nd memory to free
1004: . m3 - 3rd memory to free
1005: . m4 - 4th memory to free
1006: . m5 - 5th memory to free
1007: . m6 - 6th memory to free
1008: - m7 - 7th memory to free
1010: Level: developer
1012: Notes:
1013: Memory must have been obtained with `PetscMalloc7()`
1015: The arguments need to be in the same order as they were in the call to `PetscMalloc7()`
1017: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()`,
1018: `PetscMalloc7()`
1019: M*/
1020: #define PetscFree7(m1, m2, m3, m4, m5, m6, m7) PetscFreeA(7, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6), &(m7))
1022: PETSC_EXTERN PetscErrorCode PetscMallocA(int, PetscBool, int, const char *, const char *, size_t, void *, ...);
1023: PETSC_EXTERN PetscErrorCode PetscFreeA(int, int, const char *, const char *, void *, ...);
1024: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t, PetscBool, int, const char[], const char[], void **);
1025: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void *, int, const char[], const char[]);
1026: PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t, int, const char[], const char[], void **);
1027: PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1028: 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 **));
1029: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);
1031: /*
1032: Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1033: */
1034: PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1035: PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);
1036: #if defined(PETSC_HAVE_CUDA)
1037: PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void);
1038: PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void);
1039: #endif
1040: #if defined(PETSC_HAVE_HIP)
1041: PETSC_EXTERN PetscErrorCode PetscMallocSetHIPHost(void);
1042: PETSC_EXTERN PetscErrorCode PetscMallocResetHIPHost(void);
1043: #endif
1045: #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE
1046: #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION
1048: /*
1049: Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1050: */
1051: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1052: PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1053: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1054: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1055: PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1056: PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int, PetscLogDouble *);
1057: PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool, PetscBool);
1058: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool *, PetscBool *, PetscBool *);
1059: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int, const char[], const char[]);
1060: PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1061: PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool *);
1062: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool);
1063: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool *);
1065: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType, MPI_Datatype *);
1066: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype, PetscDataType *);
1067: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType, size_t *);
1068: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char *, PetscDataType *, PetscBool *);
1070: /*
1071: These are MPI operations for MPI_Allreduce() etc
1072: */
1073: PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1074: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1075: PETSC_EXTERN MPI_Op MPIU_SUM;
1076: PETSC_EXTERN MPI_Op MPIU_MAX;
1077: PETSC_EXTERN MPI_Op MPIU_MIN;
1078: #else
1079: #define MPIU_SUM MPI_SUM
1080: #define MPIU_MAX MPI_MAX
1081: #define MPIU_MIN MPI_MIN
1082: #endif
1083: PETSC_EXTERN MPI_Op Petsc_Garbage_SetIntersectOp;
1084: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm, const PetscInt[], PetscInt *, PetscInt *);
1086: #if (defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128)) || (defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16))
1087: /*MC
1088: MPIU_SUM___FP16___FLOAT128 - MPI_Op that acts as a replacement for `MPI_SUM` with
1089: custom `MPI_Datatype` `MPIU___FLOAT128`, `MPIU___COMPLEX128`, and `MPIU___FP16`.
1091: Level: advanced
1093: Developer Note:
1094: This should be unified with `MPIU_SUM`
1096: .seealso: `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`
1097: M*/
1098: PETSC_EXTERN MPI_Op MPIU_SUM___FP16___FLOAT128;
1099: #endif
1100: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm, const PetscInt[], PetscInt *, PetscInt *);
1102: PETSC_EXTERN PetscErrorCode MPIULong_Send(void *, PetscInt, MPI_Datatype, PetscMPIInt, PetscMPIInt, MPI_Comm) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(1, 3);
1103: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void *, PetscInt, MPI_Datatype, PetscMPIInt, PetscMPIInt, MPI_Comm) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(1, 3);
1105: /*
1106: Defines PETSc error handling.
1107: */
1108: #include <petscerror.h>
1110: PETSC_EXTERN PetscBool PetscCIEnabled; /* code is running in the PETSc test harness CI */
1111: PETSC_EXTERN PetscBool PetscCIEnabledPortableErrorOutput; /* error output is stripped to ensure portability of error messages across systems */
1112: PETSC_EXTERN const char *PetscCIFilename(const char *);
1113: PETSC_EXTERN int PetscCILinenumber(int);
1115: #define PETSC_SMALLEST_CLASSID ((PetscClassId)1211211)
1116: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1117: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1118: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[], PetscClassId *);
1119: PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject, PetscObjectId *);
1120: PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject, PetscObjectId, PetscBool *);
1122: /*
1123: Routines that get memory usage information from the OS
1124: */
1125: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1126: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1127: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1128: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);
1130: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);
1132: /*
1133: Initialization of PETSc
1134: */
1135: PETSC_EXTERN PetscErrorCode PetscInitialize(int *, char ***, const char[], const char[]);
1136: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int, char **, const char[], const char[]);
1137: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1138: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1139: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1140: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1141: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1142: PETSC_EXTERN PetscErrorCode PetscGetArgs(int *, char ***);
1143: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1144: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);
1146: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1147: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);
1148: PETSC_EXTERN PetscErrorCode PetscSysFinalizePackage(void);
1150: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[], const char[]);
1151: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1152: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1153: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject, const char[]);
1155: PETSC_EXTERN PetscErrorCode PetscMonitorCompare(PetscErrorCode (*)(void), void *, PetscErrorCode (*)(void **), PetscErrorCode (*)(void), void *, PetscErrorCode (*)(void **), PetscBool *);
1157: /*
1158: These are so that in extern C code we can caste function pointers to non-extern C
1159: function pointers. Since the regular C++ code expects its function pointers to be C++
1160: */
1161: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1162: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1163: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);
1165: /*
1166: Functions that can act on any PETSc object.
1167: */
1168: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject *);
1169: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject, MPI_Comm *);
1170: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject, PetscClassId *);
1171: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject, const char *[]);
1172: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject, const char *[]);
1173: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject, const char[]);
1174: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject, const char *[]);
1175: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject, PetscInt);
1176: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject, PetscInt *);
1177: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject, PetscObject, PetscInt);
1178: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1179: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject, PetscInt *);
1180: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1181: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject, PetscMPIInt *);
1182: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject, const char[], PetscObject);
1183: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject, const char[]);
1184: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject, const char[], PetscObject *);
1185: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject, const char[], void (*)(void));
1186: #define PetscObjectComposeFunction(a, b, ...) PetscObjectComposeFunction_Private((a), (b), (PetscVoidFunction)(__VA_ARGS__))
1187: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1188: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1189: PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1190: PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject, PetscObject);
1191: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm, PetscMPIInt *);
1193: #include <petscviewertypes.h>
1194: #include <petscoptions.h>
1196: PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer, PetscBool, PetscLogDouble);
1197: PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool *);
1199: PETSC_EXTERN PetscErrorCode PetscObjectsListGetGlobalNumbering(MPI_Comm, PetscInt, PetscObject *, PetscInt *, PetscInt *);
1201: PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer, const char[]);
1202: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject, PetscViewer);
1203: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject, PetscViewer);
1204: #define PetscObjectQueryFunction(obj, name, fptr) PetscObjectQueryFunction_Private((obj), (name), (PetscVoidFunction *)(fptr))
1205: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject, const char[], void (**)(void));
1206: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject, const char[]);
1207: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject, const char[]);
1208: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject, const char[]);
1209: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject, const char *[]);
1210: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject, const char[]);
1211: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1212: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1213: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject, PetscObject, const char[]);
1214: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1215: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject, const char[], PetscBool *);
1216: PETSC_EXTERN PetscErrorCode PetscObjectObjectTypeCompare(PetscObject, PetscObject, PetscBool *);
1217: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject, const char[], PetscBool *);
1218: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject, PetscBool *, const char[], ...);
1219: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject, PetscBool *, const char[], ...);
1220: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1221: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);
1223: #if defined(PETSC_HAVE_SAWS)
1224: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1225: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1226: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject, PetscBool);
1227: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1228: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1229: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1230: PETSC_EXTERN void PetscStackSAWsGrantAccess(void);
1231: PETSC_EXTERN void PetscStackSAWsTakeAccess(void);
1232: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1233: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);
1235: #else
1236: #define PetscSAWsBlock() PETSC_SUCCESS
1237: #define PetscObjectSAWsViewOff(obj) PETSC_SUCCESS
1238: #define PetscObjectSAWsSetBlock(obj, flg) PETSC_SUCCESS
1239: #define PetscObjectSAWsBlock(obj) PETSC_SUCCESS
1240: #define PetscObjectSAWsGrantAccess(obj) PETSC_SUCCESS
1241: #define PetscObjectSAWsTakeAccess(obj) PETSC_SUCCESS
1242: #define PetscStackViewSAWs() PETSC_SUCCESS
1243: #define PetscStackSAWsViewOff() PETSC_SUCCESS
1244: #define PetscStackSAWsTakeAccess()
1245: #define PetscStackSAWsGrantAccess()
1247: #endif
1249: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[], PetscDLMode, PetscDLHandle *);
1250: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1251: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle, const char[], void **);
1252: PETSC_EXTERN PetscErrorCode PetscDLAddr(void (*)(void), char **);
1253: #ifdef PETSC_HAVE_CXX
1254: PETSC_EXTERN PetscErrorCode PetscDemangleSymbol(const char *, char **);
1255: #endif
1257: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void *, PetscStack **);
1259: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE *, PetscBool);
1260: PETSC_EXTERN PetscErrorCode PetscObjectsView(PetscViewer);
1261: PETSC_EXTERN PetscErrorCode PetscObjectsGetObject(const char *, PetscObject *, char **);
1262: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList *);
1263: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList, const char[], PetscObject *);
1264: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList, PetscObject, char **, PetscBool *);
1265: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *, const char[], PetscObject);
1266: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *, const char[]);
1267: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList, PetscObjectList *);
1269: /*
1270: Dynamic library lists. Lists of names of routines in objects or in dynamic
1271: link libraries that will be loaded as needed.
1272: */
1274: #define PetscFunctionListAdd(list, name, fptr) PetscFunctionListAdd_Private((list), (name), (PetscVoidFunction)(fptr))
1275: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList *, const char[], PetscVoidFunction);
1276: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList *);
1277: PETSC_EXTERN PetscErrorCode PetscFunctionListClear(PetscFunctionList);
1278: #define PetscFunctionListFind(list, name, fptr) PetscFunctionListFind_Private((list), (name), (PetscVoidFunction *)(fptr))
1279: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList, const char[], PetscVoidFunction *);
1280: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm, FILE *, const char[], const char[], const char[], const char[], PetscFunctionList, const char[], const char[]);
1281: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList, PetscFunctionList *);
1282: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList, PetscViewer);
1283: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList, const char ***, int *);
1284: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintNonEmpty(PetscFunctionList);
1285: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintAll(void);
1287: PETSC_EXTERN PetscDLLibrary PetscDLLibrariesLoaded;
1288: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm, PetscDLLibrary *, const char[]);
1289: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm, PetscDLLibrary *, const char[]);
1290: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm, PetscDLLibrary *, const char[], const char[], void **);
1291: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1292: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm, const char[], char *, size_t, PetscBool *);
1293: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm, const char[], PetscDLLibrary *);
1294: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);
1296: /*
1297: Useful utility routines
1298: */
1299: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm, PetscInt *, PetscInt *);
1300: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm, PetscInt, PetscInt *, PetscInt *);
1301: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm, PetscInt *, PetscInt *);
1302: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm, PetscMPIInt);
1303: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm, PetscMPIInt);
1304: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1305: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE *);
1306: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm, const PetscInt[2], PetscInt[2]);
1307: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm, const PetscReal[2], PetscReal[2]);
1309: /*MC
1310: PetscNot - negates a logical type value and returns result as a `PetscBool`
1312: Level: beginner
1314: Note:
1315: This is useful in cases like
1316: .vb
1317: int *a;
1318: PetscBool flag = PetscNot(a)
1319: .ve
1320: where !a would not return a `PetscBool` because we cannot provide a cast from int to `PetscBool` in C.
1322: .seealso: `PetscBool`, `PETSC_TRUE`, `PETSC_FALSE`
1323: M*/
1324: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)
1326: /*MC
1327: PetscHelpPrintf - Prints help messages.
1329: Synopsis:
1330: #include <petscsys.h>
1331: PetscErrorCode (*PetscHelpPrintf)(MPI_Comm comm, const char format[],args);
1333: Not Collective, only applies on MPI rank 0; No Fortran Support
1335: Input Parameters:
1336: + comm - the MPI communicator over which the help message is printed
1337: . format - the usual printf() format string
1338: - args - arguments to be printed
1340: Level: developer
1342: Notes:
1343: You can change how help messages are printed by replacing the function pointer with a function that does not simply write to stdout.
1345: To use, write your own function, for example,
1346: .vb
1347: PetscErrorCode mypetschelpprintf(MPI_Comm comm,const char format[],....)
1348: {
1349: PetscFunctionReturn(PETSC_SUCCESS);
1350: }
1351: .ve
1352: then do the assignment
1353: .vb
1354: PetscHelpPrintf = mypetschelpprintf;
1355: .ve
1357: You can do the assignment before `PetscInitialize()`.
1359: The default routine used is called `PetscHelpPrintfDefault()`.
1361: .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscErrorPrintf()`, `PetscHelpPrintfDefault()`
1362: M*/
1363: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1365: /*
1366: Defines PETSc profiling.
1367: */
1368: #include <petsclog.h>
1370: /*
1371: Simple PETSc parallel IO for ASCII printing
1372: */
1373: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[], char[]);
1374: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm, const char[], const char[], FILE **);
1375: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm, FILE *);
1376: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1377: PETSC_EXTERN PetscErrorCode PetscFFlush(FILE *);
1378: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1379: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char *, size_t, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1380: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char *, size_t, const char[], size_t *, ...) PETSC_ATTRIBUTE_FORMAT(3, 5);
1381: PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[], size_t, const char *, PetscInt, const PetscReal[]);
1383: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1384: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1385: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1387: PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char *, size_t *);
1388: PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char *, char *);
1390: #if defined(PETSC_HAVE_POPEN)
1391: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm, const char[], const char[], const char[], FILE **);
1392: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm, FILE *);
1393: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1394: #endif
1396: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1397: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1398: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm, FILE *);
1399: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm, FILE *, size_t, char[]);
1400: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm, const char[], const char[], FILE **);
1401: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char *[]);
1403: PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1404: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer, void **);
1405: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer, void *);
1406: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer *);
1407: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm, PetscContainer *);
1408: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void *));
1409: PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void *);
1411: /*
1412: For use in debuggers
1413: */
1414: PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1415: PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1416: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt, const PetscInt[], PetscViewer);
1417: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt, const PetscReal[], PetscViewer);
1418: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt, const PetscScalar[], PetscViewer);
1420: /*
1421: Basic memory and string operations. These are usually simple wrappers
1422: around the basic Unix system calls, but a few of them have additional
1423: functionality and/or error checking.
1424: */
1425: #include <petscstring.h>
1427: #include <stddef.h>
1428: #include <stdlib.h>
1430: #if defined(PETSC_CLANG_STATIC_ANALYZER)
1431: #define PetscPrefetchBlock(a, b, c, d)
1432: #else
1433: /*MC
1434: PetscPrefetchBlock - Prefetches a block of memory
1436: Synopsis:
1437: #include <petscsys.h>
1438: void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)
1440: Not Collective
1442: Input Parameters:
1443: + a - pointer to first element to fetch (any type but usually `PetscInt` or `PetscScalar`)
1444: . n - number of elements to fetch
1445: . rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors)
1446: - t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note
1448: Level: developer
1450: Notes:
1451: The last two arguments (`rw` and `t`) must be compile-time constants.
1453: Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality. Not all architectures offer
1454: equivalent locality hints, but the following macros are always defined to their closest analogue.
1455: + `PETSC_PREFETCH_HINT_NTA` - Non-temporal. Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched).
1456: . `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.
1457: . `PETSC_PREFETCH_HINT_T1` - Fetch to level 2 and higher (not L1).
1458: - `PETSC_PREFETCH_HINT_T2` - Fetch to high-level cache only. (On many systems, T0 and T1 are equivalent.)
1460: This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid
1461: address).
1463: M*/
1464: #define PetscPrefetchBlock(a, n, rw, t) \
1465: do { \
1466: const char *_p = (const char *)(a), *_end = (const char *)((a) + (n)); \
1467: for (; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p, (rw), (t)); \
1468: } while (0)
1469: #endif
1470: /*
1471: Determine if some of the kernel computation routines use
1472: Fortran (rather than C) for the numerical calculations. On some machines
1473: and compilers (like complex numbers) the Fortran version of the routines
1474: is faster than the C/C++ versions. The flag --with-fortran-kernels
1475: should be used with ./configure to turn these on.
1476: */
1477: #if defined(PETSC_USE_FORTRAN_KERNELS)
1479: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1480: #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1481: #endif
1483: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
1484: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
1485: #endif
1487: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1488: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
1489: #endif
1491: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1492: #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
1493: #endif
1495: #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
1496: #define PETSC_USE_FORTRAN_KERNEL_NORM
1497: #endif
1499: #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
1500: #define PETSC_USE_FORTRAN_KERNEL_MAXPY
1501: #endif
1503: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1504: #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
1505: #endif
1507: #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1508: #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
1509: #endif
1511: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
1512: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
1513: #endif
1515: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1516: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
1517: #endif
1519: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
1520: #define PETSC_USE_FORTRAN_KERNEL_MDOT
1521: #endif
1523: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
1524: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
1525: #endif
1527: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
1528: #define PETSC_USE_FORTRAN_KERNEL_AYPX
1529: #endif
1531: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
1532: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
1533: #endif
1535: #endif
1537: /*
1538: Macros for indicating code that should be compiled with a C interface,
1539: rather than a C++ interface. Any routines that are dynamically loaded
1540: (such as the PCCreate_XXX() routines) must be wrapped so that the name
1541: mangler does not change the functions symbol name. This just hides the
1542: ugly extern "C" {} wrappers.
1543: */
1544: #if defined(__cplusplus)
1545: #define EXTERN_C_BEGIN extern "C" {
1546: #define EXTERN_C_END }
1547: #else
1548: #define EXTERN_C_BEGIN
1549: #define EXTERN_C_END
1550: #endif
1552: /*MC
1553: MPI_Comm - the basic object used by MPI to determine which processes are involved in a
1554: communication
1556: Level: beginner
1558: Note:
1559: This manual page is a place-holder because MPICH does not have a manual page for `MPI_Comm`
1561: .seealso: `PETSC_COMM_WORLD`, `PETSC_COMM_SELF`
1562: M*/
1564: #if defined(PETSC_HAVE_MPIIO)
1565: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4);
1566: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4);
1567: 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);
1568: 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);
1569: 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);
1570: 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);
1571: #endif
1573: /*@C
1574: PetscIntCast - casts a `PetscInt64` (which is 64 bits in size) to a `PetscInt` (which may be 32-bits in size), generates an
1575: error if the `PetscInt` is not large enough to hold the number.
1577: Not Collective; No Fortran Support
1579: Input Parameter:
1580: . a - the `PetscInt64` value
1582: Output Parameter:
1583: . b - the resulting `PetscInt` value
1585: Level: advanced
1587: Note:
1588: 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
1590: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscIntMultError()`, `PetscIntSumError()`
1591: @*/
1592: static inline PetscErrorCode PetscIntCast(PetscInt64 a, PetscInt *b)
1593: {
1594: PetscFunctionBegin;
1595: *b = 0;
1596: // if using 64-bit indices already then this comparison is tautologically true
1597: 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);
1598: *b = (PetscInt)a;
1599: PetscFunctionReturn(PETSC_SUCCESS);
1600: }
1602: /*@C
1603: PetscCountCast - casts a `PetscCount` to a `PetscInt` (which may be 32-bits in size), generates an
1604: error if the `PetscInt` is not large enough to hold the number.
1606: Not Collective; No Fortran Support
1608: Input Parameter:
1609: . a - the `PetscCount` value
1611: Output Parameter:
1612: . b - the resulting `PetscInt` value
1614: Level: advanced
1616: Note:
1617: 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
1619: .seealso: `PetscCount`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscIntMultError()`, `PetscIntSumError()`, `PetscIntCast()`
1620: @*/
1621: static inline PetscErrorCode PetscCountCast(PetscCount a, PetscInt *b)
1622: {
1623: PetscFunctionBegin;
1624: *b = 0;
1625: 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);
1626: *b = (PetscInt)a;
1627: PetscFunctionReturn(PETSC_SUCCESS);
1628: }
1630: /*@C
1631: PetscBLASIntCast - casts a `PetscInt` (which may be 64-bits in size) to a `PetscBLASInt` (which may be 32-bits in size), generates an
1632: error if the `PetscBLASInt` is not large enough to hold the number.
1634: Not Collective; No Fortran Support
1636: Input Parameter:
1637: . a - the `PetscInt` value
1639: Output Parameter:
1640: . b - the resulting `PetscBLASInt` value
1642: Level: advanced
1644: Note:
1645: Errors if the integer is negative since PETSc calls to BLAS/LAPACK never need to cast negative integer inputs
1647: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscIntCast()`, `PetscCountCast()`
1648: @*/
1649: static inline PetscErrorCode PetscBLASIntCast(PetscInt a, PetscBLASInt *b)
1650: {
1651: PetscFunctionBegin;
1652: *b = 0;
1653: if (PetscDefined(USE_64BIT_INDICES) && !PetscDefined(HAVE_64BIT_BLAS_INDICES)) {
1654: 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);
1655: }
1656: PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer to BLAS/LAPACK routine");
1657: *b = (PetscBLASInt)a;
1658: PetscFunctionReturn(PETSC_SUCCESS);
1659: }
1661: /*@C
1662: PetscCuBLASIntCast - like `PetscBLASIntCast()`, but for `PetscCuBLASInt`.
1664: Not Collective; No Fortran Support
1666: Input Parameter:
1667: . a - the `PetscInt` value
1669: Output Parameter:
1670: . b - the resulting `PetscCuBLASInt` value
1672: Level: advanced
1674: Note:
1675: Errors if the integer is negative since PETSc calls to cuBLAS and friends never need to cast negative integer inputs
1677: .seealso: `PetscCuBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscIntCast()`
1678: @*/
1679: static inline PetscErrorCode PetscCuBLASIntCast(PetscInt a, PetscCuBLASInt *b)
1680: {
1681: PetscFunctionBegin;
1682: *b = 0;
1683: 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);
1684: PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer %" PetscInt_FMT "to cuBLAS routine", a);
1685: *b = (PetscCuBLASInt)a;
1686: PetscFunctionReturn(PETSC_SUCCESS);
1687: }
1689: /*@C
1690: PetscHipBLASIntCast - like `PetscBLASIntCast()`, but for `PetscHipBLASInt`.
1692: Not Collective; No Fortran Support
1694: Input Parameter:
1695: . a - the `PetscInt` value
1697: Output Parameter:
1698: . b - the resulting `PetscHipBLASInt` value
1700: Level: advanced
1702: Note:
1703: Errors if the integer is negative since PETSc calls to hipBLAS and friends never need to cast negative integer inputs
1705: .seealso: `PetscHipBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscIntCast()`
1706: @*/
1707: static inline PetscErrorCode PetscHipBLASIntCast(PetscInt a, PetscHipBLASInt *b)
1708: {
1709: PetscFunctionBegin;
1710: *b = 0;
1711: 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);
1712: PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer %" PetscInt_FMT "to hipBLAS routine", a);
1713: *b = (PetscHipBLASInt)a;
1714: PetscFunctionReturn(PETSC_SUCCESS);
1715: }
1717: /*@C
1718: PetscMPIIntCast - casts a `PetscInt` (which may be 64-bits in size) to a `PetscMPIInt` (which may be 32-bits in size), generates an
1719: error if the `PetscMPIInt` is not large enough to hold the number.
1721: Not Collective; No Fortran Support
1723: Input Parameter:
1724: . a - the `PetscInt` value
1726: Output Parameter:
1727: . b - the resulting `PetscMPIInt` value
1729: Level: advanced
1731: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscIntCast()`
1732: @*/
1733: static inline PetscErrorCode PetscMPIIntCast(PetscInt a, PetscMPIInt *b)
1734: {
1735: PetscFunctionBegin;
1736: *b = 0;
1737: 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);
1738: *b = (PetscMPIInt)a;
1739: PetscFunctionReturn(PETSC_SUCCESS);
1740: }
1742: #define PetscInt64Mult(a, b) (((PetscInt64)(a)) * ((PetscInt64)(b)))
1744: /*@C
1745: PetscRealIntMultTruncate - Computes the product of a positive `PetscReal` and a positive
1746: `PetscInt` and truncates the value to slightly less than the maximal possible value.
1748: Not Collective; No Fortran Support
1750: Input Parameters:
1751: + a - The `PetscReal` value
1752: - b - The `PetscInt` value
1754: Level: advanced
1756: Notes:
1757: Returns the result as a `PetscInt` value.
1759: Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`.
1761: Use `PetscIntMultTruncate()` to compute the product of two positive `PetscInt` and truncate
1762: to fit a `PetscInt`.
1764: Use `PetscIntMultError()` to compute the product of two `PetscInt` if you wish to generate an
1765: error if the result will not fit in a `PetscInt`.
1767: Developer Notes:
1768: We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but
1769: requires many more checks.
1771: This is used where we compute approximate sizes for workspace and need to insure the
1772: workspace is index-able.
1774: .seealso: `PetscReal`, `PetscInt`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()`
1775: @*/
1776: static inline PetscInt PetscRealIntMultTruncate(PetscReal a, PetscInt b)
1777: {
1778: PetscInt64 r = (PetscInt64)(a * (PetscReal)b);
1779: if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
1780: return (PetscInt)r;
1781: }
1783: /*@C
1784: PetscIntMultTruncate - Computes the product of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value
1786: Not Collective; No Fortran Support
1788: Input Parameters:
1789: + a - the `PetscInt` value
1790: - b - the second value
1792: Returns:
1793: The result as a `PetscInt` value
1795: Level: advanced
1797: Notes:
1798: Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`
1800: Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt`
1802: 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`
1804: Developer Notes:
1805: We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but requires many more checks.
1807: This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
1809: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()`,
1810: `PetscIntSumTruncate()`
1811: @*/
1812: static inline PetscInt PetscIntMultTruncate(PetscInt a, PetscInt b)
1813: {
1814: PetscInt64 r = PetscInt64Mult(a, b);
1815: if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
1816: return (PetscInt)r;
1817: }
1819: /*@C
1820: PetscIntSumTruncate - Computes the sum of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value
1822: Not Collective; No Fortran Support
1824: Input Parameters:
1825: + a - the `PetscInt` value
1826: - b - the second value
1828: Returns:
1829: The result as a `PetscInt` value
1831: Level: advanced
1833: Notes:
1834: Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`
1836: Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt`
1838: 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`
1840: Developer Note:
1841: This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
1843: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`
1844: @*/
1845: static inline PetscInt PetscIntSumTruncate(PetscInt a, PetscInt b)
1846: {
1847: PetscInt64 r = ((PetscInt64)a) + ((PetscInt64)b);
1848: if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
1849: return (PetscInt)r;
1850: }
1852: /*@C
1853: PetscIntMultError - Computes the product of two positive `PetscInt` and generates an error with overflow.
1855: Not Collective; No Fortran Support
1857: Input Parameters:
1858: + a - the `PetscInt` value
1859: - b - the second value
1861: Output Parameter:
1862: . result - the result as a `PetscInt` value, or `NULL` if you do not want the result, you just want to check if it overflows
1864: Level: advanced
1866: Notes:
1867: Use `PetscInt64Mult()` to compute the product of two `PetscInt` and store in a `PetscInt64`
1869: Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt`
1871: Developer Note:
1872: 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.
1873: `PetscIntSumError()` can be used to check for this situation.
1875: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscIntMult64()`, `PetscIntSumError()`
1876: @*/
1877: static inline PetscErrorCode PetscIntMultError(PetscInt a, PetscInt b, PetscInt *result)
1878: {
1879: PetscInt64 r = PetscInt64Mult(a, b);
1881: PetscFunctionBegin;
1882: if (result) *result = (PetscInt)r;
1883: if (!PetscDefined(USE_64BIT_INDICES)) {
1884: 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);
1885: }
1886: PetscFunctionReturn(PETSC_SUCCESS);
1887: }
1889: /*@C
1891: PetscIntSumError - Computes the sum of two positive `PetscInt` and generates an error with overflow.
1893: Not Collective; No Fortran Support
1895: Input Parameters:
1896: + a - the `PetscInt` value
1897: - b - the second value
1899: Output Parameter:
1900: . c - the result as a `PetscInt` value, or `NULL` if you do not want the result, you just want to check if it overflows
1902: Level: advanced
1904: Notes:
1905: Use `PetscInt64Mult()` to compute the product of two 32-bit `PetscInt` and store in a `PetscInt64`
1907: Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt`
1909: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`
1910: @*/
1911: static inline PetscErrorCode PetscIntSumError(PetscInt a, PetscInt b, PetscInt *result)
1912: {
1913: PetscInt64 r = ((PetscInt64)a) + ((PetscInt64)b);
1915: PetscFunctionBegin;
1916: if (result) *result = (PetscInt)r;
1917: if (!PetscDefined(USE_64BIT_INDICES)) {
1918: 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);
1919: }
1920: PetscFunctionReturn(PETSC_SUCCESS);
1921: }
1923: /*
1924: The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
1925: */
1926: #if defined(hz)
1927: #undef hz
1928: #endif
1930: #if defined(PETSC_HAVE_SYS_TYPES_H)
1931: #include <sys/types.h>
1932: #endif
1934: /*MC
1936: PETSC_VERSION - This manual page provides information about how PETSc documents and uses its version information. This information is available to both C/C++
1937: and Fortran compilers when `petscsys.h` is included.
1939: The current PETSc version and the API for accessing it are defined in <A HREF="../include/petscversion.h.html">include/petscversion.html</A>
1941: The complete version number is given as the triple PETSC_VERSION_MAJOR.PETSC_VERSION_MINOR.PETSC_VERSION_SUBMINOR (in short hand x.y.z)
1943: 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
1944: where only a change in the major version number (x) indicates a change in the API.
1946: 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
1948: 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),
1949: 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
1950: version number (x.y.z).
1952: `PETSC_RELEASE_DATE` is the date the x.y version was released (i.e. the version before any patch releases)
1954: `PETSC_VERSION_DATE` is the date the x.y.z version was released
1956: `PETSC_VERSION_GIT` is the last git commit to the repository given in the form vx.y.z-wwwww
1958: `PETSC_VERSION_DATE_GIT` is the date of the last git commit to the repository
1960: `PETSC_VERSION_()` is deprecated and will eventually be removed.
1962: Level: intermediate
1963: M*/
1965: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[], size_t);
1966: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[], size_t);
1967: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[], size_t);
1968: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[], size_t);
1969: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
1970: PETSC_EXTERN PetscErrorCode PetscGetDate(char[], size_t);
1971: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
1972: PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt *, PetscInt *, PetscInt *, PetscInt *);
1974: PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt, const PetscInt[], PetscBool *);
1975: PETSC_EXTERN PetscErrorCode PetscSortedInt64(PetscInt, const PetscInt64[], PetscBool *);
1976: PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt, const PetscMPIInt[], PetscBool *);
1977: PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt, const PetscReal[], PetscBool *);
1978: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt, PetscInt[]);
1979: PETSC_EXTERN PetscErrorCode PetscSortInt64(PetscInt, PetscInt64[]);
1980: PETSC_EXTERN PetscErrorCode PetscSortCount(PetscInt, PetscCount[]);
1981: PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt, PetscInt[]);
1982: PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *, PetscInt[]);
1983: PETSC_EXTERN PetscErrorCode PetscSortedCheckDupsInt(PetscInt, const PetscInt[], PetscBool *);
1984: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt *, PetscInt[]);
1985: PETSC_EXTERN PetscErrorCode PetscCheckDupsInt(PetscInt, const PetscInt[], PetscBool *);
1986: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt *);
1987: PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt *);
1988: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt, const PetscInt[], PetscInt[]);
1989: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt, const char *[], PetscInt[]);
1990: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt, PetscInt[], PetscInt[]);
1991: PETSC_EXTERN PetscErrorCode PetscSortIntWithCountArray(PetscCount, PetscInt[], PetscCount[]);
1992: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt, PetscInt[], PetscInt[], PetscInt[]);
1993: PETSC_EXTERN PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount, PetscInt[], PetscInt[], PetscCount[]);
1994: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt, PetscMPIInt[]);
1995: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *, PetscMPIInt[]);
1996: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt, PetscMPIInt[], PetscMPIInt[]);
1997: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt, PetscMPIInt[], PetscInt[]);
1998: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt, PetscInt[], PetscScalar[]);
1999: PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt, PetscInt[], void *, size_t, void *);
2000: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt, PetscReal[]);
2001: PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt, PetscReal[], PetscInt[]);
2002: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt, const PetscReal[], PetscInt[]);
2003: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt *, PetscReal[]);
2004: PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal, PetscInt, const PetscReal[], PetscReal, PetscInt *);
2005: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt, PetscInt, PetscScalar[], PetscInt[]);
2006: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt, PetscInt, PetscReal[], PetscInt[]);
2007: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt, const PetscBool[], const PetscInt[], PetscInt *, PetscInt **, PetscInt **, PetscInt **, PetscInt **);
2008: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt, const PetscInt[], const PetscInt[], PetscInt, const PetscInt[], const PetscInt[], PetscInt *, PetscInt **, PetscInt **);
2009: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscInt *, PetscInt **);
2010: PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt, const PetscMPIInt[], PetscInt, const PetscMPIInt[], PetscInt *, PetscMPIInt **);
2011: PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *);
2013: PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt, void *, size_t, int (*)(const void *, const void *, void *), void *);
2014: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt, PetscInt[]);
2015: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt, PetscMPIInt[]);
2016: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt, PetscReal[]);
2017: PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt, void *, size_t, void *, size_t, int (*)(const void *, const void *, void *), void *);
2018: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt, PetscInt[], PetscInt[]);
2019: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt, PetscMPIInt[], PetscMPIInt[]);
2020: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt, PetscReal[], PetscInt[]);
2022: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2023: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[], size_t);
2025: /*J
2026: PetscRandomType - String with the name of a PETSc randomizer
2028: Level: beginner
2030: Note:
2031: To use `PETSCSPRNG` or `PETSCRANDOM123` you must have ./configure PETSc
2032: with the option `--download-sprng` or `--download-random123`. We recommend the default provided with PETSc.
2034: .seealso: `PetscRandomSetType()`, `PetscRandom`, `PetscRandomCreate()`
2035: J*/
2036: typedef const char *PetscRandomType;
2037: #define PETSCRAND "rand"
2038: #define PETSCRAND48 "rand48"
2039: #define PETSCSPRNG "sprng"
2040: #define PETSCRANDER48 "rander48"
2041: #define PETSCRANDOM123 "random123"
2042: #define PETSCCURAND "curand"
2044: /* Logging support */
2045: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;
2047: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);
2048: PETSC_EXTERN PetscErrorCode PetscRandomFinalizePackage(void);
2050: /* Dynamic creation and loading functions */
2051: PETSC_EXTERN PetscFunctionList PetscRandomList;
2053: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[], PetscErrorCode (*)(PetscRandom));
2054: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2055: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2056: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType *);
2057: PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom, PetscObject, const char[]);
2058: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom, PetscViewer);
2060: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm, PetscRandom *);
2061: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom, PetscScalar *);
2062: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom, PetscReal *);
2063: PETSC_EXTERN PetscErrorCode PetscRandomGetValues(PetscRandom, PetscInt, PetscScalar *);
2064: PETSC_EXTERN PetscErrorCode PetscRandomGetValuesReal(PetscRandom, PetscInt, PetscReal *);
2065: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom, PetscScalar *, PetscScalar *);
2066: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom, PetscScalar, PetscScalar);
2067: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom, unsigned long);
2068: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom, unsigned long *);
2069: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2070: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom *);
2072: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[], char[], size_t);
2073: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[], char[], size_t);
2074: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[], size_t);
2075: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[], char[]);
2076: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[], size_t);
2077: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[], char, PetscBool *);
2078: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[], char, PetscBool *);
2079: PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2080: PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2081: PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);
2083: static inline PetscBool PetscBinaryBigEndian(void)
2084: {
2085: long _petsc_v = 1;
2086: return ((char *)&_petsc_v)[0] ? PETSC_FALSE : PETSC_TRUE;
2087: }
2089: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int, void *, PetscInt, PetscInt *, PetscDataType);
2090: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm, int, void *, PetscInt, PetscInt *, PetscDataType);
2091: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int, const void *, PetscInt, PetscDataType);
2092: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm, int, const void *, PetscInt, PetscDataType);
2093: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[], PetscFileMode, int *);
2094: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2095: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm, PetscBool *);
2096: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm, PetscBool *);
2097: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm, char[], size_t);
2098: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm, const char[], char[], size_t, PetscBool *);
2099: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm, const char[], char[], size_t, PetscBool *);
2100: #if defined(PETSC_USE_SOCKET_VIEWER)
2101: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[], int, int *);
2102: #endif
2104: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int, off_t, PetscBinarySeekType, off_t *);
2105: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm, int, off_t, PetscBinarySeekType, off_t *);
2106: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *, PetscDataType, PetscInt);
2108: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2109: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[], PetscBool);
2110: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2111: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char *);
2112: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2113: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2114: PETSC_EXTERN PetscErrorCode PetscWaitOnError(void);
2116: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt *);
2117: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], PetscMPIInt **, PetscMPIInt **);
2118: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt **, PetscMPIInt **, PetscMPIInt **);
2119: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscInt ***, MPI_Request **);
2120: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscScalar ***, MPI_Request **);
2121: 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);
2122: 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);
2123: 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);
2125: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm, PetscBuildTwoSidedType);
2126: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm, PetscBuildTwoSidedType *);
2128: PETSC_EXTERN PetscErrorCode PetscSSEIsEnabled(MPI_Comm, PetscBool *, PetscBool *);
2130: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);
2132: struct _n_PetscSubcomm {
2133: MPI_Comm parent; /* parent communicator */
2134: MPI_Comm dupparent; /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2135: MPI_Comm child; /* the sub-communicator */
2136: PetscMPIInt n; /* num of subcommunicators under the parent communicator */
2137: PetscMPIInt color; /* color of processors belong to this communicator */
2138: PetscMPIInt *subsize; /* size of subcommunicator[color] */
2139: PetscSubcommType type;
2140: char *subcommprefix;
2141: };
2143: static inline MPI_Comm PetscSubcommParent(PetscSubcomm scomm)
2144: {
2145: return scomm->parent;
2146: }
2147: static inline MPI_Comm PetscSubcommChild(PetscSubcomm scomm)
2148: {
2149: return scomm->child;
2150: }
2151: static inline MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm)
2152: {
2153: return scomm->dupparent;
2154: }
2155: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm, PetscSubcomm *);
2156: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm *);
2157: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm, PetscInt);
2158: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm, PetscSubcommType);
2159: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm, PetscMPIInt, PetscMPIInt);
2160: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm, PetscViewer);
2161: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2162: PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm, const char[]);
2163: PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm, MPI_Comm *);
2164: PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm, MPI_Comm *);
2165: PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm, MPI_Comm *);
2167: PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt, PetscHeap *);
2168: PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap, PetscInt, PetscInt);
2169: PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap, PetscInt *, PetscInt *);
2170: PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap, PetscInt *, PetscInt *);
2171: PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap, PetscInt, PetscInt);
2172: PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2173: PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap *);
2174: PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap, PetscViewer);
2176: PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2177: PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm, PetscShmComm *);
2178: PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm, PetscMPIInt, PetscMPIInt *);
2179: PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm, PetscMPIInt, PetscMPIInt *);
2180: PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm, MPI_Comm *);
2182: /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2183: PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm, PetscInt, PetscOmpCtrl *);
2184: PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl, MPI_Comm *, MPI_Comm *, PetscBool *);
2185: PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl *);
2186: PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2187: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2188: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);
2190: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t, size_t, PetscSegBuffer *);
2191: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer *);
2192: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer, size_t, void *);
2193: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer, void *);
2194: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer, void *);
2195: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer, void *);
2196: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer, size_t *);
2197: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer, size_t);
2199: /* Type-safe wrapper to encourage use of PETSC_RESTRICT. Does not use PetscFunctionBegin because the error handling
2200: * prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as
2201: * possible. */
2202: static inline PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg, size_t count, PetscInt *PETSC_RESTRICT *slot)
2203: {
2204: return PetscSegBufferGet(seg, count, (void **)slot);
2205: }
2207: extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2208: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted *);
2209: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted *);
2210: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted, const char *, const char *, PetscBool *);
2212: #include <stdarg.h>
2213: PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char *, size_t, const char[], size_t *, va_list);
2214: PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE *, const char[], va_list);
2216: PETSC_EXTERN PetscSegBuffer PetscCitationsList;
2218: /*@C
2219: PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code.
2221: Not Collective; No Fortran Support
2223: Input Parameters:
2224: + cite - the bibtex item, formatted to displayed on multiple lines nicely
2225: - set - a boolean variable initially set to `PETSC_FALSE`; this is used to insure only a single registration of the citation
2227: Options Database Key:
2228: . -citations [filename] - print out the bibtex entries for the given computation
2230: Level: intermediate
2231: @*/
2232: static inline PetscErrorCode PetscCitationsRegister(const char cit[], PetscBool *set)
2233: {
2234: size_t len;
2235: char *vstring;
2237: PetscFunctionBegin;
2238: if (set && *set) PetscFunctionReturn(PETSC_SUCCESS);
2239: PetscCall(PetscStrlen(cit, &len));
2240: PetscCall(PetscSegBufferGet(PetscCitationsList, len, &vstring));
2241: PetscCall(PetscArraycpy(vstring, cit, len));
2242: if (set) *set = PETSC_TRUE;
2243: PetscFunctionReturn(PETSC_SUCCESS);
2244: }
2246: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm, char[], char[], size_t);
2247: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm, const char[], char[], size_t);
2248: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm, const char[], const char[]);
2250: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm, char[], char[], size_t);
2251: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm, const char[], char[], char[], size_t);
2252: PETSC_EXTERN PetscErrorCode PetscBoxUpload(MPI_Comm, const char[], const char[]);
2254: PETSC_EXTERN PetscErrorCode PetscGlobusGetTransfers(MPI_Comm, const char[], char[], size_t);
2255: PETSC_EXTERN PetscErrorCode PetscGlobusAuthorize(MPI_Comm, char[], size_t);
2256: PETSC_EXTERN PetscErrorCode PetscGlobusUpload(MPI_Comm, const char[], const char[]);
2258: PETSC_EXTERN PetscErrorCode PetscPullJSONValue(const char[], const char[], char[], size_t, PetscBool *);
2259: PETSC_EXTERN PetscErrorCode PetscPushJSONValue(char[], const char[], const char[], size_t);
2261: #if defined(PETSC_USE_DEBUG)
2262: static inline unsigned int PetscStrHash(const char *str)
2263: {
2264: unsigned int c, hash = 5381;
2266: while ((c = (unsigned int)*str++)) hash = ((hash << 5) + hash) + c; /* hash * 33 + c */
2267: return hash;
2268: }
2270: /*MC
2271: MPIU_Allreduce - a PETSc replacement for `MPI_Allreduce()` that tries to determine if the call from all the MPI ranks occur from the
2272: same place in the PETSc code. This helps to detect bugs where different MPI ranks follow different code paths
2273: resulting in inconsistent and incorrect calls to `MPI_Allreduce()`.
2275: Synopsis:
2276: #include <petscsys.h>
2277: PetscErrorCode MPIU_Allreduce(void *indata,void *outdata,PetscMPIInt count,MPI_Datatype datatype, MPI_Op op, MPI_Comm comm);
2279: Collective
2281: Input Parameters:
2282: + a - pointer to the input data to be reduced
2283: . c - the number of MPI data items in a and b
2284: . d - the MPI datatype, for example `MPI_INT`
2285: . e - the MPI operation, for example `MPI_SUM`
2286: - fcomm - the MPI communicator on which the operation occurs
2288: Output Parameter:
2289: . b - the reduced values
2291: Level: developer
2293: Notes:
2294: In optimized mode this directly calls `MPI_Allreduce()`
2296: This is defined as a macro that can return error codes internally so it cannot be used in a subroutine that returns void.
2298: 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
2300: .seealso: `MPI_Allreduce()`
2301: M*/
2302: #define MPIU_Allreduce(a, b, c, d, e, fcomm) \
2303: PetscMacroReturnStandard( \
2304: PetscMPIInt a_b1[6], a_b2[6]; \
2305: int _mpiu_allreduce_c_int = (int)(c); \
2306: a_b1[0] = -(PetscMPIInt)__LINE__; \
2307: a_b1[1] = -a_b1[0]; \
2308: a_b1[2] = -(PetscMPIInt)PetscStrHash(PETSC_FUNCTION_NAME); \
2309: a_b1[3] = -a_b1[2]; \
2310: a_b1[4] = -(PetscMPIInt)(c); \
2311: a_b1[5] = -a_b1[4]; \
2312: \
2313: PetscCallMPI(MPI_Allreduce(a_b1, a_b2, 6, MPI_INT, MPI_MAX, fcomm)); \
2314: PetscCheck(-a_b2[0] == a_b2[1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "MPI_Allreduce() called in different locations (code lines) on different processors"); \
2315: PetscCheck(-a_b2[2] == a_b2[3], PETSC_COMM_SELF, PETSC_ERR_PLIB, "MPI_Allreduce() called in different locations (functions) on different processors"); \
2316: 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); \
2317: PetscCallMPI(MPI_Allreduce((a), (b), (c), (d), (e), (fcomm)));)
2318: #else
2319: #define MPIU_Allreduce(a, b, c, d, e, fcomm) PetscMacroReturnStandard(PetscCallMPI(MPI_Allreduce((a), (b), (c), (d), (e), (fcomm))))
2320: #endif
2322: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
2323: PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint, PetscMPIInt, MPI_Info, MPI_Comm, void *, MPI_Win *);
2324: PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win, PetscMPIInt, MPI_Aint *, PetscMPIInt *, void *);
2325: #endif
2327: /* this is a vile hack */
2328: #if defined(PETSC_HAVE_NECMPI)
2329: #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)
2330: #define MPI_Type_free(a) (*(a) = MPI_DATATYPE_NULL, 0);
2331: #endif
2332: #endif
2334: /*
2335: List of external packages and queries on it
2336: */
2337: PETSC_EXTERN PetscErrorCode PetscHasExternalPackage(const char[], PetscBool *);
2339: /* this cannot go here because it may be in a different shared library */
2340: PETSC_EXTERN PetscErrorCode PCMPIServerBegin(void);
2341: PETSC_EXTERN PetscErrorCode PCMPIServerEnd(void);
2342: PETSC_EXTERN PetscErrorCode PCMPICommsDestroy(void);
2343: PETSC_EXTERN PetscBool PCMPIServerActive;
2345: #define PETSC_HAVE_FORTRAN PETSC_DEPRECATED_MACRO(3, 20, 0, "PETSC_USE_FORTRAN_BINDINGS", ) PETSC_USE_FORTRAN_BINDINGS