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 OpenMPI but now appears to be compiling using a non-OpenMPI 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 OpenMPI mpi.h version but now appears to be compiling using a different OpenMPI 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 OpenMPI, 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
222: where you wish PETSc to compute the required value.
224: Level: beginner
226: Developer Note:
227: I would like to use const `PetscInt` `PETSC_DETERMINE` = `PETSC_DECIDE`; but for
228: some reason this is not allowed by the standard even though `PETSC_DECIDE` is a constant value.
230: .seealso: `PETSC_DECIDE`, `PETSC_DEFAULT`, `PETSC_IGNORE`, `VecSetSizes()`
231: M*/
233: /*MC
234: PETSC_DEFAULT - standard way of passing in integer or floating point parameter
235: where you wish PETSc to use the default.
237: Level: beginner
239: Fortran Note:
240: You need to use `PETSC_DEFAULT_INTEGER` or `PETSC_DEFAULT_REAL`.
242: .seealso: `PETSC_DECIDE`, `PETSC_IGNORE`, `PETSC_DETERMINE`
243: M*/
245: /* These MUST be preprocessor defines! see https://gitlab.com/petsc/petsc/-/issues/1370 */
246: #define PETSC_DECIDE (-1)
247: #define PETSC_DETERMINE PETSC_DECIDE
248: #define PETSC_DEFAULT (-2)
250: /*MC
251: PETSC_COMM_WORLD - the equivalent of the `MPI_COMM_WORLD` communicator which represents
252: all the processes that PETSc knows about.
254: Level: beginner
256: Notes:
257: By default `PETSC_COMM_WORLD` and `MPI_COMM_WORLD` are identical unless you wish to
258: run PETSc on ONLY a subset of `MPI_COMM_WORLD`. In that case create your new (smaller)
259: communicator, call it, say comm, and set `PETSC_COMM_WORLD` = comm BEFORE calling
260: PetscInitialize(), but after `MPI_Init()` has been called.
262: The value of `PETSC_COMM_WORLD` should never be USED/accessed before `PetscInitialize()`
263: is called because it may not have a valid value yet.
265: .seealso: `PETSC_COMM_SELF`
266: M*/
267: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;
269: /*MC
270: PETSC_COMM_SELF - This is always `MPI_COMM_SELF`
272: Level: beginner
274: Notes:
275: Do not USE/access or set this variable before `PetscInitialize()` has been called.
277: .seealso: `PETSC_COMM_WORLD`
278: M*/
279: #define PETSC_COMM_SELF MPI_COMM_SELF
281: /*MC
282: PETSC_MPI_THREAD_REQUIRED - the required threading support used if PETSc initializes
283: MPI with `MPI_Init_thread()`.
285: Level: beginner
287: Notes:
288: By default `PETSC_MPI_THREAD_REQUIRED` equals `MPI_THREAD_FUNNELED`.
290: .seealso: `PetscInitialize()`
291: M*/
292: PETSC_EXTERN PetscMPIInt PETSC_MPI_THREAD_REQUIRED;
294: PETSC_EXTERN PetscBool PetscBeganMPI;
295: PETSC_EXTERN PetscBool PetscErrorHandlingInitialized;
296: PETSC_EXTERN PetscBool PetscInitializeCalled;
297: PETSC_EXTERN PetscBool PetscFinalizeCalled;
298: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;
300: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm), PetscErrorCode (*)(MPI_Comm));
301: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm, MPI_Comm *, int *);
302: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm *);
303: PETSC_EXTERN PetscErrorCode PetscCommGetComm(MPI_Comm, MPI_Comm *);
304: PETSC_EXTERN PetscErrorCode PetscCommRestoreComm(MPI_Comm, MPI_Comm *);
306: #if defined(PETSC_HAVE_KOKKOS)
307: PETSC_EXTERN PetscErrorCode PetscKokkosInitializeCheck(void); /* Initialize Kokkos if not yet. */
308: #endif
310: #if defined(PETSC_HAVE_NVSHMEM)
311: PETSC_EXTERN PetscBool PetscBeganNvshmem;
312: PETSC_EXTERN PetscBool PetscNvshmemInitialized;
313: PETSC_EXTERN PetscErrorCode PetscNvshmemFinalize(void);
314: #endif
316: #if defined(PETSC_HAVE_ELEMENTAL)
317: PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
318: PETSC_EXTERN PetscErrorCode PetscElementalInitialized(PetscBool *);
319: PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
320: #endif
322: /*MC
323: PetscMalloc - Allocates memory, One should use `PetscNew()`, `PetscMalloc1()` or `PetscCalloc1()` usually instead of this
325: Synopsis:
326: #include <petscsys.h>
327: PetscErrorCode PetscMalloc(size_t m,void **result)
329: Not Collective
331: Input Parameter:
332: . m - number of bytes to allocate
334: Output Parameter:
335: . result - memory allocated
337: Level: beginner
339: Notes:
340: Memory is always allocated at least double aligned
342: It is safe to allocate size 0 and pass the resulting pointer (which may or may not be `NULL`) to `PetscFree()`.
344: .seealso: `PetscFree()`, `PetscNew()`
345: M*/
346: #define PetscMalloc(a, b) ((*PetscTrMalloc)((a), PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (void **)(b)))
348: /*MC
349: PetscRealloc - Reallocates memory
351: Synopsis:
352: #include <petscsys.h>
353: PetscErrorCode PetscRealloc(size_t m,void **result)
355: Not Collective
357: Input Parameters:
358: + m - number of bytes to allocate
359: - result - previous memory
361: Output Parameter:
362: . result - new memory allocated
364: Level: developer
366: Notes:
367: Memory is always allocated at least double aligned
369: .seealso: `PetscMalloc()`, `PetscFree()`, `PetscNew()`
370: M*/
371: #define PetscRealloc(a, b) ((*PetscTrRealloc)((a), __LINE__, PETSC_FUNCTION_NAME, __FILE__, (void **)(b)))
373: /*MC
374: PetscAddrAlign - Rounds up an address to `PETSC_MEMALIGN` alignment
376: Synopsis:
377: #include <petscsys.h>
378: void *PetscAddrAlign(void *addr)
380: Not Collective
382: Input Parameter:
383: . addr - address to align (any pointer type)
385: Level: developer
387: .seealso: `PetscMallocAlign()`
388: M*/
389: #define PetscAddrAlign(a) ((void *)((((PETSC_UINTPTR_T)(a)) + (PETSC_MEMALIGN - 1)) & ~(PETSC_MEMALIGN - 1)))
391: /*MC
392: PetscCalloc - Allocates a cleared (zeroed) memory region aligned to `PETSC_MEMALIGN`
394: Synopsis:
395: #include <petscsys.h>
396: PetscErrorCode PetscCalloc(size_t m,void **result)
398: Not Collective
400: Input Parameter:
401: . m - number of bytes to allocate
403: Output Parameter:
404: . result - memory allocated
406: Level: beginner
408: Notes:
409: Memory is always allocated at least double aligned. This macro is useful in allocating memory pointed by void pointers
411: It is safe to allocate size 0 and pass the resulting pointer (which may or may not be `NULL`) to `PetscFree()`.
413: .seealso: `PetscFree()`, `PetscNew()`
414: M*/
415: #define PetscCalloc(m, result) PetscMallocA(1, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)m), (result))
417: /*MC
418: PetscMalloc1 - Allocates an array of memory aligned to `PETSC_MEMALIGN`
420: Synopsis:
421: #include <petscsys.h>
422: PetscErrorCode PetscMalloc1(size_t m1,type **r1)
424: Not Collective
426: Input Parameter:
427: . m1 - number of elements to allocate (may be zero)
429: Output Parameter:
430: . r1 - memory allocated
432: Level: beginner
434: Note:
435: This uses the sizeof() of the memory type requested to determine the total memory to be allocated, therefore you should not
436: multiply the number of elements requested by the `sizeof()` the type. For example use
437: .vb
438: PetscInt *id;
439: PetscMalloc1(10,&id);
440: .ve
441: not
442: .vb
443: PetscInt *id;
444: PetscMalloc1(10*sizeof(PetscInt),&id);
445: .ve
447: Does not zero the memory allocated, use `PetscCalloc1()` to obtain memory that has been zeroed.
449: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc2()`
450: M*/
451: #define PetscMalloc1(m1, r1) PetscMallocA(1, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1))
453: /*MC
454: PetscCalloc1 - Allocates a cleared (zeroed) array of memory aligned to `PETSC_MEMALIGN`
456: Synopsis:
457: #include <petscsys.h>
458: PetscErrorCode PetscCalloc1(size_t m1,type **r1)
460: Not Collective
462: Input Parameter:
463: . m1 - number of elements to allocate in 1st chunk (may be zero)
465: Output Parameter:
466: . r1 - memory allocated
468: Level: beginner
470: Notes:
471: See `PetsMalloc1()` for more details on usage.
473: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc2()`
474: M*/
475: #define PetscCalloc1(m1, r1) PetscMallocA(1, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1))
477: /*MC
478: PetscMalloc2 - Allocates 2 arrays of memory both aligned to `PETSC_MEMALIGN`
480: Synopsis:
481: #include <petscsys.h>
482: PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)
484: Not Collective
486: Input Parameters:
487: + m1 - number of elements to allocate in 1st chunk (may be zero)
488: - m2 - number of elements to allocate in 2nd chunk (may be zero)
490: Output Parameters:
491: + r1 - memory allocated in first chunk
492: - r2 - memory allocated in second chunk
494: Level: developer
496: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc2()`
497: M*/
498: #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))
500: /*MC
501: PetscCalloc2 - Allocates 2 cleared (zeroed) arrays of memory both aligned to `PETSC_MEMALIGN`
503: Synopsis:
504: #include <petscsys.h>
505: PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)
507: Not Collective
509: Input Parameters:
510: + m1 - number of elements to allocate in 1st chunk (may be zero)
511: - m2 - number of elements to allocate in 2nd chunk (may be zero)
513: Output Parameters:
514: + r1 - memory allocated in first chunk
515: - r2 - memory allocated in second chunk
517: Level: developer
519: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc2()`
520: M*/
521: #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))
523: /*MC
524: PetscMalloc3 - Allocates 3 arrays of memory, all aligned to `PETSC_MEMALIGN`
526: Synopsis:
527: #include <petscsys.h>
528: PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
530: Not Collective
532: Input Parameters:
533: + m1 - number of elements to allocate in 1st chunk (may be zero)
534: . m2 - number of elements to allocate in 2nd chunk (may be zero)
535: - m3 - number of elements to allocate in 3rd chunk (may be zero)
537: Output Parameters:
538: + r1 - memory allocated in first chunk
539: . r2 - memory allocated in second chunk
540: - r3 - memory allocated in third chunk
542: Level: developer
544: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc3()`, `PetscFree3()`
545: M*/
546: #define PetscMalloc3(m1, r1, m2, r2, m3, r3) \
547: 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))
549: /*MC
550: PetscCalloc3 - Allocates 3 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
552: Synopsis:
553: #include <petscsys.h>
554: PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
556: Not Collective
558: Input Parameters:
559: + m1 - number of elements to allocate in 1st chunk (may be zero)
560: . m2 - number of elements to allocate in 2nd chunk (may be zero)
561: - m3 - number of elements to allocate in 3rd chunk (may be zero)
563: Output Parameters:
564: + r1 - memory allocated in first chunk
565: . r2 - memory allocated in second chunk
566: - r3 - memory allocated in third chunk
568: Level: developer
570: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc2()`, `PetscMalloc3()`, `PetscFree3()`
571: M*/
572: #define PetscCalloc3(m1, r1, m2, r2, m3, r3) \
573: 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))
575: /*MC
576: PetscMalloc4 - Allocates 4 arrays of memory, all aligned to `PETSC_MEMALIGN`
578: Synopsis:
579: #include <petscsys.h>
580: PetscErrorCode PetscMalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
582: Not Collective
584: Input Parameters:
585: + m1 - number of elements to allocate in 1st chunk (may be zero)
586: . m2 - number of elements to allocate in 2nd chunk (may be zero)
587: . m3 - number of elements to allocate in 3rd chunk (may be zero)
588: - m4 - number of elements to allocate in 4th chunk (may be zero)
590: Output Parameters:
591: + r1 - memory allocated in first chunk
592: . r2 - memory allocated in second chunk
593: . r3 - memory allocated in third chunk
594: - r4 - memory allocated in fourth chunk
596: Level: developer
598: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc4()`, `PetscFree4()`
599: M*/
600: #define PetscMalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \
601: 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))
603: /*MC
604: PetscCalloc4 - Allocates 4 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
606: Synopsis:
607: #include <petscsys.h>
608: PetscErrorCode PetscCalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
610: Not Collective
612: Input Parameters:
613: + m1 - number of elements to allocate in 1st chunk (may be zero)
614: . m2 - number of elements to allocate in 2nd chunk (may be zero)
615: . m3 - number of elements to allocate in 3rd chunk (may be zero)
616: - m4 - number of elements to allocate in 4th chunk (may be zero)
618: Output Parameters:
619: + r1 - memory allocated in first chunk
620: . r2 - memory allocated in second chunk
621: . r3 - memory allocated in third chunk
622: - r4 - memory allocated in fourth chunk
624: Level: developer
626: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc4()`, `PetscFree4()`
627: M*/
628: #define PetscCalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \
629: 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))
631: /*MC
632: PetscMalloc5 - Allocates 5 arrays of memory, all aligned to `PETSC_MEMALIGN`
634: Synopsis:
635: #include <petscsys.h>
636: 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)
638: Not Collective
640: Input Parameters:
641: + m1 - number of elements to allocate in 1st chunk (may be zero)
642: . m2 - number of elements to allocate in 2nd chunk (may be zero)
643: . m3 - number of elements to allocate in 3rd chunk (may be zero)
644: . m4 - number of elements to allocate in 4th chunk (may be zero)
645: - m5 - number of elements to allocate in 5th chunk (may be zero)
647: Output Parameters:
648: + r1 - memory allocated in first chunk
649: . r2 - memory allocated in second chunk
650: . r3 - memory allocated in third chunk
651: . r4 - memory allocated in fourth chunk
652: - r5 - memory allocated in fifth chunk
654: Level: developer
656: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc5()`, `PetscFree5()`
657: M*/
658: #define PetscMalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \
659: 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))
661: /*MC
662: PetscCalloc5 - Allocates 5 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
664: Synopsis:
665: #include <petscsys.h>
666: 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)
668: Not Collective
670: Input Parameters:
671: + m1 - number of elements to allocate in 1st chunk (may be zero)
672: . m2 - number of elements to allocate in 2nd chunk (may be zero)
673: . m3 - number of elements to allocate in 3rd chunk (may be zero)
674: . m4 - number of elements to allocate in 4th chunk (may be zero)
675: - m5 - number of elements to allocate in 5th chunk (may be zero)
677: Output Parameters:
678: + r1 - memory allocated in first chunk
679: . r2 - memory allocated in second chunk
680: . r3 - memory allocated in third chunk
681: . r4 - memory allocated in fourth chunk
682: - r5 - memory allocated in fifth chunk
684: Level: developer
686: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc5()`, `PetscFree5()`
687: M*/
688: #define PetscCalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \
689: 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))
691: /*MC
692: PetscMalloc6 - Allocates 6 arrays of memory, all aligned to `PETSC_MEMALIGN`
694: Synopsis:
695: #include <petscsys.h>
696: 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)
698: Not Collective
700: Input Parameters:
701: + m1 - number of elements to allocate in 1st chunk (may be zero)
702: . m2 - number of elements to allocate in 2nd chunk (may be zero)
703: . m3 - number of elements to allocate in 3rd chunk (may be zero)
704: . m4 - number of elements to allocate in 4th chunk (may be zero)
705: . m5 - number of elements to allocate in 5th chunk (may be zero)
706: - m6 - number of elements to allocate in 6th chunk (may be zero)
708: Output Parameteasr:
709: + r1 - memory allocated in first chunk
710: . r2 - memory allocated in second chunk
711: . r3 - memory allocated in third chunk
712: . r4 - memory allocated in fourth chunk
713: . r5 - memory allocated in fifth chunk
714: - r6 - memory allocated in sixth chunk
716: Level: developer
718: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc6()`, `PetscFree3()`, `PetscFree4()`, `PetscFree5()`, `PetscFree6()`
719: M*/
720: #define PetscMalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \
721: 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))
723: /*MC
724: PetscCalloc6 - Allocates 6 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
726: Synopsis:
727: #include <petscsys.h>
728: 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)
730: Not Collective
732: Input Parameters:
733: + m1 - number of elements to allocate in 1st chunk (may be zero)
734: . m2 - number of elements to allocate in 2nd chunk (may be zero)
735: . m3 - number of elements to allocate in 3rd chunk (may be zero)
736: . m4 - number of elements to allocate in 4th chunk (may be zero)
737: . m5 - number of elements to allocate in 5th chunk (may be zero)
738: - m6 - number of elements to allocate in 6th chunk (may be zero)
740: Output Parameters:
741: + r1 - memory allocated in first chunk
742: . r2 - memory allocated in second chunk
743: . r3 - memory allocated in third chunk
744: . r4 - memory allocated in fourth chunk
745: . r5 - memory allocated in fifth chunk
746: - r6 - memory allocated in sixth chunk
748: Level: developer
750: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscMalloc6()`, `PetscFree6()`
751: M*/
752: #define PetscCalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \
753: 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))
755: /*MC
756: PetscMalloc7 - Allocates 7 arrays of memory, all aligned to `PETSC_MEMALIGN`
758: Synopsis:
759: #include <petscsys.h>
760: 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)
762: Not Collective
764: Input Parameters:
765: + m1 - number of elements to allocate in 1st chunk (may be zero)
766: . m2 - number of elements to allocate in 2nd chunk (may be zero)
767: . m3 - number of elements to allocate in 3rd chunk (may be zero)
768: . m4 - number of elements to allocate in 4th chunk (may be zero)
769: . m5 - number of elements to allocate in 5th chunk (may be zero)
770: . m6 - number of elements to allocate in 6th chunk (may be zero)
771: - m7 - number of elements to allocate in 7th chunk (may be zero)
773: Output Parameters:
774: + r1 - memory allocated in first chunk
775: . r2 - memory allocated in second chunk
776: . r3 - memory allocated in third chunk
777: . r4 - memory allocated in fourth chunk
778: . r5 - memory allocated in fifth chunk
779: . r6 - memory allocated in sixth chunk
780: - r7 - memory allocated in seventh chunk
782: Level: developer
784: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc7()`, `PetscFree7()`
785: M*/
786: #define PetscMalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \
787: 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))
789: /*MC
790: PetscCalloc7 - Allocates 7 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
792: Synopsis:
793: #include <petscsys.h>
794: 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)
796: Not Collective
798: Input Parameters:
799: + m1 - number of elements to allocate in 1st chunk (may be zero)
800: . m2 - number of elements to allocate in 2nd chunk (may be zero)
801: . m3 - number of elements to allocate in 3rd chunk (may be zero)
802: . m4 - number of elements to allocate in 4th chunk (may be zero)
803: . m5 - number of elements to allocate in 5th chunk (may be zero)
804: . m6 - number of elements to allocate in 6th chunk (may be zero)
805: - m7 - number of elements to allocate in 7th chunk (may be zero)
807: Output Parameters:
808: + r1 - memory allocated in first chunk
809: . r2 - memory allocated in second chunk
810: . r3 - memory allocated in third chunk
811: . r4 - memory allocated in fourth chunk
812: . r5 - memory allocated in fifth chunk
813: . r6 - memory allocated in sixth chunk
814: - r7 - memory allocated in seventh chunk
816: Level: developer
818: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscMalloc7()`, `PetscFree7()`
819: M*/
820: #define PetscCalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \
821: 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))
823: /*MC
824: PetscNew - Allocates memory of a particular type, zeros the memory! Aligned to `PETSC_MEMALIGN`
826: Synopsis:
827: #include <petscsys.h>
828: PetscErrorCode PetscNew(type **result)
830: Not Collective
832: Output Parameter:
833: . result - memory allocated, sized to match pointer type
835: Level: beginner
837: .seealso: `PetscFree()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc1()`
838: M*/
839: #define PetscNew(b) PetscCalloc1(1, (b))
841: #define PetscNewLog(o, b) PETSC_DEPRECATED_MACRO(3, 18, 0, "PetscNew()", ) PetscNew((b))
843: /*MC
844: PetscFree - Frees memory
846: Synopsis:
847: #include <petscsys.h>
848: PetscErrorCode PetscFree(void *memory)
850: Not Collective
852: Input Parameter:
853: . memory - memory to free (the pointer is ALWAYS set to `NULL` upon success)
855: Level: beginner
857: Note:
858: Do not free memory obtained with `PetscMalloc2()`, `PetscCalloc2()` etc, they must be freed with `PetscFree2()` etc.
860: It is safe to call `PetscFree()` on a `NULL` pointer.
862: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc1()`
863: M*/
864: #define PetscFree(a) ((PetscErrorCode)((*PetscTrFree)((void *)(a), __LINE__, PETSC_FUNCTION_NAME, __FILE__) || ((a) = PETSC_NULLPTR, PETSC_SUCCESS)))
866: /*MC
867: PetscFree2 - Frees 2 chunks of memory obtained with `PetscMalloc2()`
869: Synopsis:
870: #include <petscsys.h>
871: PetscErrorCode PetscFree2(void *memory1,void *memory2)
873: Not Collective
875: Input Parameters:
876: + memory1 - memory to free
877: - memory2 - 2nd memory to free
879: Level: developer
881: Note:
882: Memory must have been obtained with `PetscMalloc2()`
884: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`
885: M*/
886: #define PetscFree2(m1, m2) PetscFreeA(2, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2))
888: /*MC
889: PetscFree3 - Frees 3 chunks of memory obtained with `PetscMalloc3()`
891: Synopsis:
892: #include <petscsys.h>
893: PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)
895: Not Collective
897: Input Parameters:
898: + memory1 - memory to free
899: . memory2 - 2nd memory to free
900: - memory3 - 3rd memory to free
902: Level: developer
904: Note:
905: Memory must have been obtained with `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: Note:
929: Memory must have been obtained with `PetscMalloc4()`
931: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`
932: M*/
933: #define PetscFree4(m1, m2, m3, m4) PetscFreeA(4, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4))
935: /*MC
936: PetscFree5 - Frees 5 chunks of memory obtained with `PetscMalloc5()`
938: Synopsis:
939: #include <petscsys.h>
940: PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)
942: Not Collective
944: Input Parameters:
945: + m1 - memory to free
946: . m2 - 2nd memory to free
947: . m3 - 3rd memory to free
948: . m4 - 4th memory to free
949: - m5 - 5th memory to free
951: Level: developer
953: Note:
954: Memory must have been obtained with `PetscMalloc5()`
956: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`
957: M*/
958: #define PetscFree5(m1, m2, m3, m4, m5) PetscFreeA(5, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5))
960: /*MC
961: PetscFree6 - Frees 6 chunks of memory obtained with `PetscMalloc6()`
963: Synopsis:
964: #include <petscsys.h>
965: PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)
967: Not Collective
969: Input Parameters:
970: + m1 - memory to free
971: . m2 - 2nd memory to free
972: . m3 - 3rd memory to free
973: . m4 - 4th memory to free
974: . m5 - 5th memory to free
975: - m6 - 6th memory to free
977: Level: developer
979: Note:
980: Memory must have been obtained with `PetscMalloc6()`
982: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()`
983: M*/
984: #define PetscFree6(m1, m2, m3, m4, m5, m6) PetscFreeA(6, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6))
986: /*MC
987: PetscFree7 - Frees 7 chunks of memory obtained with `PetscMalloc7()`
989: Synopsis:
990: #include <petscsys.h>
991: PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)
993: Not Collective
995: Input Parameters:
996: + m1 - memory to free
997: . m2 - 2nd memory to free
998: . m3 - 3rd memory to free
999: . m4 - 4th memory to free
1000: . m5 - 5th memory to free
1001: . m6 - 6th memory to free
1002: - m7 - 7th memory to free
1004: Level: developer
1006: Note:
1007: Memory must have been obtained with `PetscMalloc7()`
1009: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()`,
1010: `PetscMalloc7()`
1011: M*/
1012: #define PetscFree7(m1, m2, m3, m4, m5, m6, m7) PetscFreeA(7, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6), &(m7))
1014: PETSC_EXTERN PetscErrorCode PetscMallocA(int, PetscBool, int, const char *, const char *, size_t, void *, ...);
1015: PETSC_EXTERN PetscErrorCode PetscFreeA(int, int, const char *, const char *, void *, ...);
1016: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t, PetscBool, int, const char[], const char[], void **);
1017: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void *, int, const char[], const char[]);
1018: PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t, int, const char[], const char[], void **);
1019: PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1020: 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 **));
1021: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);
1023: /*
1024: Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1025: */
1026: PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1027: PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);
1028: #if defined(PETSC_HAVE_CUDA)
1029: PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void);
1030: PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void);
1031: #endif
1032: #if defined(PETSC_HAVE_HIP)
1033: PETSC_EXTERN PetscErrorCode PetscMallocSetHIPHost(void);
1034: PETSC_EXTERN PetscErrorCode PetscMallocResetHIPHost(void);
1035: #endif
1037: #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE
1038: #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION
1040: /*
1041: Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1042: */
1043: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1044: PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1045: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1046: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1047: PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1048: PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int, PetscLogDouble *);
1049: PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool, PetscBool);
1050: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool *, PetscBool *, PetscBool *);
1051: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int, const char[], const char[]);
1052: PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1053: PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool *);
1054: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool);
1055: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool *);
1057: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType, MPI_Datatype *);
1058: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype, PetscDataType *);
1059: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType, size_t *);
1060: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char *, PetscDataType *, PetscBool *);
1062: /*
1063: These are MPI operations for MPI_Allreduce() etc
1064: */
1065: PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1066: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1067: PETSC_EXTERN MPI_Op MPIU_SUM;
1068: PETSC_EXTERN MPI_Op MPIU_MAX;
1069: PETSC_EXTERN MPI_Op MPIU_MIN;
1070: #else
1071: #define MPIU_SUM MPI_SUM
1072: #define MPIU_MAX MPI_MAX
1073: #define MPIU_MIN MPI_MIN
1074: #endif
1075: PETSC_EXTERN MPI_Op Petsc_Garbage_SetIntersectOp;
1076: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm, const PetscInt[], PetscInt *, PetscInt *);
1078: #if (defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128)) || (defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16))
1079: /*MC
1080: MPIU_SUM___FP16___FLOAT128 - MPI_Op that acts as a replacement for `MPI_SUM` with
1081: custom `MPI_Datatype` `MPIU___FLOAT128`, `MPIU___COMPLEX128`, and `MPIU___FP16`.
1083: Level: advanced
1085: Developer Note:
1086: This should be unified with `MPIU_SUM`
1088: .seealso: `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`
1089: M*/
1090: PETSC_EXTERN MPI_Op MPIU_SUM___FP16___FLOAT128;
1091: #endif
1092: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm, const PetscInt[], PetscInt *, PetscInt *);
1094: PETSC_EXTERN PetscErrorCode MPIULong_Send(void *, PetscInt, MPI_Datatype, PetscMPIInt, PetscMPIInt, MPI_Comm) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(1, 3);
1095: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void *, PetscInt, MPI_Datatype, PetscMPIInt, PetscMPIInt, MPI_Comm) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(1, 3);
1097: /*
1098: Defines PETSc error handling.
1099: */
1100: #include <petscerror.h>
1102: PETSC_EXTERN PetscBool PetscCIEnabled; /* code is running in the PETSc test harness CI */
1103: PETSC_EXTERN PetscBool PetscCIEnabledPortableErrorOutput; /* error output is stripped to ensure portability of error messages across systems */
1104: PETSC_EXTERN const char *PetscCIFilename(const char *);
1105: PETSC_EXTERN int PetscCILinenumber(int);
1107: #define PETSC_SMALLEST_CLASSID ((PetscClassId)1211211)
1108: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1109: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1110: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[], PetscClassId *);
1111: PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject, PetscObjectId *);
1112: PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject, PetscObjectId, PetscBool *);
1114: /*
1115: Routines that get memory usage information from the OS
1116: */
1117: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1118: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1119: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1120: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);
1122: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);
1124: /*
1125: Initialization of PETSc
1126: */
1127: PETSC_EXTERN PetscErrorCode PetscInitialize(int *, char ***, const char[], const char[]);
1128: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int, char **, const char[], const char[]);
1129: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1130: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1131: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1132: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1133: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1134: PETSC_EXTERN PetscErrorCode PetscGetArgs(int *, char ***);
1135: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1136: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);
1138: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1139: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);
1140: PETSC_EXTERN PetscErrorCode PetscSysFinalizePackage(void);
1142: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[], const char[]);
1143: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1144: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1145: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject, const char[]);
1147: PETSC_EXTERN PetscErrorCode PetscMonitorCompare(PetscErrorCode (*)(void), void *, PetscErrorCode (*)(void **), PetscErrorCode (*)(void), void *, PetscErrorCode (*)(void **), PetscBool *);
1149: /*
1150: These are so that in extern C code we can caste function pointers to non-extern C
1151: function pointers. Since the regular C++ code expects its function pointers to be C++
1152: */
1153: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1154: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1155: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);
1157: /*
1158: Functions that can act on any PETSc object.
1159: */
1160: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject *);
1161: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject, MPI_Comm *);
1162: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject, PetscClassId *);
1163: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject, const char *[]);
1164: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject, const char *[]);
1165: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject, const char[]);
1166: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject, const char *[]);
1167: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject, PetscInt);
1168: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject, PetscInt *);
1169: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject, PetscObject, PetscInt);
1170: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1171: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject, PetscInt *);
1172: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1173: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject, PetscMPIInt *);
1174: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject, const char[], PetscObject);
1175: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject, const char[]);
1176: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject, const char[], PetscObject *);
1177: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject, const char[], void (*)(void));
1178: #define PetscObjectComposeFunction(a, b, ...) PetscObjectComposeFunction_Private((a), (b), (PetscVoidFunction)(__VA_ARGS__))
1179: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1180: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1181: PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1182: PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject, PetscObject);
1183: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm, PetscMPIInt *);
1185: #include <petscviewertypes.h>
1186: #include <petscoptions.h>
1188: PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer, PetscBool, PetscLogDouble);
1189: PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool *);
1191: PETSC_EXTERN PetscErrorCode PetscObjectsListGetGlobalNumbering(MPI_Comm, PetscInt, PetscObject *, PetscInt *, PetscInt *);
1193: PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer, const char[]);
1194: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject, PetscViewer);
1195: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject, PetscViewer);
1196: #define PetscObjectQueryFunction(obj, name, fptr) PetscObjectQueryFunction_Private((obj), (name), (PetscVoidFunction *)(fptr))
1197: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject, const char[], void (**)(void));
1198: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject, const char[]);
1199: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject, const char[]);
1200: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject, const char[]);
1201: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject, const char *[]);
1202: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject, const char[]);
1203: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1204: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1205: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject, PetscObject, const char[]);
1206: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1207: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject, const char[], PetscBool *);
1208: PETSC_EXTERN PetscErrorCode PetscObjectObjectTypeCompare(PetscObject, PetscObject, PetscBool *);
1209: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject, const char[], PetscBool *);
1210: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject, PetscBool *, const char[], ...);
1211: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject, PetscBool *, const char[], ...);
1212: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1213: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);
1215: #if defined(PETSC_HAVE_SAWS)
1216: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1217: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1218: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject, PetscBool);
1219: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1220: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1221: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1222: PETSC_EXTERN void PetscStackSAWsGrantAccess(void);
1223: PETSC_EXTERN void PetscStackSAWsTakeAccess(void);
1224: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1225: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);
1227: #else
1228: #define PetscSAWsBlock() PETSC_SUCCESS
1229: #define PetscObjectSAWsViewOff(obj) PETSC_SUCCESS
1230: #define PetscObjectSAWsSetBlock(obj, flg) PETSC_SUCCESS
1231: #define PetscObjectSAWsBlock(obj) PETSC_SUCCESS
1232: #define PetscObjectSAWsGrantAccess(obj) PETSC_SUCCESS
1233: #define PetscObjectSAWsTakeAccess(obj) PETSC_SUCCESS
1234: #define PetscStackViewSAWs() PETSC_SUCCESS
1235: #define PetscStackSAWsViewOff() PETSC_SUCCESS
1236: #define PetscStackSAWsTakeAccess()
1237: #define PetscStackSAWsGrantAccess()
1239: #endif
1241: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[], PetscDLMode, PetscDLHandle *);
1242: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1243: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle, const char[], void **);
1244: PETSC_EXTERN PetscErrorCode PetscDLAddr(void (*)(void), char **);
1245: #ifdef PETSC_HAVE_CXX
1246: PETSC_EXTERN PetscErrorCode PetscDemangleSymbol(const char *, char **);
1247: #endif
1249: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void *, PetscStack **);
1251: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE *, PetscBool);
1252: PETSC_EXTERN PetscErrorCode PetscObjectsView(PetscViewer);
1253: PETSC_EXTERN PetscErrorCode PetscObjectsGetObject(const char *, PetscObject *, char **);
1254: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList *);
1255: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList, const char[], PetscObject *);
1256: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList, PetscObject, char **, PetscBool *);
1257: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *, const char[], PetscObject);
1258: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *, const char[]);
1259: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList, PetscObjectList *);
1261: /*
1262: Dynamic library lists. Lists of names of routines in objects or in dynamic
1263: link libraries that will be loaded as needed.
1264: */
1266: #define PetscFunctionListAdd(list, name, fptr) PetscFunctionListAdd_Private((list), (name), (PetscVoidFunction)(fptr))
1267: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList *, const char[], PetscVoidFunction);
1268: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList *);
1269: PETSC_EXTERN PetscErrorCode PetscFunctionListClear(PetscFunctionList);
1270: #define PetscFunctionListFind(list, name, fptr) PetscFunctionListFind_Private((list), (name), (PetscVoidFunction *)(fptr))
1271: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList, const char[], PetscVoidFunction *);
1272: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm, FILE *, const char[], const char[], const char[], const char[], PetscFunctionList, const char[], const char[]);
1273: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList, PetscFunctionList *);
1274: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList, PetscViewer);
1275: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList, const char ***, int *);
1276: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintNonEmpty(PetscFunctionList);
1277: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintAll(void);
1279: PETSC_EXTERN PetscDLLibrary PetscDLLibrariesLoaded;
1280: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm, PetscDLLibrary *, const char[]);
1281: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm, PetscDLLibrary *, const char[]);
1282: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm, PetscDLLibrary *, const char[], const char[], void **);
1283: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1284: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm, const char[], char *, size_t, PetscBool *);
1285: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm, const char[], PetscDLLibrary *);
1286: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);
1288: /*
1289: Useful utility routines
1290: */
1291: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm, PetscInt *, PetscInt *);
1292: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm, PetscInt, PetscInt *, PetscInt *);
1293: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm, PetscInt *, PetscInt *);
1294: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm, PetscMPIInt);
1295: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm, PetscMPIInt);
1296: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1297: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE *);
1298: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm, const PetscInt[2], PetscInt[2]);
1299: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm, const PetscReal[2], PetscReal[2]);
1301: /*MC
1302: PetscNot - negates a logical type value and returns result as a `PetscBool`
1304: Level: beginner
1306: Note:
1307: This is useful in cases like
1308: .vb
1309: int *a;
1310: PetscBool flag = PetscNot(a)
1311: .ve
1312: where !a would not return a `PetscBool` because we cannot provide a cast from int to `PetscBool` in C.
1314: .seealso: `PetscBool`, `PETSC_TRUE`, `PETSC_FALSE`
1315: M*/
1316: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)
1318: /*MC
1319: PetscHelpPrintf - Prints help messages.
1321: Synopsis:
1322: #include <petscsys.h>
1323: PetscErrorCode (*PetscHelpPrintf)(MPI_Comm comm, const char format[],args);
1325: Not Collective, only applies on MPI rank 0; No Fortran Support
1327: Input Parameters:
1328: + comm - the MPI communicator over which the help message is printed
1329: . format - the usual printf() format string
1330: - args - arguments to be printed
1332: Level: developer
1334: Note:
1335: You can change how help messages are printed by replacing the function pointer with a function that does not simply write to stdout.
1337: To use, write your own function, for example,
1338: .vb
1339: PetscErrorCode mypetschelpprintf(MPI_Comm comm,const char format[],....)
1340: {
1341: PetscFunctionReturn(PETSC_SUCCESS);
1342: }
1343: .ve
1344: then do the assignment
1345: $ PetscHelpPrintf = mypetschelpprintf;
1346: You can do the assignment before `PetscInitialize()`.
1348: The default routine used is called `PetscHelpPrintfDefault()`.
1350: .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscErrorPrintf()`, `PetscHelpPrintfDefault()`
1351: M*/
1352: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1354: /*
1355: Defines PETSc profiling.
1356: */
1357: #include <petsclog.h>
1359: /*
1360: Simple PETSc parallel IO for ASCII printing
1361: */
1362: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[], char[]);
1363: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm, const char[], const char[], FILE **);
1364: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm, FILE *);
1365: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1366: PETSC_EXTERN PetscErrorCode PetscFFlush(FILE *);
1367: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1368: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char *, size_t, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1369: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char *, size_t, const char[], size_t *, ...) PETSC_ATTRIBUTE_FORMAT(3, 5);
1370: PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[], size_t, const char *, PetscInt, const PetscReal[]);
1372: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1373: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1374: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1376: PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char *, size_t *);
1377: PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char *, char *);
1379: #if defined(PETSC_HAVE_POPEN)
1380: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm, const char[], const char[], const char[], FILE **);
1381: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm, FILE *);
1382: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1383: #endif
1385: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1386: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1387: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm, FILE *);
1388: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm, FILE *, size_t, char[]);
1389: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm, const char[], const char[], FILE **);
1390: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char *[]);
1392: PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1393: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer, void **);
1394: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer, void *);
1395: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer *);
1396: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm, PetscContainer *);
1397: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void *));
1398: PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void *);
1400: /*
1401: For use in debuggers
1402: */
1403: PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1404: PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1405: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt, const PetscInt[], PetscViewer);
1406: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt, const PetscReal[], PetscViewer);
1407: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt, const PetscScalar[], PetscViewer);
1409: /*
1410: Basic memory and string operations. These are usually simple wrappers
1411: around the basic Unix system calls, but a few of them have additional
1412: functionality and/or error checking.
1413: */
1414: #include <petscstring.h>
1416: #include <stddef.h>
1417: #include <stdlib.h>
1419: #if defined(PETSC_CLANG_STATIC_ANALYZER)
1420: #define PetscPrefetchBlock(a, b, c, d)
1421: #else
1422: /*MC
1423: PetscPrefetchBlock - Prefetches a block of memory
1425: Synopsis:
1426: #include <petscsys.h>
1427: void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)
1429: Not Collective
1431: Input Parameters:
1432: + a - pointer to first element to fetch (any type but usually `PetscInt` or `PetscScalar`)
1433: . n - number of elements to fetch
1434: . rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors)
1435: - t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note
1437: Level: developer
1439: Notes:
1440: The last two arguments (`rw` and `t`) must be compile-time constants.
1442: Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality. Not all architectures offer
1443: equivalent locality hints, but the following macros are always defined to their closest analogue.
1444: + `PETSC_PREFETCH_HINT_NTA` - Non-temporal. Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched).
1445: . `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.
1446: . `PETSC_PREFETCH_HINT_T1` - Fetch to level 2 and higher (not L1).
1447: - `PETSC_PREFETCH_HINT_T2` - Fetch to high-level cache only. (On many systems, T0 and T1 are equivalent.)
1449: This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid
1450: address).
1452: M*/
1453: #define PetscPrefetchBlock(a, n, rw, t) \
1454: do { \
1455: const char *_p = (const char *)(a), *_end = (const char *)((a) + (n)); \
1456: for (; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p, (rw), (t)); \
1457: } while (0)
1458: #endif
1459: /*
1460: Determine if some of the kernel computation routines use
1461: Fortran (rather than C) for the numerical calculations. On some machines
1462: and compilers (like complex numbers) the Fortran version of the routines
1463: is faster than the C/C++ versions. The flag --with-fortran-kernels
1464: should be used with ./configure to turn these on.
1465: */
1466: #if defined(PETSC_USE_FORTRAN_KERNELS)
1468: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1469: #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1470: #endif
1472: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
1473: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
1474: #endif
1476: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1477: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
1478: #endif
1480: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1481: #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
1482: #endif
1484: #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
1485: #define PETSC_USE_FORTRAN_KERNEL_NORM
1486: #endif
1488: #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
1489: #define PETSC_USE_FORTRAN_KERNEL_MAXPY
1490: #endif
1492: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1493: #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
1494: #endif
1496: #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1497: #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
1498: #endif
1500: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
1501: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
1502: #endif
1504: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1505: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
1506: #endif
1508: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
1509: #define PETSC_USE_FORTRAN_KERNEL_MDOT
1510: #endif
1512: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
1513: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
1514: #endif
1516: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
1517: #define PETSC_USE_FORTRAN_KERNEL_AYPX
1518: #endif
1520: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
1521: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
1522: #endif
1524: #endif
1526: /*
1527: Macros for indicating code that should be compiled with a C interface,
1528: rather than a C++ interface. Any routines that are dynamically loaded
1529: (such as the PCCreate_XXX() routines) must be wrapped so that the name
1530: mangler does not change the functions symbol name. This just hides the
1531: ugly extern "C" {} wrappers.
1532: */
1533: #if defined(__cplusplus)
1534: #define EXTERN_C_BEGIN extern "C" {
1535: #define EXTERN_C_END }
1536: #else
1537: #define EXTERN_C_BEGIN
1538: #define EXTERN_C_END
1539: #endif
1541: /*MC
1542: MPI_Comm - the basic object used by MPI to determine which processes are involved in a
1543: communication
1545: Level: beginner
1547: Note:
1548: This manual page is a place-holder because MPICH does not have a manual page for `MPI_Comm`
1550: .seealso: `PETSC_COMM_WORLD`, `PETSC_COMM_SELF`
1551: M*/
1553: #if defined(PETSC_HAVE_MPIIO)
1554: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4);
1555: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4);
1556: 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);
1557: 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);
1558: 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);
1559: 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);
1560: #endif
1562: /*@C
1563: PetscIntCast - casts a `PetscInt64` (which is 64 bits in size) to a `PetscInt` (which may be 32-bits in size), generates an
1564: error if the `PetscInt` is not large enough to hold the number.
1566: Not Collective; No Fortran Support
1568: Input Parameter:
1569: . a - the `PetscInt64` value
1571: Output Parameter:
1572: . b - the resulting `PetscInt` value
1574: Level: advanced
1576: Notes:
1577: 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
1579: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscIntMultError()`, `PetscIntSumError()`
1580: @*/
1581: static inline PetscErrorCode PetscIntCast(PetscInt64 a, PetscInt *b)
1582: {
1583: PetscFunctionBegin;
1584: *b = 0;
1585: // if using 64-bit indices already then this comparison is tautologically true
1586: 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);
1587: *b = (PetscInt)a;
1588: PetscFunctionReturn(PETSC_SUCCESS);
1589: }
1591: /*@C
1592: PetscCountCast - casts a `PetscCount` to a `PetscInt` (which may be 32-bits in size), generates an
1593: error if the `PetscInt` is not large enough to hold the number.
1595: Not Collective; No Fortran Support
1597: Input Parameter:
1598: . a - the `PetscCount` value
1600: Output Parameter:
1601: . b - the resulting `PetscInt` value
1603: Level: advanced
1605: Note:
1606: 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
1608: .seealso: `PetscCount`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscIntMultError()`, `PetscIntSumError()`, `PetscIntCast()`
1609: @*/
1610: static inline PetscErrorCode PetscCountCast(PetscCount a, PetscInt *b)
1611: {
1612: PetscFunctionBegin;
1613: *b = 0;
1614: 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);
1615: *b = (PetscInt)a;
1616: PetscFunctionReturn(PETSC_SUCCESS);
1617: }
1619: /*@C
1620: PetscBLASIntCast - casts a `PetscInt` (which may be 64-bits in size) to a `PetscBLASInt` (which may be 32-bits in size), generates an
1621: error if the `PetscBLASInt` is not large enough to hold the number.
1623: Not Collective; No Fortran Support
1625: Input Parameter:
1626: . a - the `PetscInt` value
1628: Output Parameter:
1629: . b - the resulting `PetscBLASInt` value
1631: Level: advanced
1633: Note:
1634: Errors if the integer is negative since PETSc calls to BLAS/LAPACK never need to cast negative integer inputs
1636: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscIntCast()`, `PetscCountCast()`
1637: @*/
1638: static inline PetscErrorCode PetscBLASIntCast(PetscInt a, PetscBLASInt *b)
1639: {
1640: PetscFunctionBegin;
1641: *b = 0;
1642: if (PetscDefined(USE_64BIT_INDICES) && !PetscDefined(HAVE_64BIT_BLAS_INDICES)) {
1643: 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);
1644: }
1645: PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer to BLAS/LAPACK routine");
1646: *b = (PetscBLASInt)a;
1647: PetscFunctionReturn(PETSC_SUCCESS);
1648: }
1650: /*@C
1651: PetscCuBLASIntCast - like `PetscBLASIntCast()`, but for `PetscCuBLASInt`.
1653: Not Collective; No Fortran Support
1655: Input Parameter:
1656: . a - the `PetscInt` value
1658: Output Parameter:
1659: . b - the resulting `PetscCuBLASInt` value
1661: Level: advanced
1663: Note:
1664: Errors if the integer is negative since PETSc calls to cuBLAS and friends never need to cast negative integer inputs
1666: .seealso: `PetscCuBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscIntCast()`
1667: @*/
1668: static inline PetscErrorCode PetscCuBLASIntCast(PetscInt a, PetscCuBLASInt *b)
1669: {
1670: PetscFunctionBegin;
1671: *b = 0;
1672: 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);
1673: PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer %" PetscInt_FMT "to cuBLAS routine", a);
1674: *b = (PetscCuBLASInt)a;
1675: PetscFunctionReturn(PETSC_SUCCESS);
1676: }
1678: /*@C
1679: PetscHipBLASIntCast - like `PetscBLASIntCast()`, but for `PetscHipBLASInt`.
1681: Not Collective; No Fortran Support
1683: Input Parameter:
1684: . a - the `PetscInt` value
1686: Output Parameter:
1687: . b - the resulting `PetscHipBLASInt` value
1689: Level: advanced
1691: Note:
1692: Errors if the integer is negative since PETSc calls to hipBLAS and friends never need to cast negative integer inputs
1694: .seealso: `PetscHipBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscIntCast()`
1695: @*/
1696: static inline PetscErrorCode PetscHipBLASIntCast(PetscInt a, PetscHipBLASInt *b)
1697: {
1698: PetscFunctionBegin;
1699: *b = 0;
1700: 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);
1701: PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer %" PetscInt_FMT "to hipBLAS routine", a);
1702: *b = (PetscHipBLASInt)a;
1703: PetscFunctionReturn(PETSC_SUCCESS);
1704: }
1706: /*@C
1707: PetscMPIIntCast - casts a `PetscInt` (which may be 64-bits in size) to a PetscMPIInt (which may be 32-bits in size), generates an
1708: error if the `PetscMPIInt` is not large enough to hold the number.
1710: Not Collective; No Fortran Support
1712: Input Parameter:
1713: . a - the `PetscInt` value
1715: Output Parameter:
1716: . b - the resulting `PetscMPIInt` value
1718: Level: advanced
1720: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscIntCast()`
1721: @*/
1722: static inline PetscErrorCode PetscMPIIntCast(PetscInt a, PetscMPIInt *b)
1723: {
1724: PetscFunctionBegin;
1725: *b = 0;
1726: 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);
1727: *b = (PetscMPIInt)a;
1728: PetscFunctionReturn(PETSC_SUCCESS);
1729: }
1731: #define PetscInt64Mult(a, b) (((PetscInt64)(a)) * ((PetscInt64)(b)))
1733: /*@C
1734: PetscRealIntMultTruncate - Computes the product of a positive `PetscReal` and a positive
1735: `PetscInt` and truncates the value to slightly less than the maximal possible value.
1737: Not Collective; No Fortran Support
1739: Input Parameters:
1740: + a - The `PetscReal` value
1741: - b - The `PetscInt` value
1743: Level: advanced
1745: Notes:
1746: Returns the result as a `PetscInt` value.
1748: Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`.
1750: Use `PetscIntMultTruncate()` to compute the product of two positive `PetscInt` and truncate
1751: to fit a `PetscInt`.
1753: Use `PetscIntMultError()` to compute the product of two `PetscInt` if you wish to generate an
1754: error if the result will not fit in a `PetscInt`.
1756: Developers Notes:
1757: We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but
1758: requires many more checks.
1760: This is used where we compute approximate sizes for workspace and need to insure the
1761: workspace is index-able.
1763: .seealso: `PetscReal`, `PetscInt`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()`
1764: @*/
1765: static inline PetscInt PetscRealIntMultTruncate(PetscReal a, PetscInt b)
1766: {
1767: PetscInt64 r = (PetscInt64)(a * (PetscReal)b);
1768: if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
1769: return (PetscInt)r;
1770: }
1772: /*@C
1774: PetscIntMultTruncate - Computes the product of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value
1776: Not Collective; No Fortran Support
1778: Input Parameters:
1779: + a - the PetscInt value
1780: - b - the second value
1782: Returns:
1783: the result as a `PetscInt` value
1785: Level: advanced
1787: Notes:
1788: Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`
1790: Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt`
1792: 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`
1794: Developers Notes:
1795: We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but requires many more checks.
1797: This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
1799: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()`,
1800: `PetscIntSumTruncate()`
1801: @*/
1802: static inline PetscInt PetscIntMultTruncate(PetscInt a, PetscInt b)
1803: {
1804: PetscInt64 r = PetscInt64Mult(a, b);
1805: if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
1806: return (PetscInt)r;
1807: }
1809: /*@C
1811: PetscIntSumTruncate - Computes the sum of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value
1813: Not Collective; No Fortran Support
1815: Input Parameters:
1816: + a - the `PetscInt` value
1817: - b - the second value
1819: Returns:
1820: the result as a `PetscInt` value
1822: Level: advanced
1824: Notes:
1825: Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`
1827: Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt`
1829: 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`
1831: Developers Notes:
1832: This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
1834: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`
1835: @*/
1836: static inline PetscInt PetscIntSumTruncate(PetscInt a, PetscInt b)
1837: {
1838: PetscInt64 r = ((PetscInt64)a) + ((PetscInt64)b);
1839: if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
1840: return (PetscInt)r;
1841: }
1843: /*@C
1845: PetscIntMultError - Computes the product of two positive `PetscInt` and generates an error with overflow.
1847: Not Collective; No Fortran Support
1849: Input Parameters:
1850: + a - the `PetscInt` value
1851: - b - the second value
1853: Output Parameter:
1854: . result - the result as a `PetscInt` value, or `NULL` if you do not want the result, you just want to check if it overflows
1856: Level: advanced
1858: Notes:
1859: Use `PetscInt64Mult()` to compute the product of two `PetscInt` and store in a `PetscInt64`
1861: Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt`
1863: Developers Note:
1864: 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.
1865: `PetscIntSumError()` can be used to check for this situation.
1867: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscIntMult64()`, `PetscIntSumError()`
1868: @*/
1869: static inline PetscErrorCode PetscIntMultError(PetscInt a, PetscInt b, PetscInt *result)
1870: {
1871: PetscInt64 r = PetscInt64Mult(a, b);
1873: PetscFunctionBegin;
1874: if (result) *result = (PetscInt)r;
1875: if (!PetscDefined(USE_64BIT_INDICES)) {
1876: 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);
1877: }
1878: PetscFunctionReturn(PETSC_SUCCESS);
1879: }
1881: /*@C
1883: PetscIntSumError - Computes the sum of two positive `PetscInt` and generates an error with overflow.
1885: Not Collective; No Fortran Support
1887: Input Parameters:
1888: + a - the `PetscInt` value
1889: - b - the second value
1891: Output Parameter:
1892: . c - the result as a `PetscInt` value, or `NULL` if you do not want the result, you just want to check if it overflows
1894: Level: advanced
1896: Notes:
1897: Use `PetscInt64Mult()` to compute the product of two 32-bit `PetscInt` and store in a `PetscInt64`
1899: Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt`
1901: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`
1902: @*/
1903: static inline PetscErrorCode PetscIntSumError(PetscInt a, PetscInt b, PetscInt *result)
1904: {
1905: PetscInt64 r = ((PetscInt64)a) + ((PetscInt64)b);
1907: PetscFunctionBegin;
1908: if (result) *result = (PetscInt)r;
1909: if (!PetscDefined(USE_64BIT_INDICES)) {
1910: 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);
1911: }
1912: PetscFunctionReturn(PETSC_SUCCESS);
1913: }
1915: /*
1916: The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
1917: */
1918: #if defined(hz)
1919: #undef hz
1920: #endif
1922: #if defined(PETSC_HAVE_SYS_TYPES_H)
1923: #include <sys/types.h>
1924: #endif
1926: /*MC
1928: PETSC_VERSION - This manual page provides information about how PETSc documents and uses its version information. This information is available to both C/C++
1929: and Fortran compilers when `petscsys.h` is included.
1931: The current PETSc version and the API for accessing it are defined in <A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/include/petscversion.h.html">include/petscversion.html</A>
1933: The complete version number is given as the triple PETSC_VERSION_MAJOR.PETSC_VERSION_MINOR.PETSC_VERSION_SUBMINOR (in short hand x.y.z)
1935: 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
1936: where only a change in the major version number (x) indicates a change in the API.
1938: 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
1940: 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),
1941: 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
1942: version number (x.y.z).
1944: `PETSC_RELEASE_DATE` is the date the x.y version was released (i.e. the version before any patch releases)
1946: `PETSC_VERSION_DATE` is the date the x.y.z version was released
1948: `PETSC_VERSION_GIT` is the last git commit to the repository given in the form vx.y.z-wwwww
1950: `PETSC_VERSION_DATE_GIT` is the date of the last git commit to the repository
1952: `PETSC_VERSION_()` is deprecated and will eventually be removed.
1954: Level: intermediate
1955: M*/
1957: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[], size_t);
1958: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[], size_t);
1959: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[], size_t);
1960: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[], size_t);
1961: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
1962: PETSC_EXTERN PetscErrorCode PetscGetDate(char[], size_t);
1963: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
1964: PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt *, PetscInt *, PetscInt *, PetscInt *);
1966: PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt, const PetscInt[], PetscBool *);
1967: PETSC_EXTERN PetscErrorCode PetscSortedInt64(PetscInt, const PetscInt64[], PetscBool *);
1968: PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt, const PetscMPIInt[], PetscBool *);
1969: PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt, const PetscReal[], PetscBool *);
1970: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt, PetscInt[]);
1971: PETSC_EXTERN PetscErrorCode PetscSortInt64(PetscInt, PetscInt64[]);
1972: PETSC_EXTERN PetscErrorCode PetscSortCount(PetscInt, PetscCount[]);
1973: PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt, PetscInt[]);
1974: PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *, PetscInt[]);
1975: PETSC_EXTERN PetscErrorCode PetscSortedCheckDupsInt(PetscInt, const PetscInt[], PetscBool *);
1976: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt *, PetscInt[]);
1977: PETSC_EXTERN PetscErrorCode PetscCheckDupsInt(PetscInt, const PetscInt[], PetscBool *);
1978: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt *);
1979: PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt *);
1980: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt, const PetscInt[], PetscInt[]);
1981: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt, const char *[], PetscInt[]);
1982: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt, PetscInt[], PetscInt[]);
1983: PETSC_EXTERN PetscErrorCode PetscSortIntWithCountArray(PetscCount, PetscInt[], PetscCount[]);
1984: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt, PetscInt[], PetscInt[], PetscInt[]);
1985: PETSC_EXTERN PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount, PetscInt[], PetscInt[], PetscCount[]);
1986: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt, PetscMPIInt[]);
1987: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *, PetscMPIInt[]);
1988: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt, PetscMPIInt[], PetscMPIInt[]);
1989: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt, PetscMPIInt[], PetscInt[]);
1990: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt, PetscInt[], PetscScalar[]);
1991: PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt, PetscInt[], void *, size_t, void *);
1992: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt, PetscReal[]);
1993: PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt, PetscReal[], PetscInt[]);
1994: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt, const PetscReal[], PetscInt[]);
1995: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt *, PetscReal[]);
1996: PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal, PetscInt, const PetscReal[], PetscReal, PetscInt *);
1997: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt, PetscInt, PetscScalar[], PetscInt[]);
1998: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt, PetscInt, PetscReal[], PetscInt[]);
1999: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt, const PetscBool[], const PetscInt[], PetscInt *, PetscInt **, PetscInt **, PetscInt **, PetscInt **);
2000: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt, const PetscInt[], const PetscInt[], PetscInt, const PetscInt[], const PetscInt[], PetscInt *, PetscInt **, PetscInt **);
2001: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscInt *, PetscInt **);
2002: PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt, const PetscMPIInt[], PetscInt, const PetscMPIInt[], PetscInt *, PetscMPIInt **);
2003: PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *);
2005: PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt, void *, size_t, int (*)(const void *, const void *, void *), void *);
2006: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt, PetscInt[]);
2007: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt, PetscMPIInt[]);
2008: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt, PetscReal[]);
2009: PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt, void *, size_t, void *, size_t, int (*)(const void *, const void *, void *), void *);
2010: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt, PetscInt[], PetscInt[]);
2011: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt, PetscMPIInt[], PetscMPIInt[]);
2012: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt, PetscReal[], PetscInt[]);
2014: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2015: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[], size_t);
2017: /*J
2018: PetscRandomType - String with the name of a PETSc randomizer
2020: Level: beginner
2022: Note:
2023: To use `PETSCSPRNG` or `PETSCRANDOM123` you must have ./configure PETSc
2024: with the option `--download-sprng` or `--download-random123`. We recommend the default provided with PETSc.
2026: .seealso: `PetscRandomSetType()`, `PetscRandom`, `PetscRandomCreate()`
2027: J*/
2028: typedef const char *PetscRandomType;
2029: #define PETSCRAND "rand"
2030: #define PETSCRAND48 "rand48"
2031: #define PETSCSPRNG "sprng"
2032: #define PETSCRANDER48 "rander48"
2033: #define PETSCRANDOM123 "random123"
2034: #define PETSCCURAND "curand"
2036: /* Logging support */
2037: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;
2039: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);
2040: PETSC_EXTERN PetscErrorCode PetscRandomFinalizePackage(void);
2042: /* Dynamic creation and loading functions */
2043: PETSC_EXTERN PetscFunctionList PetscRandomList;
2045: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[], PetscErrorCode (*)(PetscRandom));
2046: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2047: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2048: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType *);
2049: PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom, PetscObject, const char[]);
2050: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom, PetscViewer);
2052: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm, PetscRandom *);
2053: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom, PetscScalar *);
2054: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom, PetscReal *);
2055: PETSC_EXTERN PetscErrorCode PetscRandomGetValues(PetscRandom, PetscInt, PetscScalar *);
2056: PETSC_EXTERN PetscErrorCode PetscRandomGetValuesReal(PetscRandom, PetscInt, PetscReal *);
2057: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom, PetscScalar *, PetscScalar *);
2058: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom, PetscScalar, PetscScalar);
2059: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom, unsigned long);
2060: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom, unsigned long *);
2061: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2062: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom *);
2064: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[], char[], size_t);
2065: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[], char[], size_t);
2066: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[], size_t);
2067: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[], char[]);
2068: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[], size_t);
2069: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[], char, PetscBool *);
2070: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[], char, PetscBool *);
2071: PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2072: PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2073: PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);
2075: static inline PetscBool PetscBinaryBigEndian(void)
2076: {
2077: long _petsc_v = 1;
2078: return ((char *)&_petsc_v)[0] ? PETSC_FALSE : PETSC_TRUE;
2079: }
2081: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int, void *, PetscInt, PetscInt *, PetscDataType);
2082: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm, int, void *, PetscInt, PetscInt *, PetscDataType);
2083: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int, const void *, PetscInt, PetscDataType);
2084: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm, int, const void *, PetscInt, PetscDataType);
2085: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[], PetscFileMode, int *);
2086: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2087: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm, PetscBool *);
2088: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm, PetscBool *);
2089: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm, char[], size_t);
2090: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm, const char[], char[], size_t, PetscBool *);
2091: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm, const char[], char[], size_t, PetscBool *);
2092: #if defined(PETSC_USE_SOCKET_VIEWER)
2093: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[], int, int *);
2094: #endif
2096: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int, off_t, PetscBinarySeekType, off_t *);
2097: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm, int, off_t, PetscBinarySeekType, off_t *);
2098: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *, PetscDataType, PetscInt);
2100: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2101: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[], PetscBool);
2102: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2103: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char *);
2104: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2105: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2106: PETSC_EXTERN PetscErrorCode PetscWaitOnError(void);
2108: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt *);
2109: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], PetscMPIInt **, PetscMPIInt **);
2110: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt **, PetscMPIInt **, PetscMPIInt **);
2111: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscInt ***, MPI_Request **);
2112: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscScalar ***, MPI_Request **);
2113: 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);
2114: 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);
2115: 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);
2117: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm, PetscBuildTwoSidedType);
2118: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm, PetscBuildTwoSidedType *);
2120: PETSC_EXTERN PetscErrorCode PetscSSEIsEnabled(MPI_Comm, PetscBool *, PetscBool *);
2122: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);
2124: struct _n_PetscSubcomm {
2125: MPI_Comm parent; /* parent communicator */
2126: MPI_Comm dupparent; /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2127: MPI_Comm child; /* the sub-communicator */
2128: PetscMPIInt n; /* num of subcommunicators under the parent communicator */
2129: PetscMPIInt color; /* color of processors belong to this communicator */
2130: PetscMPIInt *subsize; /* size of subcommunicator[color] */
2131: PetscSubcommType type;
2132: char *subcommprefix;
2133: };
2135: static inline MPI_Comm PetscSubcommParent(PetscSubcomm scomm)
2136: {
2137: return scomm->parent;
2138: }
2139: static inline MPI_Comm PetscSubcommChild(PetscSubcomm scomm)
2140: {
2141: return scomm->child;
2142: }
2143: static inline MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm)
2144: {
2145: return scomm->dupparent;
2146: }
2147: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm, PetscSubcomm *);
2148: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm *);
2149: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm, PetscInt);
2150: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm, PetscSubcommType);
2151: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm, PetscMPIInt, PetscMPIInt);
2152: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm, PetscViewer);
2153: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2154: PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm, const char[]);
2155: PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm, MPI_Comm *);
2156: PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm, MPI_Comm *);
2157: PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm, MPI_Comm *);
2159: PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt, PetscHeap *);
2160: PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap, PetscInt, PetscInt);
2161: PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap, PetscInt *, PetscInt *);
2162: PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap, PetscInt *, PetscInt *);
2163: PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap, PetscInt, PetscInt);
2164: PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2165: PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap *);
2166: PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap, PetscViewer);
2168: PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2169: PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm, PetscShmComm *);
2170: PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm, PetscMPIInt, PetscMPIInt *);
2171: PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm, PetscMPIInt, PetscMPIInt *);
2172: PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm, MPI_Comm *);
2174: /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2175: PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm, PetscInt, PetscOmpCtrl *);
2176: PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl, MPI_Comm *, MPI_Comm *, PetscBool *);
2177: PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl *);
2178: PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2179: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2180: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);
2182: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t, size_t, PetscSegBuffer *);
2183: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer *);
2184: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer, size_t, void *);
2185: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer, void *);
2186: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer, void *);
2187: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer, void *);
2188: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer, size_t *);
2189: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer, size_t);
2191: /* Type-safe wrapper to encourage use of PETSC_RESTRICT. Does not use PetscFunctionBegin because the error handling
2192: * prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as
2193: * possible. */
2194: static inline PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg, size_t count, PetscInt *PETSC_RESTRICT *slot)
2195: {
2196: return PetscSegBufferGet(seg, count, (void **)slot);
2197: }
2199: extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2200: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted *);
2201: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted *);
2202: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted, const char *, const char *, PetscBool *);
2204: #include <stdarg.h>
2205: PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char *, size_t, const char[], size_t *, va_list);
2206: PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE *, const char[], va_list);
2208: PETSC_EXTERN PetscSegBuffer PetscCitationsList;
2210: /*@C
2211: PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code.
2213: Not Collective; No Fortran Support
2215: Input Parameters:
2216: + cite - the bibtex item, formatted to displayed on multiple lines nicely
2217: - set - a boolean variable initially set to `PETSC_FALSE`; this is used to insure only a single registration of the citation
2219: Options Database Key:
2220: . -citations [filename] - print out the bibtex entries for the given computation
2222: Level: intermediate
2223: @*/
2224: static inline PetscErrorCode PetscCitationsRegister(const char cit[], PetscBool *set)
2225: {
2226: size_t len;
2227: char *vstring;
2229: PetscFunctionBegin;
2230: if (set && *set) PetscFunctionReturn(PETSC_SUCCESS);
2231: PetscCall(PetscStrlen(cit, &len));
2232: PetscCall(PetscSegBufferGet(PetscCitationsList, len, &vstring));
2233: PetscCall(PetscArraycpy(vstring, cit, len));
2234: if (set) *set = PETSC_TRUE;
2235: PetscFunctionReturn(PETSC_SUCCESS);
2236: }
2238: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm, char[], char[], size_t);
2239: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm, const char[], char[], size_t);
2240: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm, const char[], const char[]);
2242: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm, char[], char[], size_t);
2243: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm, const char[], char[], char[], size_t);
2244: PETSC_EXTERN PetscErrorCode PetscBoxUpload(MPI_Comm, const char[], const char[]);
2246: PETSC_EXTERN PetscErrorCode PetscGlobusGetTransfers(MPI_Comm, const char[], char[], size_t);
2247: PETSC_EXTERN PetscErrorCode PetscGlobusAuthorize(MPI_Comm, char[], size_t);
2248: PETSC_EXTERN PetscErrorCode PetscGlobusUpload(MPI_Comm, const char[], const char[]);
2250: PETSC_EXTERN PetscErrorCode PetscPullJSONValue(const char[], const char[], char[], size_t, PetscBool *);
2251: PETSC_EXTERN PetscErrorCode PetscPushJSONValue(char[], const char[], const char[], size_t);
2253: #if defined(PETSC_USE_DEBUG)
2254: static inline unsigned int PetscStrHash(const char *str)
2255: {
2256: unsigned int c, hash = 5381;
2258: while ((c = (unsigned int)*str++)) hash = ((hash << 5) + hash) + c; /* hash * 33 + c */
2259: return hash;
2260: }
2262: /*MC
2263: MPIU_Allreduce - a PETSc replacement for `MPI_Allreduce()` that tries to determine if the call from all the MPI ranks occur from the
2264: same place in the PETSc code. This helps to detect bugs where different MPI ranks follow different code paths
2265: resulting in inconsistent and incorrect calls to `MPI_Allreduce()`.
2267: Synopsis:
2268: #include <petscsys.h>
2269: PetscErrorCode MPIU_Allreduce(void *indata,void *outdata,PetscMPIInt count,MPI_Datatype datatype, MPI_Op op, MPI_Comm comm);
2271: Collective
2273: Input Parameters:
2274: + a - pointer to the input data to be reduced
2275: . c - the number of MPI data items in a and b
2276: . d - the MPI datatype, for example `MPI_INT`
2277: . e - the MPI operation, for example `MPI_SUM`
2278: - fcomm - the MPI communicator on which the operation occurs
2280: Output Parameter:
2281: . b - the reduced values
2283: Level: developer
2285: Notes:
2286: In optimized mode this directly calls `MPI_Allreduce()`
2288: This is defined as a macro that can return error codes internally so it cannot be used in a subroutine that returns void.
2290: 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
2292: .seealso: `MPI_Allreduce()`
2293: M*/
2294: #define MPIU_Allreduce(a, b, c, d, e, fcomm) \
2295: PetscMacroReturnStandard( \
2296: PetscMPIInt a_b1[6], a_b2[6]; \
2297: int _mpiu_allreduce_c_int = (int)(c); \
2298: a_b1[0] = -(PetscMPIInt)__LINE__; \
2299: a_b1[1] = -a_b1[0]; \
2300: a_b1[2] = -(PetscMPIInt)PetscStrHash(PETSC_FUNCTION_NAME); \
2301: a_b1[3] = -a_b1[2]; \
2302: a_b1[4] = -(PetscMPIInt)(c); \
2303: a_b1[5] = -a_b1[4]; \
2304: \
2305: PetscCallMPI(MPI_Allreduce(a_b1, a_b2, 6, MPI_INT, MPI_MAX, fcomm)); \
2306: PetscCheck(-a_b2[0] == a_b2[1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "MPI_Allreduce() called in different locations (code lines) on different processors"); \
2307: PetscCheck(-a_b2[2] == a_b2[3], PETSC_COMM_SELF, PETSC_ERR_PLIB, "MPI_Allreduce() called in different locations (functions) on different processors"); \
2308: 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); \
2309: PetscCallMPI(MPI_Allreduce((a), (b), (c), (d), (e), (fcomm)));)
2310: #else
2311: #define MPIU_Allreduce(a, b, c, d, e, fcomm) PetscMacroReturnStandard(PetscCallMPI(MPI_Allreduce((a), (b), (c), (d), (e), (fcomm))))
2312: #endif
2314: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
2315: PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint, PetscMPIInt, MPI_Info, MPI_Comm, void *, MPI_Win *);
2316: PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win, PetscMPIInt, MPI_Aint *, PetscMPIInt *, void *);
2317: #endif
2319: /* this is a vile hack */
2320: #if defined(PETSC_HAVE_NECMPI)
2321: #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)
2322: #define MPI_Type_free(a) (*(a) = MPI_DATATYPE_NULL, 0);
2323: #endif
2324: #endif
2326: /*
2327: List of external packages and queries on it
2328: */
2329: PETSC_EXTERN PetscErrorCode PetscHasExternalPackage(const char[], PetscBool *);
2331: /* this cannot go here because it may be in a different shared library */
2332: PETSC_EXTERN PetscErrorCode PCMPIServerBegin(void);
2333: PETSC_EXTERN PetscErrorCode PCMPIServerEnd(void);
2334: PETSC_EXTERN PetscErrorCode PCMPICommsDestroy(void);
2335: PETSC_EXTERN PetscBool PCMPIServerActive;
2337: #define PETSC_HAVE_FORTRAN PETSC_DEPRECATED_MACRO(3, 20, 0, "PETSC_USE_FORTRAN_BINDINGS", ) PETSC_USE_FORTRAN_BINDINGS