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: /*MC
10: PeOP - indicates an argument to a PETSc function is optional and one can pass `NULL` instead. This is used by the Fortran API generator
12: Level: developer
14: Example:
15: .vb
16: PetscErrorCode XXXX(Vec v, PeOp PetscObject obj, PeOp PetscInt *idx, PeOp PetscInt *array[])
17: .ve
19: Notes:
20: This is not part of the PETSc public API and should only be used in PETSc source code.
22: Put this in the function declaration in front of each variable that is optional
24: Developer Note:
25: Shortened form of PETSc optional
27: .seealso: `PeNS`, `PeNSS`, `PetscCtxRt`, `PetscInitialize()`
28: M*/
29: #define PeOp
31: /*MC
32: PeNS - indicates a function that does not use the PETSc standard arguments which make it easy to generate automatic language stubs for other languages
34: Level: developer
36: Notes:
37: This is not part of the PETSc public API and should only be used in PETSc source code.
39: Put this at the end of the function declaration closing parenthesis
41: Developer Note:
42: Shortened form of PETSc non-standard
44: .seealso: `PeOp`, `PeNSS`, `PetscCtxRt`, `PetscInitialize()`
45: M*/
46: #define PeNS
48: /*MC
49: PeNSS - indicates a function that needs a special treatment in the C-side stub when generating the binding for other languages
51: Level: developer
53: Notes:
54: This is not part of the PETSc public API and should only be used in PETSc source code.
56: Put this at the end of the function declaration closing parenthesis
58: It is similar to PeNS; in Fortran it will generate the Fortran interface definition automatically but not the C stub, which should be added manually under the appropriate `ftn-custom` directory
60: Developer Note:
61: Shortened form of PETSc non-standard stub
63: .seealso: `PeOp`, `PeNS`, `PetscCtxRt`, `PetscInitialize()`
64: M*/
65: #define PeNSS
67: /* ========================================================================== */
68: /*
69: petscconf.h is contained in ${PETSC_ARCH}/include/petscconf.h it is
70: found automatically by the compiler due to the -I${PETSC_DIR}/${PETSC_ARCH}/include that
71: PETSc's makefiles add to the compiler rules.
72: For --prefix installs the directory ${PETSC_ARCH} does not exist and petscconf.h is in the same
73: directory as the other PETSc include files.
74: */
75: #include <petscconf.h>
76: #include <petscpkg_version.h>
77: #include <petscconf_poison.h>
78: #include <petscfix.h>
79: #include <petscmacros.h>
81: /* SUBMANSEC = Sys */
83: #if defined(PETSC_DESIRE_FEATURE_TEST_MACROS)
84: /*
85: Feature test macros must be included before headers defined by IEEE Std 1003.1-2001
86: We only turn these in PETSc source files that require them by setting PETSC_DESIRE_FEATURE_TEST_MACROS
87: */
88: #if defined(PETSC__POSIX_C_SOURCE_200112L) && !defined(_POSIX_C_SOURCE)
89: #define _POSIX_C_SOURCE 200112L
90: #endif
91: #if defined(PETSC__BSD_SOURCE) && !defined(_BSD_SOURCE)
92: #define _BSD_SOURCE
93: #endif
94: #if defined(PETSC__DEFAULT_SOURCE) && !defined(_DEFAULT_SOURCE)
95: #define _DEFAULT_SOURCE
96: #endif
97: #if defined(PETSC__GNU_SOURCE) && !defined(_GNU_SOURCE)
98: #define _GNU_SOURCE
99: #endif
100: #endif
102: #include <petscsystypes.h>
104: /* ========================================================================== */
106: /*
107: Defines the interface to MPI allowing the use of all MPI functions.
109: PETSc does not use the C++ binding of MPI at ALL. The following flag
110: makes sure the C++ bindings are not included. The C++ bindings REQUIRE
111: putting mpi.h before ANY C++ include files, we cannot control this
112: with all PETSc users. Users who want to use the MPI C++ bindings can include
113: mpicxx.h directly in their code
114: */
115: #if !defined(MPICH_SKIP_MPICXX)
116: #define MPICH_SKIP_MPICXX 1
117: #endif
118: #if !defined(OMPI_SKIP_MPICXX)
119: #define OMPI_SKIP_MPICXX 1
120: #endif
121: #if defined(PETSC_HAVE_MPIUNI)
122: #include <petsc/mpiuni/mpi.h>
123: #else
124: #include <mpi.h>
125: #endif
127: /*
128: Perform various sanity checks that the correct mpi.h is being included at compile time.
129: This usually happens because
130: * either an unexpected mpi.h is in the default compiler path (i.e. in /usr/include) or
131: * an extra include path -I/something (which contains the unexpected mpi.h) is being passed to the compiler
132: Note: with MPICH and OpenMPI, accept versions [x.y.z, x+1.0.0) as compatible
133: */
134: #if defined(PETSC_HAVE_MPIUNI)
135: #if !defined(MPIUNI_H)
136: #error "PETSc was configured with --with-mpi=0 but now appears to be compiling using a different mpi.h"
137: #endif
138: #elif defined(PETSC_HAVE_I_MPI)
139: #if !defined(I_MPI_NUMVERSION)
140: #error "PETSc was configured with I_MPI but now appears to be compiling using a non-I_MPI mpi.h"
141: #elif I_MPI_NUMVERSION != PETSC_PKG_I_MPI_NUMVERSION
142: #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"
143: #endif
144: #elif defined(PETSC_HAVE_MVAPICH2)
145: #if !defined(MVAPICH2_NUMVERSION)
146: #error "PETSc was configured with MVAPICH2 but now appears to be compiling using a non-MVAPICH2 mpi.h"
147: #elif MVAPICH2_NUMVERSION != PETSC_PKG_MVAPICH2_NUMVERSION
148: #error "PETSc was configured with one MVAPICH2 mpi.h version but now appears to be compiling using a different MVAPICH2 mpi.h version"
149: #endif
150: #elif defined(PETSC_HAVE_MPICH)
151: #if !defined(MPICH_NUMVERSION) || defined(MVAPICH2_NUMVERSION) || defined(I_MPI_NUMVERSION)
152: #error "PETSc was configured with MPICH but now appears to be compiling using a non-MPICH mpi.h"
153: #elif PETSC_PKG_MPICH_VERSION_GT(MPICH_NUMVERSION / 10000000, MPICH_NUMVERSION / 100000 % 100, MPICH_NUMVERSION / 1000 % 100)
154: #error "PETSc was configured with one MPICH mpi.h version but now appears to be compiling using an older MPICH mpi.h version"
155: #elif PETSC_PKG_MPICH_VERSION_LT(MPICH_NUMVERSION / 10000000, 0, 0)
156: #error "PETSc was configured with one MPICH mpi.h version but now appears to be compiling using a newer major MPICH mpi.h version"
157: #endif
158: #elif defined(PETSC_HAVE_OPENMPI)
159: #if !defined(OMPI_MAJOR_VERSION)
160: #error "PETSc was configured with Open MPI but now appears to be compiling using a non-Open MPI mpi.h"
161: #elif PETSC_PKG_OPENMPI_VERSION_GT(OMPI_MAJOR_VERSION, OMPI_MINOR_VERSION, OMPI_RELEASE_VERSION)
162: #error "PETSc was configured with one Open MPI mpi.h version but now appears to be compiling using an older Open MPI mpi.h version"
163: #elif PETSC_PKG_OPENMPI_VERSION_LT(OMPI_MAJOR_VERSION, 0, 0)
164: #error "PETSc was configured with one Open MPI mpi.h version but now appears to be compiling using a newer major Open MPI mpi.h version"
165: #endif
166: #elif defined(PETSC_HAVE_MSMPI_VERSION)
167: #if !defined(MSMPI_VER)
168: #error "PETSc was configured with MSMPI but now appears to be compiling using a non-MSMPI mpi.h"
169: #elif (MSMPI_VER != PETSC_HAVE_MSMPI_VERSION)
170: #error "PETSc was configured with one MSMPI mpi.h version but now appears to be compiling using a different MSMPI mpi.h version"
171: #endif
172: #elif defined(OMPI_MAJOR_VERSION) || defined(MPICH_NUMVERSION) || defined(MSMPI_VER)
173: #error "PETSc was configured with undetermined MPI - but now appears to be compiling using any of Open MPI, MS-MPI or a MPICH variant"
174: #endif
176: /*
177: Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler
178: see the top of mpicxx.h in the MPICH2 distribution.
179: */
180: #include <stdio.h>
182: /* MSMPI on 32-bit Microsoft Windows requires this yukky hack - that breaks MPI standard compliance */
183: #if !defined(MPIAPI)
184: #define MPIAPI
185: #endif
187: PETSC_EXTERN MPI_Datatype MPIU_ENUM PETSC_ATTRIBUTE_MPI_TYPE_TAG(PetscEnum);
188: #define MPIU_BOOL MPI_C_BOOL PETSC_DEPRECATED_MACRO(3, 24, 0, "MPI_C_BOOL", )
190: /*MC
191: MPIU_INT - Portable MPI datatype corresponding to `PetscInt` independent of the precision of `PetscInt`
193: Level: beginner
195: Note:
196: In MPI calls that require an MPI datatype that matches a `PetscInt` or array of `PetscInt` values, pass this value.
198: .seealso: `PetscReal`, `PetscScalar`, `PetscComplex`, `PetscInt`, `MPIU_COUNT`, `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`
199: M*/
201: PETSC_EXTERN MPI_Datatype MPIU_FORTRANADDR;
203: #if defined(PETSC_USE_64BIT_INDICES)
204: #define MPIU_INT MPIU_INT64
205: #else
206: #define MPIU_INT MPI_INT
207: #endif
209: /*MC
210: MPIU_COUNT - Portable MPI datatype corresponding to `PetscCount` independent of the precision of `PetscCount`
212: Level: beginner
214: Note:
215: In MPI calls that require an MPI datatype that matches a `PetscCount` or array of `PetscCount` values, pass this value.
217: Developer Note:
218: It seems `MPI_AINT` is unsigned so this may be the wrong choice here since `PetscCount` is signed
220: .seealso: `PetscReal`, `PetscScalar`, `PetscComplex`, `PetscInt`, `MPIU_INT`, `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`
221: M*/
222: #define MPIU_COUNT MPI_AINT
224: /*
225: For the rare cases when one needs to send a size_t object with MPI
226: */
227: PETSC_EXTERN MPI_Datatype MPIU_SIZE_T PETSC_ATTRIBUTE_MPI_TYPE_TAG(size_t);
229: /*
230: You can use PETSC_STDOUT as a replacement of stdout. You can also change
231: the value of PETSC_STDOUT to redirect all standard output elsewhere
232: */
233: PETSC_EXTERN FILE *PETSC_STDOUT;
235: /*
236: You can use PETSC_STDERR as a replacement of stderr. You can also change
237: the value of PETSC_STDERR to redirect all standard error elsewhere
238: */
239: PETSC_EXTERN FILE *PETSC_STDERR;
241: /*
242: Handle inclusion when using clang compiler with CUDA support
243: __float128 is not available for the device
244: */
245: #if defined(__clang__) && (defined(__CUDA_ARCH__) || defined(__HIPCC__))
246: #define PETSC_SKIP_REAL___FLOAT128
247: #endif
249: /*
250: Declare extern C stuff after including external header files
251: */
253: PETSC_EXTERN PetscBool PETSC_RUNNING_ON_VALGRIND;
254: /*
255: Defines elementary mathematics functions and constants.
256: */
257: #include <petscmath.h>
259: /*MC
260: PETSC_IGNORE - same as `NULL`, means PETSc will ignore this argument
262: Level: beginner
264: Note:
265: Accepted by many PETSc functions to not set a parameter and instead use a default value
267: Fortran Note:
268: Use `PETSC_NULL_INTEGER`, `PETSC_NULL_SCALAR` etc
270: .seealso: `PETSC_DECIDE`, `PETSC_DEFAULT`, `PETSC_DETERMINE`
271: M*/
272: #define PETSC_IGNORE PETSC_NULLPTR
273: #define PETSC_NULL PETSC_DEPRECATED_MACRO(3, 19, 0, "PETSC_NULLPTR", ) PETSC_NULLPTR
275: /*MC
276: PETSC_UNLIMITED - standard way of passing an integer or floating point parameter to indicate PETSc there is no bound on the value allowed
278: Level: beginner
280: Example Usage:
281: .vb
282: KSPSetTolerances(ksp, PETSC_CURRENT, PETSC_CURRENT, PETSC_UNLIMITED, PETSC_UNLIMITED);
283: .ve
284: indicates that the solver is allowed to take any number of iterations and will not stop early no matter how the residual gets.
286: Fortran Note:
287: Use `PETSC_UNLIMITED_INTEGER` or `PETSC_UNLIMITED_REAL`.
289: .seealso: `PETSC_DEFAULT`, `PETSC_IGNORE`, `PETSC_DETERMINE`, `PETSC_DECIDE`
290: M*/
292: /*MC
293: PETSC_DECIDE - standard way of passing an integer or floating point parameter to indicate PETSc should determine an appropriate value
295: Level: beginner
297: Example Usage:
298: .vb
299: VecSetSizes(ksp, PETSC_DECIDE, 10);
300: .ve
301: indicates that the global size of the vector is 10 and the local size will be automatically determined so that the sum of the
302: local sizes is the global size, see `PetscSplitOwnership()`.
304: Fortran Note:
305: Use `PETSC_DECIDE_INTEGER` or `PETSC_DECIDE_REAL`.
307: .seealso: `PETSC_DEFAULT`, `PETSC_IGNORE`, `PETSC_DETERMINE`, `PETSC_UNLIMITED`
308: M*/
310: /*MC
311: PETSC_DETERMINE - standard way of passing an integer or floating point parameter to indicate PETSc should determine an appropriate value
313: Level: beginner
315: Example Usage:
316: .vb
317: VecSetSizes(ksp, 10, PETSC_DETERMINE);
318: .ve
319: indicates that the local size of the vector is 10 and the global size will be automatically summing up all the local sizes.
321: Note:
322: Same as `PETSC_DECIDE`
324: Fortran Note:
325: Use `PETSC_DETERMINE_INTEGER` or `PETSC_DETERMINE_REAL`.
327: Developer Note:
328: I would like to use const `PetscInt` `PETSC_DETERMINE` = `PETSC_DECIDE`; but for
329: some reason this is not allowed by the standard even though `PETSC_DECIDE` is a constant value.
331: .seealso: `PETSC_DECIDE`, `PETSC_DEFAULT`, `PETSC_IGNORE`, `VecSetSizes()`, `PETSC_UNLIMITED`
332: M*/
334: /*MC
335: PETSC_CURRENT - standard way of indicating to an object not to change the current value of the parameter in the object
337: Level: beginner
339: Note:
340: Use `PETSC_DECIDE` to use the value that was set by PETSc when the object's type was set
342: Fortran Note:
343: Use `PETSC_CURRENT_INTEGER` or `PETSC_CURRENT_REAL`.
345: .seealso: `PETSC_DECIDE`, `PETSC_IGNORE`, `PETSC_DETERMINE`, `PETSC_DEFAULT`, `PETSC_UNLIMITED`
346: M*/
348: /*MC
349: PETSC_DEFAULT - deprecated, see `PETSC_CURRENT` and `PETSC_DETERMINE`
351: Level: beginner
353: Note:
354: The name is confusing since it tells the object to continue to use the value it is using, not the default value when the object's type was set.
356: Developer Note:
357: Unfortunately this was used for two different purposes in the past, to actually trigger the use of a default value or to continue the
358: use of currently set value (in, for example, `KSPSetTolerances()`.
360: .seealso: `PETSC_DECIDE`, `PETSC_IGNORE`, `PETSC_DETERMINE`, `PETSC_CURRENT`, `PETSC_UNLIMITED`
361: M*/
363: /* These MUST be preprocessor defines! see https://gitlab.com/petsc/petsc/-/issues/1370 */
364: #define PETSC_DECIDE (-1)
365: #define PETSC_DETERMINE PETSC_DECIDE
366: #define PETSC_CURRENT (-2)
367: #define PETSC_UNLIMITED (-3)
368: /* PETSC_DEFAULT is deprecated in favor of PETSC_CURRENT for use in KSPSetTolerances() and similar functions */
369: #define PETSC_DEFAULT PETSC_CURRENT
371: /*MC
372: PETSC_COMM_WORLD - the equivalent of the `MPI_COMM_WORLD` communicator which represents all the processes that PETSc knows about.
374: Level: beginner
376: Notes:
377: By default `PETSC_COMM_WORLD` and `MPI_COMM_WORLD` are identical unless you wish to
378: run PETSc on ONLY a subset of `MPI_COMM_WORLD`. In that case create your new (smaller)
379: communicator, call it, say comm, and set `PETSC_COMM_WORLD` = comm BEFORE calling
380: `PetscInitialize()`, but after `MPI_Init()` has been called.
382: The value of `PETSC_COMM_WORLD` should never be used or accessed before `PetscInitialize()`
383: is called because it may not have a valid value yet.
385: .seealso: `PETSC_COMM_SELF`
386: M*/
387: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;
389: /*MC
390: PETSC_COMM_SELF - This is always `MPI_COMM_SELF`
392: Level: beginner
394: Note:
395: Do not USE/access or set this variable before `PetscInitialize()` has been called.
397: .seealso: `PETSC_COMM_WORLD`
398: M*/
399: #define PETSC_COMM_SELF MPI_COMM_SELF
401: /*MC
402: PETSC_MPI_THREAD_REQUIRED - the required threading support used if PETSc initializes MPI with `MPI_Init_thread()`.
404: No Fortran Support
406: Level: beginner
408: Note:
409: By default `PETSC_MPI_THREAD_REQUIRED` equals `MPI_THREAD_FUNNELED` when the MPI implementation provides `MPI_Init_thread()`, otherwise it equals `MPI_THREAD_SINGLE`
411: .seealso: `PetscInitialize()`
412: M*/
413: PETSC_EXTERN PetscMPIInt PETSC_MPI_THREAD_REQUIRED;
415: /*MC
416: PetscBeganMPI - indicates if PETSc initialized MPI using `MPI_Init()` during `PetscInitialize()` or if MPI was already initialized with `MPI_Init()`
418: Synopsis:
419: #include <petscsys.h>
420: PetscBool PetscBeganMPI;
422: No Fortran Support
424: Level: developer
426: Note:
427: `MPI_Init()` can never be called after `PetscInitialize()`
429: .seealso: `PetscInitialize()`, `PetscInitializeCalled`
430: M*/
431: PETSC_EXTERN PetscBool PetscBeganMPI;
433: PETSC_EXTERN PetscBool PetscErrorHandlingInitialized;
434: PETSC_EXTERN PetscBool PetscInitializeCalled;
435: PETSC_EXTERN PetscBool PetscFinalizeCalled;
436: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;
438: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm), PetscErrorCode (*)(MPI_Comm));
439: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm, MPI_Comm *, int *);
440: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm *);
441: PETSC_EXTERN PetscErrorCode PetscCommGetComm(MPI_Comm, MPI_Comm *);
442: PETSC_EXTERN PetscErrorCode PetscCommRestoreComm(MPI_Comm, MPI_Comm *);
444: #if defined(PETSC_HAVE_KOKKOS)
445: PETSC_EXTERN PetscErrorCode PetscKokkosInitializeCheck(void); /* Initialize Kokkos if not yet. */
446: #endif
448: #if defined(PETSC_HAVE_NVSHMEM)
449: PETSC_EXTERN PetscBool PetscBeganNvshmem;
450: PETSC_EXTERN PetscBool PetscNvshmemInitialized;
451: PETSC_EXTERN PetscErrorCode PetscNvshmemFinalize(void);
452: #endif
454: #if defined(PETSC_HAVE_ELEMENTAL)
455: PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
456: PETSC_EXTERN PetscErrorCode PetscElementalInitialized(PetscBool *);
457: PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
458: #endif
460: /*MC
461: PetscMalloc - Allocates memory for use with PETSc. One should use `PetscNew()`, `PetscMalloc1()` or `PetscCalloc1()` usually instead of `PetscMalloc()`
463: Synopsis:
464: #include <petscsys.h>
465: PetscErrorCode PetscMalloc(size_t m,void **result)
467: Not Collective
469: Input Parameter:
470: . m - number of bytes to allocate
472: Output Parameter:
473: . result - memory allocated
475: Level: beginner
477: Notes:
478: Memory is always allocated at least double aligned
480: It is safe to allocate with an m of 0 and pass the resulting pointer to `PetscFree()`.
481: However, the pointer should never be dereferenced or the program will crash.
483: Developer Note:
484: All the `PetscMallocN()` routines actually call `PetscMalloc()` behind the scenes.
486: Except for data structures that store information about the PETSc options database all memory allocated by PETSc is
487: obtained with `PetscMalloc()` or `PetscCalloc()`
489: .seealso: `PetscFree()`, `PetscNew()`, `PetscCalloc()`
490: M*/
491: #define PetscMalloc(a, b) ((*PetscTrMalloc)((a), PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (void **)(b)))
493: /*MC
494: PetscRealloc - Reallocates memory
496: Synopsis:
497: #include <petscsys.h>
498: PetscErrorCode PetscRealloc(size_t m,void **result)
500: Not Collective
502: Input Parameters:
503: + m - number of bytes to allocate
504: - result - previous memory
506: Output Parameter:
507: . result - new memory allocated
509: Level: developer
511: Notes:
512: `results` must have already been obtained with `PetscMalloc()`
514: Memory is always allocated at least double aligned
516: .seealso: `PetscMalloc()`, `PetscFree()`, `PetscNew()`
517: M*/
518: #define PetscRealloc(a, b) ((*PetscTrRealloc)((a), __LINE__, PETSC_FUNCTION_NAME, __FILE__, (void **)(b)))
520: /*MC
521: PetscAddrAlign - Rounds up an address to `PETSC_MEMALIGN` alignment
523: Synopsis:
524: #include <petscsys.h>
525: void *PetscAddrAlign(void *addr)
527: Not Collective
529: Input Parameter:
530: . addr - address to align (any pointer type)
532: Level: developer
534: .seealso: `PetscMallocAlign()`
535: M*/
536: #define PetscAddrAlign(a) ((void *)((((PETSC_UINTPTR_T)(a)) + (PETSC_MEMALIGN - 1)) & ~(PETSC_MEMALIGN - 1)))
538: /*MC
539: PetscCalloc - Allocates a cleared (zeroed) memory region aligned to `PETSC_MEMALIGN`, similar to `PetscMalloc()`
541: Synopsis:
542: #include <petscsys.h>
543: PetscErrorCode PetscCalloc(size_t m,void **result)
545: Not Collective
547: Input Parameter:
548: . m - number of bytes to allocate
550: Output Parameter:
551: . result - memory allocated
553: Level: beginner
555: Notes:
556: Memory is always allocated at least double aligned. This macro is useful in allocating memory pointed by void pointers
558: It is safe to allocate with an m of 0 and pass the resulting pointer to `PetscFree()`.
560: However, the pointer should never be dereferenced or the program will crash.
562: Developer Note:
563: All `PetscCallocN()` routines call `PetscCalloc()` behind the scenes.
565: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`
566: M*/
567: #define PetscCalloc(m, result) PetscMallocA(1, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)m), (result))
569: /*MC
570: PetscMalloc1 - Allocates an array of memory aligned to `PETSC_MEMALIGN`
572: Synopsis:
573: #include <petscsys.h>
574: PetscErrorCode PetscMalloc1(size_t m1,type **r1)
576: Not Collective
578: Input Parameter:
579: . m1 - number of elements to allocate (may be zero)
581: Output Parameter:
582: . r1 - memory allocated
584: Level: beginner
586: Note:
587: This uses `sizeof()` of the memory type requested to determine the total memory to be allocated; therefore, you should not
588: multiply the number of elements requested by the `sizeof()` the type. For example, use
589: .vb
590: PetscInt *id;
591: PetscMalloc1(10,&id);
592: .ve
593: not
594: .vb
595: PetscInt *id;
596: PetscMalloc1(10*sizeof(PetscInt),&id);
597: .ve
599: Does not zero the memory allocated, use `PetscCalloc1()` to obtain memory that has been zeroed.
601: The `PetscMalloc[N]()` and `PetscCalloc[N]()` take an argument of type `size_t`! However, most codes use `value`, computed via `int` or `PetscInt` variables. This can overflow in
602: 32bit `int` computation - while computation in 64bit `size_t` would not overflow!
603: It's best if any arithmetic that is done for size computations is done with `size_t` type - avoiding arithmetic overflow!
605: `PetscMalloc[N]()` and `PetscCalloc[N]()` attempt to work-around this by casting the first variable to `size_t`.
606: This works for most expressions, but not all, such as
607: .vb
608: PetscInt *id, a, b;
609: PetscMalloc1(use_a_squared ? a * a * b : a * b, &id); // use_a_squared is cast to size_t, but a and b are still PetscInt
610: PetscMalloc1(a + b * b, &id); // a is cast to size_t, but b * b is performed at PetscInt precision first due to order-of-operations
611: .ve
613: These expressions should either be avoided, or appropriately cast variables to `size_t`:
614: .vb
615: PetscInt *id, a, b;
616: PetscMalloc1(use_a_squared ? (size_t)a * a * b : (size_t)a * b, &id); // Cast a to size_t before multiplication
617: PetscMalloc1(b * b + a, &id); // b is automatically cast to size_t and order-of-operations ensures size_t precision is maintained
618: .ve
620: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc2()`
621: M*/
622: #define PetscMalloc1(m1, r1) PetscMallocA(1, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1))
624: /*MC
625: PetscCalloc1 - Allocates a cleared (zeroed) array of memory aligned to `PETSC_MEMALIGN`
627: Synopsis:
628: #include <petscsys.h>
629: PetscErrorCode PetscCalloc1(size_t m1,type **r1)
631: Not Collective
633: Input Parameter:
634: . m1 - number of elements to allocate in 1st chunk (may be zero)
636: Output Parameter:
637: . r1 - memory allocated
639: Level: beginner
641: Note:
642: See `PetscMalloc1()` for more details on usage.
644: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc2()`
645: M*/
646: #define PetscCalloc1(m1, r1) PetscMallocA(1, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1))
648: /*MC
649: PetscMalloc2 - Allocates 2 arrays of memory both aligned to `PETSC_MEMALIGN`
651: Synopsis:
652: #include <petscsys.h>
653: PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)
655: Not Collective
657: Input Parameters:
658: + m1 - number of elements to allocate in 1st chunk (may be zero)
659: - m2 - number of elements to allocate in 2nd chunk (may be zero)
661: Output Parameters:
662: + r1 - memory allocated in first chunk
663: - r2 - memory allocated in second chunk
665: Level: developer
667: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc2()`
668: M*/
669: #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))
671: /*MC
672: PetscCalloc2 - Allocates 2 cleared (zeroed) arrays of memory both aligned to `PETSC_MEMALIGN`
674: Synopsis:
675: #include <petscsys.h>
676: PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)
678: Not Collective
680: Input Parameters:
681: + m1 - number of elements to allocate in 1st chunk (may be zero)
682: - m2 - number of elements to allocate in 2nd chunk (may be zero)
684: Output Parameters:
685: + r1 - memory allocated in first chunk
686: - r2 - memory allocated in second chunk
688: Level: developer
690: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc2()`
691: M*/
692: #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))
694: /*MC
695: PetscMalloc3 - Allocates 3 arrays of memory, all aligned to `PETSC_MEMALIGN`
697: Synopsis:
698: #include <petscsys.h>
699: PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
701: Not Collective
703: Input Parameters:
704: + m1 - number of elements to allocate in 1st chunk (may be zero)
705: . m2 - number of elements to allocate in 2nd chunk (may be zero)
706: - m3 - number of elements to allocate in 3rd chunk (may be zero)
708: Output Parameters:
709: + r1 - memory allocated in first chunk
710: . r2 - memory allocated in second chunk
711: - r3 - memory allocated in third chunk
713: Level: developer
715: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc3()`, `PetscFree3()`
716: M*/
717: #define PetscMalloc3(m1, r1, m2, r2, m3, r3) \
718: 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))
720: /*MC
721: PetscCalloc3 - Allocates 3 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
723: Synopsis:
724: #include <petscsys.h>
725: PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
727: Not Collective
729: Input Parameters:
730: + m1 - number of elements to allocate in 1st chunk (may be zero)
731: . m2 - number of elements to allocate in 2nd chunk (may be zero)
732: - m3 - number of elements to allocate in 3rd chunk (may be zero)
734: Output Parameters:
735: + r1 - memory allocated in first chunk
736: . r2 - memory allocated in second chunk
737: - r3 - memory allocated in third chunk
739: Level: developer
741: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc2()`, `PetscMalloc3()`, `PetscFree3()`
742: M*/
743: #define PetscCalloc3(m1, r1, m2, r2, m3, r3) \
744: 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))
746: /*MC
747: PetscMalloc4 - Allocates 4 arrays of memory, all aligned to `PETSC_MEMALIGN`
749: Synopsis:
750: #include <petscsys.h>
751: PetscErrorCode PetscMalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
753: Not Collective
755: Input Parameters:
756: + m1 - number of elements to allocate in 1st chunk (may be zero)
757: . m2 - number of elements to allocate in 2nd chunk (may be zero)
758: . m3 - number of elements to allocate in 3rd chunk (may be zero)
759: - m4 - number of elements to allocate in 4th chunk (may be zero)
761: Output Parameters:
762: + r1 - memory allocated in first chunk
763: . r2 - memory allocated in second chunk
764: . r3 - memory allocated in third chunk
765: - r4 - memory allocated in fourth chunk
767: Level: developer
769: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc4()`, `PetscFree4()`
770: M*/
771: #define PetscMalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \
772: 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))
774: /*MC
775: PetscCalloc4 - Allocates 4 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
777: Synopsis:
778: #include <petscsys.h>
779: PetscErrorCode PetscCalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
781: Not Collective
783: Input Parameters:
784: + m1 - number of elements to allocate in 1st chunk (may be zero)
785: . m2 - number of elements to allocate in 2nd chunk (may be zero)
786: . m3 - number of elements to allocate in 3rd chunk (may be zero)
787: - m4 - number of elements to allocate in 4th chunk (may be zero)
789: Output Parameters:
790: + r1 - memory allocated in first chunk
791: . r2 - memory allocated in second chunk
792: . r3 - memory allocated in third chunk
793: - r4 - memory allocated in fourth chunk
795: Level: developer
797: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc4()`, `PetscFree4()`
798: M*/
799: #define PetscCalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \
800: 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))
802: /*MC
803: PetscMalloc5 - Allocates 5 arrays of memory, all aligned to `PETSC_MEMALIGN`
805: Synopsis:
806: #include <petscsys.h>
807: 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)
809: Not Collective
811: Input Parameters:
812: + m1 - number of elements to allocate in 1st chunk (may be zero)
813: . m2 - number of elements to allocate in 2nd chunk (may be zero)
814: . m3 - number of elements to allocate in 3rd chunk (may be zero)
815: . m4 - number of elements to allocate in 4th chunk (may be zero)
816: - m5 - number of elements to allocate in 5th chunk (may be zero)
818: Output Parameters:
819: + r1 - memory allocated in first chunk
820: . r2 - memory allocated in second chunk
821: . r3 - memory allocated in third chunk
822: . r4 - memory allocated in fourth chunk
823: - r5 - memory allocated in fifth chunk
825: Level: developer
827: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc5()`, `PetscFree5()`
828: M*/
829: #define PetscMalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \
830: 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))
832: /*MC
833: PetscCalloc5 - Allocates 5 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
835: Synopsis:
836: #include <petscsys.h>
837: 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)
839: Not Collective
841: Input Parameters:
842: + m1 - number of elements to allocate in 1st chunk (may be zero)
843: . m2 - number of elements to allocate in 2nd chunk (may be zero)
844: . m3 - number of elements to allocate in 3rd chunk (may be zero)
845: . m4 - number of elements to allocate in 4th chunk (may be zero)
846: - m5 - number of elements to allocate in 5th chunk (may be zero)
848: Output Parameters:
849: + r1 - memory allocated in first chunk
850: . r2 - memory allocated in second chunk
851: . r3 - memory allocated in third chunk
852: . r4 - memory allocated in fourth chunk
853: - r5 - memory allocated in fifth chunk
855: Level: developer
857: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc5()`, `PetscFree5()`
858: M*/
859: #define PetscCalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \
860: 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))
862: /*MC
863: PetscMalloc6 - Allocates 6 arrays of memory, all aligned to `PETSC_MEMALIGN`
865: Synopsis:
866: #include <petscsys.h>
867: 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)
869: Not Collective
871: Input Parameters:
872: + m1 - number of elements to allocate in 1st chunk (may be zero)
873: . m2 - number of elements to allocate in 2nd chunk (may be zero)
874: . m3 - number of elements to allocate in 3rd chunk (may be zero)
875: . m4 - number of elements to allocate in 4th chunk (may be zero)
876: . m5 - number of elements to allocate in 5th chunk (may be zero)
877: - m6 - number of elements to allocate in 6th chunk (may be zero)
879: Output Parameteasr:
880: + r1 - memory allocated in first chunk
881: . r2 - memory allocated in second chunk
882: . r3 - memory allocated in third chunk
883: . r4 - memory allocated in fourth chunk
884: . r5 - memory allocated in fifth chunk
885: - r6 - memory allocated in sixth chunk
887: Level: developer
889: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc6()`, `PetscFree3()`, `PetscFree4()`, `PetscFree5()`, `PetscFree6()`
890: M*/
891: #define PetscMalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \
892: 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))
894: /*MC
895: PetscCalloc6 - Allocates 6 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
897: Synopsis:
898: #include <petscsys.h>
899: 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)
901: Not Collective
903: Input Parameters:
904: + m1 - number of elements to allocate in 1st chunk (may be zero)
905: . m2 - number of elements to allocate in 2nd chunk (may be zero)
906: . m3 - number of elements to allocate in 3rd chunk (may be zero)
907: . m4 - number of elements to allocate in 4th chunk (may be zero)
908: . m5 - number of elements to allocate in 5th chunk (may be zero)
909: - m6 - number of elements to allocate in 6th chunk (may be zero)
911: Output Parameters:
912: + r1 - memory allocated in first chunk
913: . r2 - memory allocated in second chunk
914: . r3 - memory allocated in third chunk
915: . r4 - memory allocated in fourth chunk
916: . r5 - memory allocated in fifth chunk
917: - r6 - memory allocated in sixth chunk
919: Level: developer
921: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscMalloc6()`, `PetscFree6()`
922: M*/
923: #define PetscCalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \
924: 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))
926: /*MC
927: PetscMalloc7 - Allocates 7 arrays of memory, all aligned to `PETSC_MEMALIGN`
929: Synopsis:
930: #include <petscsys.h>
931: 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)
933: Not Collective
935: Input Parameters:
936: + m1 - number of elements to allocate in 1st chunk (may be zero)
937: . m2 - number of elements to allocate in 2nd chunk (may be zero)
938: . m3 - number of elements to allocate in 3rd chunk (may be zero)
939: . m4 - number of elements to allocate in 4th chunk (may be zero)
940: . m5 - number of elements to allocate in 5th chunk (may be zero)
941: . m6 - number of elements to allocate in 6th chunk (may be zero)
942: - m7 - number of elements to allocate in 7th chunk (may be zero)
944: Output Parameters:
945: + r1 - memory allocated in first chunk
946: . r2 - memory allocated in second chunk
947: . r3 - memory allocated in third chunk
948: . r4 - memory allocated in fourth chunk
949: . r5 - memory allocated in fifth chunk
950: . r6 - memory allocated in sixth chunk
951: - r7 - memory allocated in seventh chunk
953: Level: developer
955: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc7()`, `PetscFree7()`
956: M*/
957: #define PetscMalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \
958: 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))
960: /*MC
961: PetscCalloc7 - Allocates 7 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN`
963: Synopsis:
964: #include <petscsys.h>
965: 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)
967: Not Collective
969: Input Parameters:
970: + m1 - number of elements to allocate in 1st chunk (may be zero)
971: . m2 - number of elements to allocate in 2nd chunk (may be zero)
972: . m3 - number of elements to allocate in 3rd chunk (may be zero)
973: . m4 - number of elements to allocate in 4th chunk (may be zero)
974: . m5 - number of elements to allocate in 5th chunk (may be zero)
975: . m6 - number of elements to allocate in 6th chunk (may be zero)
976: - m7 - number of elements to allocate in 7th chunk (may be zero)
978: Output Parameters:
979: + r1 - memory allocated in first chunk
980: . r2 - memory allocated in second chunk
981: . r3 - memory allocated in third chunk
982: . r4 - memory allocated in fourth chunk
983: . r5 - memory allocated in fifth chunk
984: . r6 - memory allocated in sixth chunk
985: - r7 - memory allocated in seventh chunk
987: Level: developer
989: .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscMalloc7()`, `PetscFree7()`
990: M*/
991: #define PetscCalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \
992: 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))
994: /*MC
995: PetscNew - Allocates memory of a particular type, zeros the memory! Aligned to `PETSC_MEMALIGN`
997: Synopsis:
998: #include <petscsys.h>
999: PetscErrorCode PetscNew(type **result)
1001: Not Collective
1003: Output Parameter:
1004: . result - memory allocated, sized to match pointer `type`
1006: Level: beginner
1008: Developer Note:
1009: Calls `PetscCalloc()` with the appropriate memory size obtained from `type`
1011: .seealso: `PetscFree()`, `PetscMalloc()`, `PetscCall()`, `PetscCalloc1()`, `PetscMalloc1()`
1012: M*/
1013: #define PetscNew(b) PetscCalloc1(1, (b))
1015: #define PetscNewLog(o, b) PETSC_DEPRECATED_MACRO(3, 18, 0, "PetscNew()", ) PetscNew(b)
1017: /*MC
1018: PetscFree - Frees memory
1020: Synopsis:
1021: #include <petscsys.h>
1022: PetscErrorCode PetscFree(void *memory)
1024: Not Collective
1026: Input Parameter:
1027: . memory - memory to free (the pointer is ALWAYS set to `NULL` upon success)
1029: Level: beginner
1031: Notes:
1032: Do not free memory obtained with `PetscMalloc2()`, `PetscCalloc2()` etc, they must be freed with `PetscFree2()` etc.
1034: It is safe to call `PetscFree()` on a `NULL` pointer.
1036: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc1()`
1037: M*/
1038: #define PetscFree(a) ((PetscErrorCode)((*PetscTrFree)((void *)(a), __LINE__, PETSC_FUNCTION_NAME, __FILE__) || ((a) = PETSC_NULLPTR, PETSC_SUCCESS)))
1040: /*MC
1041: PetscFree2 - Frees 2 chunks of memory obtained with `PetscMalloc2()`
1043: Synopsis:
1044: #include <petscsys.h>
1045: PetscErrorCode PetscFree2(void *memory1,void *memory2)
1047: Not Collective
1049: Input Parameters:
1050: + memory1 - memory to free
1051: - memory2 - 2nd memory to free
1053: Level: developer
1055: Notes:
1056: Memory must have been obtained with `PetscMalloc2()`
1058: The arguments need to be in the same order as they were in the call to `PetscMalloc2()`
1060: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`
1061: M*/
1062: #define PetscFree2(m1, m2) PetscFreeA(2, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2))
1064: /*MC
1065: PetscFree3 - Frees 3 chunks of memory obtained with `PetscMalloc3()`
1067: Synopsis:
1068: #include <petscsys.h>
1069: PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)
1071: Not Collective
1073: Input Parameters:
1074: + memory1 - memory to free
1075: . memory2 - 2nd memory to free
1076: - memory3 - 3rd memory to free
1078: Level: developer
1080: Notes:
1081: Memory must have been obtained with `PetscMalloc3()`
1083: The arguments need to be in the same order as they were in the call to `PetscMalloc3()`
1085: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`
1086: M*/
1087: #define PetscFree3(m1, m2, m3) PetscFreeA(3, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3))
1089: /*MC
1090: PetscFree4 - Frees 4 chunks of memory obtained with `PetscMalloc4()`
1092: Synopsis:
1093: #include <petscsys.h>
1094: PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)
1096: Not Collective
1098: Input Parameters:
1099: + m1 - memory to free
1100: . m2 - 2nd memory to free
1101: . m3 - 3rd memory to free
1102: - m4 - 4th memory to free
1104: Level: developer
1106: Notes:
1107: Memory must have been obtained with `PetscMalloc4()`
1109: The arguments need to be in the same order as they were in the call to `PetscMalloc4()`
1111: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`
1112: M*/
1113: #define PetscFree4(m1, m2, m3, m4) PetscFreeA(4, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4))
1115: /*MC
1116: PetscFree5 - Frees 5 chunks of memory obtained with `PetscMalloc5()`
1118: Synopsis:
1119: #include <petscsys.h>
1120: PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)
1122: Not Collective
1124: Input Parameters:
1125: + m1 - memory to free
1126: . m2 - 2nd memory to free
1127: . m3 - 3rd memory to free
1128: . m4 - 4th memory to free
1129: - m5 - 5th memory to free
1131: Level: developer
1133: Notes:
1134: Memory must have been obtained with `PetscMalloc5()`
1136: The arguments need to be in the same order as they were in the call to `PetscMalloc5()`
1138: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`
1139: M*/
1140: #define PetscFree5(m1, m2, m3, m4, m5) PetscFreeA(5, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5))
1142: /*MC
1143: PetscFree6 - Frees 6 chunks of memory obtained with `PetscMalloc6()`
1145: Synopsis:
1146: #include <petscsys.h>
1147: PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)
1149: Not Collective
1151: Input Parameters:
1152: + m1 - memory to free
1153: . m2 - 2nd memory to free
1154: . m3 - 3rd memory to free
1155: . m4 - 4th memory to free
1156: . m5 - 5th memory to free
1157: - m6 - 6th memory to free
1159: Level: developer
1161: Notes:
1162: Memory must have been obtained with `PetscMalloc6()`
1164: The arguments need to be in the same order as they were in the call to `PetscMalloc6()`
1166: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()`
1167: M*/
1168: #define PetscFree6(m1, m2, m3, m4, m5, m6) PetscFreeA(6, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6))
1170: /*MC
1171: PetscFree7 - Frees 7 chunks of memory obtained with `PetscMalloc7()`
1173: Synopsis:
1174: #include <petscsys.h>
1175: PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)
1177: Not Collective
1179: Input Parameters:
1180: + m1 - memory to free
1181: . m2 - 2nd memory to free
1182: . m3 - 3rd memory to free
1183: . m4 - 4th memory to free
1184: . m5 - 5th memory to free
1185: . m6 - 6th memory to free
1186: - m7 - 7th memory to free
1188: Level: developer
1190: Notes:
1191: Memory must have been obtained with `PetscMalloc7()`
1193: The arguments need to be in the same order as they were in the call to `PetscMalloc7()`
1195: .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()`,
1196: `PetscMalloc7()`
1197: M*/
1198: #define PetscFree7(m1, m2, m3, m4, m5, m6, m7) PetscFreeA(7, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6), &(m7))
1200: PETSC_EXTERN PetscErrorCode PetscMallocA(int, PetscBool, int, const char *, const char *, size_t, void *, ...);
1201: PETSC_EXTERN PetscErrorCode PetscFreeA(int, int, const char *, const char *, void *, ...);
1202: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t, PetscBool, int, const char[], const char[], void **);
1203: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void *, int, const char[], const char[]);
1204: PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t, int, const char[], const char[], void **);
1205: PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1206: 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 **));
1207: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);
1209: /*
1210: Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1211: */
1212: PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1213: PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);
1214: #if defined(PETSC_HAVE_CUDA)
1215: PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void);
1216: PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void);
1217: #endif
1218: #if defined(PETSC_HAVE_HIP)
1219: PETSC_EXTERN PetscErrorCode PetscMallocSetHIPHost(void);
1220: PETSC_EXTERN PetscErrorCode PetscMallocResetHIPHost(void);
1221: #endif
1223: #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE
1224: #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION
1226: /*
1227: Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1228: */
1229: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1230: PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1231: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1232: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1233: PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1234: PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int, PetscLogDouble *);
1235: PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool, PetscBool);
1236: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool *, PetscBool *, PetscBool *);
1237: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int, const char[], const char[]);
1238: PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1239: PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool *);
1240: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool);
1241: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool *);
1243: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType, MPI_Datatype *);
1244: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype, PetscDataType *);
1245: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType, size_t *);
1246: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char *, PetscDataType *, PetscBool *);
1248: /*
1249: These are MPI operations for MPI_Allreduce() etc
1250: */
1251: PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1252: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1253: PETSC_EXTERN MPI_Op MPIU_SUM;
1254: PETSC_EXTERN MPI_Op MPIU_MAX;
1255: PETSC_EXTERN MPI_Op MPIU_MIN;
1256: #else
1257: #define MPIU_SUM MPI_SUM
1258: #define MPIU_MAX MPI_MAX
1259: #define MPIU_MIN MPI_MIN
1260: #endif
1261: PETSC_EXTERN MPI_Op Petsc_Garbage_SetIntersectOp;
1262: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm, const PetscInt[], PetscInt *, PetscInt *);
1264: #if (defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128)) || (defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16))
1265: /*MC
1266: MPIU_SUM___FP16___FLOAT128 - MPI_Op that acts as a replacement for `MPI_SUM` with
1267: custom `MPI_Datatype` `MPIU___FLOAT128`, `MPIU___COMPLEX128`, and `MPIU___FP16`.
1269: Level: advanced
1271: Developer Note:
1272: This should be unified with `MPIU_SUM`
1274: .seealso: `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`
1275: M*/
1276: PETSC_EXTERN MPI_Op MPIU_SUM___FP16___FLOAT128;
1277: #endif
1279: /*
1280: These are so that in extern C code we can cast function pointers to non-extern C
1281: function pointers. Since the regular C++ code expects its function pointers to be C++
1282: */
1284: /*S
1285: PetscVoidFn - A prototype of a `void fn(void)` function
1287: Level: advanced
1289: Notes:
1290: `PetscVoidFn *` plays the role of `void *` for function pointers in the PETSc API that do not return an error code.
1291: It is used where a function pointer is needed but it is not possible to use the full prototype of the function.
1293: `PetscErrorCodeFn` is similar to `PetscVoidFn` but should be used when the function returns a `PetscErrorCode`
1295: The deprecated `PetscVoidFunction` works as a replacement for `PetscVoidFn` *.
1297: The deprecated `PetscVoidStarFunction` works as a replacement for `PetscVoidFn` **.
1299: .seealso: `PetscErrorCodeFn`, `PetscObject`, `PetscObjectDestroy()`
1300: S*/
1301: PETSC_EXTERN_TYPEDEF typedef void PetscVoidFn(void);
1303: PETSC_EXTERN_TYPEDEF typedef PetscVoidFn *PetscVoidFunction;
1304: PETSC_EXTERN_TYPEDEF typedef PetscVoidFn **PetscVoidStarFunction;
1306: /*S
1307: PetscErrorCodeFn - a function typedef that represents abstractly a function that returns a PETSc error code
1308: and takes any number of arguments. Since C/C++ has no way to express this concept, it is implemented as `void (fn)(void)`.
1310: Level: advanced
1312: Notes:
1313: `PetscErrorCodeFn *` plays the role of `void *` for function pointers in the PETSc API that return an error code.
1314: It is used where a function pointer is needed but it is not possible to use the full prototype of the function,
1315: for example `VecSetOperation()`.
1317: `PetscVoidFn` is similar to `PetscErrorCodeFn` but should be used when the function does not return a `PetscErrorCode`.
1319: Developer Notes:
1320: This function type is equivalent to `PetscVoidFn`*.
1322: At the C/C++ syntax level this construct adds nothing of value to the PETSc source code. It provides a way, at the abstract
1323: PETSc API level, to indicate specifically functions that return PETSc error codes as opposed to any C/C++ function.
1325: .seealso: `PetscVoidFn`, `PetscObject`, `PetscObjectDestroy()`, `VecSetOperation()`
1326: S*/
1327: PETSC_EXTERN_TYPEDEF typedef void PetscErrorCodeFn(void);
1329: /*
1330: Defines PETSc error handling.
1331: */
1332: #include <petscerror.h>
1334: PETSC_EXTERN PetscBool PetscCIEnabled; /* code is running in the PETSc test harness CI */
1335: PETSC_EXTERN PetscBool PetscCIEnabledPortableErrorOutput; /* error output is stripped to ensure portability of error messages across systems */
1336: PETSC_EXTERN const char *PetscCIFilename(const char *);
1337: PETSC_EXTERN int PetscCILinenumber(int);
1339: #define PETSC_SMALLEST_CLASSID 1211211
1340: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1341: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1342: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[], PetscClassId *);
1343: PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject, PetscObjectId *);
1344: PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject, PetscObjectId, PetscBool *);
1346: /*
1347: Routines that get memory usage information from the OS
1348: */
1349: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1350: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1351: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1352: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);
1354: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);
1356: /*
1357: Initialization of PETSc
1358: */
1359: PETSC_EXTERN PetscErrorCode PetscInitialize(int *, char ***, const char[], const char[]);
1360: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int, char *[], const char[], const char[]);
1361: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1362: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1363: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1364: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1365: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1366: PETSC_EXTERN PetscErrorCode PetscGetArgs(int *, char ***);
1367: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1368: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);
1370: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1371: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);
1372: PETSC_EXTERN PetscErrorCode PetscSysFinalizePackage(void);
1374: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[], const char[]);
1375: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1376: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1377: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject, const char[]);
1379: /*
1380: Functions that can act on any PETSc object.
1381: */
1382: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject *);
1383: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject, MPI_Comm *);
1384: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject, PetscClassId *);
1385: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject, const char *[]);
1386: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject, const char *[]);
1387: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject, const char[]);
1388: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject, const char *[]);
1389: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject, PetscInt);
1390: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject, PetscInt *);
1391: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject, PetscObject, PetscInt);
1392: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1393: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject, PetscInt *);
1394: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1395: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject, PetscMPIInt *);
1396: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject, const char[], PetscObject);
1397: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject, const char[]);
1398: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject, const char[], PetscObject *);
1399: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject, const char[], PetscErrorCodeFn *);
1400: #define PetscObjectComposeFunction(a, b, ...) PetscObjectComposeFunction_Private((a), (b), (PetscErrorCodeFn *)(__VA_ARGS__))
1401: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1402: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1403: PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1404: PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject, PetscObject);
1405: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm, PetscMPIInt *);
1407: /*MC
1408: PetscObjectParameterSetDefault - sets a parameter default value in a `PetscObject` to a new default value.
1409: If the current value matches the old default value, then the current value is also set to the new value.
1411: No Fortran Support
1413: Synopsis:
1414: #include <petscsys.h>
1415: PetscBool PetscObjectParameterSetDefault(PetscObject obj, char* NAME, PetscReal value);
1417: Input Parameters:
1418: + obj - the `PetscObject`
1419: . NAME - the name of the parameter, unquoted
1420: - value - the new value
1422: Level: developer
1424: Notes:
1425: The defaults for an object are the values set when the object's type is set.
1427: This should only be used in object constructors, such as, `SNESCreate_NGS()`.
1429: This only works for parameters that are declared in the struct with `PetscObjectParameterDeclare()`
1431: .seealso: `PetscObjectParameterDeclare()`, `PetscInitialize()`, `PetscFinalize()`, `PetscObject`, `SNESParametersInitialize()`
1432: M*/
1433: #define PetscObjectParameterSetDefault(obj, NAME, value) \
1434: do { \
1435: if (obj->NAME == obj->default_##NAME) obj->NAME = value; \
1436: obj->default_##NAME = value; \
1437: } while (0)
1439: /*MC
1440: PetscObjectParameterDeclare - declares a parameter in a `PetscObject` and a location to store its default
1442: No Fortran Support
1444: Synopsis:
1445: #include <petscsys.h>
1446: PetscBool PetscObjectParameterDeclare(type, char* NAME)
1448: Input Parameters:
1449: + type - the type of the parameter, for example `PetscInt`
1450: - NAME - the name of the parameter, unquoted
1452: Level: developer.
1454: .seealso: `PetscObjectParameterSetDefault()`, `PetscInitialize()`, `PetscFinalize()`, `PetscObject`, `SNESParametersInitialize()`
1455: M*/
1456: #define PetscObjectParameterDeclare(type, NAME) type NAME, default_##NAME
1457: #define PetscObjectParameterDeclarePtr(type, NAME) type *NAME, *default_##NAME
1459: /*MC
1460: PetscCtx - indicates an argument that can be a pointer to any C struct (or Fortran derived type).
1462: Level: developer
1464: Notes:
1465: This should not be used for arrays of unknown type.
1467: Fortran Notes:
1468: A Fortran code that calls a function with a `PetscCtx` argument would declare the variable `ctx` with
1469: .vb
1470: type(AppType) :: ctx
1471: .ve
1472: where `AppType` is a Fortran derived type. Or the argument can be a `PetscObject`.
1474: Developer Note:
1475: `PetscCtx` is used instead of `void *` in PETSc code to enhance the clarity of the PETSc source code since `void *` serves so many different roles.
1476: The getAPI() code processor also uses the variable type to generate correct bindings for other languages.
1478: .seealso: [](sec_fortran_context), `PetscCtxRt`, PetscCtxDestroyFn()`, `PeOp`, `PeNS`, `PetscInitialize()`, `DMGetApplicationContext()`,
1479: `DMSetApplicationContextDestroy()`
1480: M*/
1481: typedef void *PetscCtx;
1483: /*MC
1484: PetscCtxRt - indicates an argument that returns a pointer to a C struct (or Fortran derived type) which is generally an application context
1486: Level: developer
1488: Notes:
1489: A PETSc object (in C or Fortran) can be used as a PETSc context
1491: This should not be used for functions that return pointers to arrays of unknown type. Thus it is used for, for example,
1492: `KSPGetApplicationContext()` but not used for `DMNetworkGetComponent()`
1494: A PETSc object (in C or Fortran) can be used as a PETSc context
1496: It is also used for functions that destroy an application context. For example, the destroy function passed to `DMSetApplicationContextDestroy()`
1497: which has a prototype of `PetscCtxDestroyFn()`
1499: This typedef is not part of the PETSc public API and should only be used in PETSc source code.
1501: For pointers to arrays of unknown type and for functions that return PETSc internal objects that are opaque to users, such
1502: as `KSPMonitorDynamicToleranceCreate()` a `void **` should be used.
1504: Fortran Notes:
1505: A Fortran code that calls a function with a `PetscCtxRt` argument must declare the variable `ctx` with
1506: .vb
1507: type(AppType), pointer :: ctx
1508: .ve
1509: where `AppType` is a Fortran derived type.
1511: If one passes a PETSc function with a `PetscCtxRt` argument as an argument in Fortran one must use the function named suffixed with `Cptr`,
1512: for example `KSPConvergedDefaultDestroyCptr`, see src/ksp/ksp/tutorials/ex1f.F90.
1514: Developer Notes:
1515: C++ compilers generate a warning or error if one passes a pointer to a pointer to a specific type (instead of `void`), for example,
1516: .vb
1517: extern calledfunction(void **);
1518: SomeCtx *ctx;
1519: calledfunction(&ctx); << warning that it is passing a pointer to a pointer to a SomeCtx instead of a void **
1520: .ve
1521: By using the common practice of prototyping the function as
1522: .vb
1523: extern calledfunction(void *);
1524: .ve
1525: the warning message is averted.
1527: `PetscCtxRt` is used instead of `void *` in PETSc code to enhance the clarity of the PETSc source code since `void *` serves so many different roles.
1528: The getAPI() code processor also uses the variable type to generate correct bindings for other languages.
1530: The Fortran C stub and Fortran interface definition generated for functions with a `PetscCtxRt` argument are the C function name suffixed with
1531: `Cptr`, for example `KSPConvergedDefaultDestroyCptr`. The Fortran user API is a macro with the original C funtion name, for example,
1532: `KSPConvergedDefaultDestroy` that calls the `KSPConvergedDefaultDestroyCptr` version and then calls `c_f_pointer()` to handle the equivalent of a `void**` cast
1533: to the users Fortran derived type argument.
1535: .seealso: [](sec_fortran_context), `PetscCtx`, `PetscCtxDestroyFn()`, `PeOp`, `PeNS`, `PetscInitialize()`, `DMGetApplicationContext()`,
1536: `DMSetApplicationContextDestroy()`
1537: M*/
1538: typedef void *PetscCtxRt;
1540: /*S
1541: PetscCtxDestroyFn - A prototype of a `PetscErrorCode (*)(PetscCtxRt)` function that is used to free application contexts
1543: Level: intermediate
1545: Notes:
1546: Used in the prototype of functions such as `DMSetApplicationContextDestroy()`
1548: The function argument is a `PetscCtxRt` which is psychologically equivalent to a `void **` meaning that this function is called with a pointer to
1549: the application context (which is itself a pointer) thus the destroy implementation must first reference the context via, for example,
1550: `*(AppCtx **)arg`. Note that syntactically `PetscCtxRt` is defined as a `void *`, this is because C++ does
1551: not accept passing a pointer to a pointer to a `void**` but it does accept passing a pointer to a pointer to `void *`.
1553: PETSc destroy functions take the address of the context (rather than just the context) so that that the destroy function can "zero the pointer" when
1554: appropriate, preventing accidental later use of a dangling pointer.
1556: .seealso: `PetscObject`, `PetscCtxDestroyDefault()`, `PetscObjectDestroy()`, `DMSetApplicationContextDestroy()`
1557: S*/
1558: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PetscCtxDestroyFn(PetscCtxRt ctx);
1560: PETSC_EXTERN PetscCtxDestroyFn PetscCtxDestroyDefault;
1561: PETSC_DEPRECATED_FUNCTION(3, 23, 0, "PetscCtxDestroyDefault()", ) static inline PetscErrorCode PetscContainerCtxDestroyDefault(PetscCtxRt a)
1562: {
1563: return PetscCtxDestroyDefault(a);
1564: }
1566: PETSC_EXTERN PetscErrorCode PetscMonitorCompare(PetscErrorCode (*)(void), void *, PetscCtxDestroyFn *, PetscErrorCode (*)(void), void *, PetscCtxDestroyFn *, PetscBool *);
1568: #include <petscviewertypes.h>
1569: #include <petscoptions.h>
1571: PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer, PetscBool, PetscLogDouble);
1572: PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool *);
1574: PETSC_EXTERN PetscErrorCode PetscObjectsListGetGlobalNumbering(MPI_Comm, PetscInt, PetscObject[], PetscInt *, PetscInt *);
1576: PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer, const char[]);
1577: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject, PetscViewer);
1578: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject, PetscViewer);
1579: #define PetscObjectQueryFunction(obj, name, fptr) PetscObjectQueryFunction_Private((obj), (name), (PetscErrorCodeFn **)(fptr))
1580: PETSC_EXTERN PetscErrorCode PetscObjectHasFunction(PetscObject, const char[], PetscBool *);
1581: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject, const char[], PetscErrorCodeFn **);
1582: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject, const char[]);
1583: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject, const char[]);
1584: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject, const char[]);
1585: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject, const char *[]);
1586: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject, const char[]);
1587: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1588: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1589: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject, PetscObject, const char[]);
1590: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1591: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject, const char[], PetscBool *);
1592: PETSC_EXTERN PetscErrorCode PetscObjectObjectTypeCompare(PetscObject, PetscObject, PetscBool *);
1593: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject, const char[], PetscBool *);
1594: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject, PetscBool *, const char[], ...);
1595: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject, PetscBool *, const char[], ...);
1596: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1597: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);
1599: #if defined(PETSC_HAVE_SAWS)
1600: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1601: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1602: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject, PetscBool);
1603: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1604: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1605: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1606: PETSC_EXTERN void PetscStackSAWsGrantAccess(void);
1607: PETSC_EXTERN void PetscStackSAWsTakeAccess(void);
1608: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1609: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);
1611: #else
1612: #define PetscSAWsBlock() PETSC_SUCCESS
1613: #define PetscObjectSAWsViewOff(obj) PETSC_SUCCESS
1614: #define PetscObjectSAWsSetBlock(obj, flg) PETSC_SUCCESS
1615: #define PetscObjectSAWsBlock(obj) PETSC_SUCCESS
1616: #define PetscObjectSAWsGrantAccess(obj) PETSC_SUCCESS
1617: #define PetscObjectSAWsTakeAccess(obj) PETSC_SUCCESS
1618: #define PetscStackViewSAWs() PETSC_SUCCESS
1619: #define PetscStackSAWsViewOff() PETSC_SUCCESS
1620: #define PetscStackSAWsTakeAccess()
1621: #define PetscStackSAWsGrantAccess()
1623: #endif
1625: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[], PetscDLMode, PetscDLHandle *);
1626: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1627: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle, const char[], void **);
1628: PETSC_EXTERN PetscErrorCode PetscDLAddr(PetscVoidFn *, char *[]);
1629: PETSC_EXTERN PetscErrorCode PetscDemangleSymbol(const char *, char *[]);
1631: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void *, PetscStack **);
1633: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE *, PetscBool);
1634: PETSC_EXTERN PetscErrorCode PetscObjectsView(PetscViewer);
1635: PETSC_EXTERN PetscErrorCode PetscObjectsGetObject(const char *, PetscObject *, const char *[]);
1636: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList *);
1637: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList, const char[], PetscObject *);
1638: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList, PetscObject, const char *[], PetscBool *);
1639: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *, const char[], PetscObject);
1640: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *, const char[]);
1641: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList, PetscObjectList *);
1643: /*
1644: Dynamic library lists. Lists of names of routines in objects or in dynamic
1645: link libraries that will be loaded as needed.
1646: */
1648: #define PetscFunctionListAdd(list, name, fptr) PetscFunctionListAdd_Private((list), (name), (PetscErrorCodeFn *)(fptr))
1649: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList *, const char[], PetscErrorCodeFn *);
1650: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList *);
1651: PETSC_EXTERN PetscErrorCode PetscFunctionListClear(PetscFunctionList);
1652: #define PetscFunctionListFind(list, name, fptr) PetscFunctionListFind_Private((list), (name), (PetscErrorCodeFn **)(fptr))
1653: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList, const char[], PetscErrorCodeFn **);
1654: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm, FILE *, const char[], const char[], const char[], const char[], PetscFunctionList, const char[], const char[]);
1655: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList, PetscFunctionList *);
1656: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList, PetscViewer);
1657: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList, const char ***, int *);
1658: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintNonEmpty(PetscFunctionList);
1659: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintAll(void);
1661: PETSC_EXTERN PetscDLLibrary PetscDLLibrariesLoaded;
1662: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm, PetscDLLibrary *, const char[]);
1663: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm, PetscDLLibrary *, const char[]);
1664: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm, PetscDLLibrary *, const char[], const char[], void **);
1665: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1666: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm, const char[], char *, size_t, PetscBool *);
1667: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm, const char[], PetscDLLibrary *);
1668: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);
1670: /*
1671: Useful utility routines
1672: */
1673: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm, PetscInt *, PetscInt *);
1674: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm, PetscInt, PetscInt *, PetscInt *);
1675: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm, PetscInt *, PetscInt *);
1676: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm, PetscMPIInt);
1677: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm, PetscMPIInt);
1678: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1679: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE *);
1680: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm, const PetscInt[2], PetscInt[2]);
1681: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm, const PetscReal[2], PetscReal[2]);
1683: /*MC
1684: PetscNot - negates a logical type value and returns result as a `PetscBool`
1686: Level: beginner
1688: Note:
1689: This is useful in cases like
1690: .vb
1691: int *a;
1692: PetscBool flag = PetscNot(a)
1693: .ve
1694: where !a would not return a `PetscBool` because we cannot provide a cast from int to `PetscBool` in C.
1696: .seealso: `PetscBool`, `PETSC_TRUE`, `PETSC_FALSE`
1697: M*/
1698: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)
1700: /*MC
1701: PetscHelpPrintf - Prints help messages.
1703: Synopsis:
1704: #include <petscsys.h>
1705: PetscErrorCode (*PetscHelpPrintf)(MPI_Comm comm, const char format[],args);
1707: Not Collective, only applies on MPI rank 0; No Fortran Support
1709: Input Parameters:
1710: + comm - the MPI communicator over which the help message is printed
1711: . format - the usual printf() format string
1712: - args - arguments to be printed
1714: Level: developer
1716: Notes:
1717: You can change how help messages are printed by replacing the function pointer with a function that does not simply write to stdout.
1719: To use, write your own function, for example,
1720: .vb
1721: PetscErrorCode mypetschelpprintf(MPI_Comm comm,const char format[],....)
1722: {
1723: PetscFunctionReturn(PETSC_SUCCESS);
1724: }
1725: .ve
1726: then do the assignment
1727: .vb
1728: PetscHelpPrintf = mypetschelpprintf;
1729: .ve
1731: You can do the assignment before `PetscInitialize()`.
1733: The default routine used is called `PetscHelpPrintfDefault()`.
1735: .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscErrorPrintf()`, `PetscHelpPrintfDefault()`
1736: M*/
1737: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1739: /*
1740: Defines PETSc profiling.
1741: */
1742: #include <petsclog.h>
1744: /*
1745: Simple PETSc parallel IO for ASCII printing
1746: */
1747: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[], char[]);
1748: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm, const char[], const char[], FILE **);
1749: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm, FILE *);
1750: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1751: PETSC_EXTERN PetscErrorCode PetscFFlush(FILE *);
1752: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1753: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char *, size_t, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1754: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char *, size_t, const char[], size_t *, ...) PETSC_ATTRIBUTE_FORMAT(3, 5);
1755: PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[], size_t, const char *, PetscInt, const PetscReal[]);
1757: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1758: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2);
1759: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1761: PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char *, size_t *);
1762: PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char *, char *);
1764: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm, const char[], const char[], const char[], FILE **);
1765: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm, FILE *);
1766: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1768: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3);
1769: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4);
1770: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm, FILE *);
1771: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm, FILE *, size_t, char[]);
1772: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm, const char[], const char[], FILE **);
1773: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char *[]);
1775: PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1776: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer, void *);
1777: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer, void *);
1778: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer *);
1779: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm, PetscContainer *);
1780: PETSC_EXTERN PetscErrorCode PetscContainerSetCtxDestroy(PetscContainer, PetscCtxDestroyFn *);
1781: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 23, 0, "PetscContainerSetCtxDestroy()", ) PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void *));
1782: PETSC_EXTERN PetscErrorCode PetscObjectContainerCompose(PetscObject, const char *name, void *, PetscCtxDestroyFn *);
1783: PETSC_EXTERN PetscErrorCode PetscObjectContainerQuery(PetscObject, const char *, PetscCtxRt);
1785: /*
1786: For use in debuggers
1787: */
1788: PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1789: PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1790: PETSC_EXTERN PetscErrorCode PetscIntViewNumColumns(PetscInt, PetscInt, const PetscInt[], PetscViewer);
1791: PETSC_EXTERN PetscErrorCode PetscRealViewNumColumns(PetscInt, PetscInt, const PetscReal[], PetscViewer);
1792: PETSC_EXTERN PetscErrorCode PetscScalarViewNumColumns(PetscInt, PetscInt, const PetscScalar[], PetscViewer);
1793: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt, const PetscInt[], PetscViewer);
1794: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt, const PetscReal[], PetscViewer);
1795: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt, const PetscScalar[], PetscViewer);
1797: /*
1798: Basic memory and string operations. These are usually simple wrappers
1799: around the basic Unix system calls, but a few of them have additional
1800: functionality and/or error checking.
1801: */
1802: #include <petscstring.h>
1804: #include <stddef.h>
1805: #include <stdlib.h>
1807: #if defined(PETSC_CLANG_STATIC_ANALYZER)
1808: #define PetscPrefetchBlock(a, b, c, d)
1809: #else
1810: /*MC
1811: PetscPrefetchBlock - Prefetches a block of memory
1813: Synopsis:
1814: #include <petscsys.h>
1815: void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)
1817: Not Collective
1819: Input Parameters:
1820: + a - pointer to first element to fetch (any type but usually `PetscInt` or `PetscScalar`)
1821: . n - number of elements to fetch
1822: . rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors)
1823: - t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note
1825: Level: developer
1827: Notes:
1828: The last two arguments (`rw` and `t`) must be compile-time constants.
1830: Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality. Not all architectures offer
1831: equivalent locality hints, but the following macros are always defined to their closest analogue.
1832: + `PETSC_PREFETCH_HINT_NTA` - Non-temporal. Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched).
1833: . `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.
1834: . `PETSC_PREFETCH_HINT_T1` - Fetch to level 2 and higher (not L1).
1835: - `PETSC_PREFETCH_HINT_T2` - Fetch to high-level cache only. (On many systems, T0 and T1 are equivalent.)
1837: This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid
1838: address).
1840: M*/
1841: #define PetscPrefetchBlock(a, n, rw, t) \
1842: do { \
1843: const char *_p = (const char *)(a), *_end = (const char *)((a) + (n)); \
1844: for (; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p, (rw), (t)); \
1845: } while (0)
1846: #endif
1847: /*
1848: Determine if some of the kernel computation routines use
1849: Fortran (rather than C) for the numerical calculations. On some machines
1850: and compilers (like complex numbers) the Fortran version of the routines
1851: is faster than the C/C++ versions. The flag --with-fortran-kernels
1852: should be used with ./configure to turn these on.
1853: */
1854: #if defined(PETSC_USE_FORTRAN_KERNELS)
1856: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1857: #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1858: #endif
1860: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1861: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
1862: #endif
1864: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1865: #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
1866: #endif
1868: #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
1869: #define PETSC_USE_FORTRAN_KERNEL_MAXPY
1870: #endif
1872: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1873: #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
1874: #endif
1876: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
1877: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
1878: #endif
1880: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1881: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
1882: #endif
1884: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
1885: #define PETSC_USE_FORTRAN_KERNEL_MDOT
1886: #endif
1888: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
1889: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
1890: #endif
1892: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
1893: #define PETSC_USE_FORTRAN_KERNEL_AYPX
1894: #endif
1896: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
1897: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
1898: #endif
1900: #endif
1902: /*
1903: Macros for indicating code that should be compiled with a C interface,
1904: rather than a C++ interface. Any routines that are dynamically loaded
1905: (such as the PCCreate_XXX() routines) must be wrapped so that the name
1906: mangler does not change the functions symbol name. This just hides the
1907: ugly extern "C" {} wrappers.
1908: */
1909: #if defined(__cplusplus)
1910: #define EXTERN_C_BEGIN extern "C" {
1911: #define EXTERN_C_END }
1912: #else
1913: #define EXTERN_C_BEGIN
1914: #define EXTERN_C_END
1915: #endif
1917: /*MC
1918: MPI_Comm - the basic object used by MPI to determine which processes are involved in a
1919: communication
1921: Level: beginner
1923: Note:
1924: This manual page is a place-holder because MPICH does not have a manual page for `MPI_Comm`
1926: .seealso: `PETSC_COMM_WORLD`, `PETSC_COMM_SELF`
1927: M*/
1929: #if defined(PETSC_HAVE_MPIIO)
1930: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4);
1931: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4);
1932: 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);
1933: 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);
1934: 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);
1935: 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);
1936: #endif
1938: #if defined(PETSC_HAVE_MPI_COUNT)
1939: typedef MPI_Count MPIU_Count;
1940: #else
1941: typedef PetscInt64 MPIU_Count;
1942: #endif
1944: /*@C
1945: PetscIntCast - casts a `MPI_Count`, `PetscInt64`, `PetscCount`, or `size_t` to a `PetscInt` (which may be 32-bits in size), generates an
1946: error if the `PetscInt` is not large enough to hold the number.
1948: Not Collective; No Fortran Support
1950: Input Parameter:
1951: . a - the `PetscInt64` value
1953: Output Parameter:
1954: . b - the resulting `PetscInt` value, or `NULL` if the result is not needed
1956: Level: advanced
1958: Note:
1959: 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
1961: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscCIntCast()`, `PetscIntMultError()`, `PetscIntSumError()`
1962: @*/
1963: static inline PetscErrorCode PetscIntCast(MPIU_Count a, PetscInt *b)
1964: {
1965: PetscFunctionBegin;
1966: if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */
1967: PetscCheck(sizeof(MPIU_Count) <= sizeof(PetscInt) || (a <= (MPIU_Count)PETSC_INT_MAX && a >= (MPIU_Count)PETSC_INT_MIN), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt64_FMT " is too big for PetscInt, you may need to ./configure using --with-64-bit-indices", (PetscInt64)a);
1968: if (b) *b = (PetscInt)a;
1969: PetscFunctionReturn(PETSC_SUCCESS);
1970: }
1972: /*@C
1973: PetscBLASIntCast - casts a `MPI_Count`, `PetscInt`, `PetscCount` or `PetscInt64` to a `PetscBLASInt` (which may be 32-bits in size), generates an
1974: error if the `PetscBLASInt` is not large enough to hold the number.
1976: Not Collective; No Fortran Support
1978: Input Parameter:
1979: . a - the `PetscInt` value
1981: Output Parameter:
1982: . b - the resulting `PetscBLASInt` value, or `NULL` if the result is not needed
1984: Level: advanced
1986: Note:
1987: Errors if the integer is negative since PETSc calls to BLAS/LAPACK never need to cast negative integer inputs
1989: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscCIntCast()`, `PetscIntCast()`
1990: @*/
1991: static inline PetscErrorCode PetscBLASIntCast(MPIU_Count a, PetscBLASInt *b)
1992: {
1993: PetscFunctionBegin;
1994: if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */
1995: PetscCheck(sizeof(MPIU_Count) <= sizeof(PetscBLASInt) || a <= (MPIU_Count)PETSC_BLAS_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt64_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", (PetscInt64)a);
1996: PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer to BLAS/LAPACK routine");
1997: if (b) *b = (PetscBLASInt)a;
1998: PetscFunctionReturn(PETSC_SUCCESS);
1999: }
2001: /*@C
2002: PetscCuBLASIntCast - like `PetscBLASIntCast()`, but for `PetscCuBLASInt`.
2004: Not Collective; No Fortran Support
2006: Input Parameter:
2007: . a - the `PetscInt` value
2009: Output Parameter:
2010: . b - the resulting `PetscCuBLASInt` value, or `NULL` if the result is not needed
2012: Level: advanced
2014: Note:
2015: Errors if the integer is negative since PETSc calls to cuBLAS and friends never need to cast negative integer inputs
2017: .seealso: `PetscCuBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscCIntCast()`, `PetscIntCast()`
2018: @*/
2019: static inline PetscErrorCode PetscCuBLASIntCast(MPIU_Count a, PetscCuBLASInt *b)
2020: {
2021: PetscFunctionBegin;
2022: if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */
2023: PetscCheck(sizeof(MPIU_Count) <= sizeof(PetscCuBLASInt) || a <= (MPIU_Count)PETSC_CUBLAS_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt64_FMT " is too big for cuBLAS, which is restricted to 32-bit integers.", (PetscInt64)a);
2024: PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer %" PetscInt64_FMT "to cuBLAS routine", (PetscInt64)a);
2025: if (b) *b = (PetscCuBLASInt)a;
2026: PetscFunctionReturn(PETSC_SUCCESS);
2027: }
2029: /*@C
2030: PetscHipBLASIntCast - like `PetscBLASIntCast()`, but for `PetscHipBLASInt`.
2032: Not Collective; No Fortran Support
2034: Input Parameter:
2035: . a - the `PetscInt` value
2037: Output Parameter:
2038: . b - the resulting `PetscHipBLASInt` value, or `NULL` if the result is not needed
2040: Level: advanced
2042: Note:
2043: Errors if the integer is negative since PETSc calls to hipBLAS and friends never need to cast negative integer inputs
2045: .seealso: `PetscHipBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscCIntCast()`, `PetscIntCast()`
2046: @*/
2047: static inline PetscErrorCode PetscHipBLASIntCast(MPIU_Count a, PetscHipBLASInt *b)
2048: {
2049: PetscFunctionBegin;
2050: if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */
2051: PetscCheck(sizeof(MPIU_Count) <= sizeof(PetscHipBLASInt) || a <= (MPIU_Count)PETSC_HIPBLAS_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt64_FMT " is too big for hipBLAS, which is restricted to 32-bit integers.", (PetscInt64)a);
2052: PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer %" PetscInt64_FMT "to hipBLAS routine", (PetscInt64)a);
2053: if (b) *b = (PetscHipBLASInt)a;
2054: PetscFunctionReturn(PETSC_SUCCESS);
2055: }
2057: /*@C
2058: PetscMPIIntCast - casts a `MPI_Count`, `PetscInt`, `PetscCount`, or `PetscInt64` to a `PetscMPIInt` (which is always 32-bits in size), generates an
2059: error if the `PetscMPIInt` is not large enough to hold the number.
2061: Not Collective; No Fortran Support
2063: Input Parameter:
2064: . a - the `PetscInt` value
2066: Output Parameter:
2067: . b - the resulting `PetscMPIInt` value, or `NULL` if the result is not needed
2069: Level: advanced
2071: .seealso: [](stylePetscCount), `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscIntCast()`
2072: @*/
2073: static inline PetscErrorCode PetscMPIIntCast(MPIU_Count a, PetscMPIInt *b)
2074: {
2075: PetscFunctionBegin;
2076: if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */
2077: PetscCheck(a <= (MPIU_Count)PETSC_MPI_INT_MAX && a >= (MPIU_Count)PETSC_MPI_INT_MIN, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt64_FMT " is too big for MPI buffer length. Maximum supported value is %d", (PetscInt64)a, PETSC_MPI_INT_MAX);
2078: if (b) *b = (PetscMPIInt)a;
2079: PetscFunctionReturn(PETSC_SUCCESS);
2080: }
2082: /*@C
2083: PetscCIntCast - casts a `MPI_Count`, `PetscInt`, `PetscCount`, or `PetscInt64` to a `int`, generates an error if the `int` is not large enough to hold the number.
2085: Not Collective; No Fortran Support
2087: Input Parameter:
2088: . a - the `PetscInt` value
2090: Output Parameter:
2091: . b - the resulting `int` value, or `NULL` if the result is not needed
2093: Level: advanced
2095: .seealso: [](stylePetscCount), `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscIntCast()`
2096: @*/
2097: static inline PetscErrorCode PetscCIntCast(MPIU_Count a, int *b)
2098: {
2099: PetscFunctionBegin;
2100: if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */
2101: PetscCheck(a <= INT_MAX && a >= INT_MIN, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt64_FMT " is too big to be casted to an int. Maximum supported value is %d", (PetscInt64)a, INT_MAX);
2102: if (b) *b = (int)a;
2103: PetscFunctionReturn(PETSC_SUCCESS);
2104: }
2106: /*MC
2107: PetscInt64Mult - Computes the product of two variables after casting them to `PetscInt64`.
2109: Not Collective; No Fortran Support
2111: Input Parameters:
2112: + a - the first variable
2113: - b - the second variable
2115: Level: advanced
2117: .seealso: [](stylePetscCount), `PetscIntMultError()`, `PetscIntMultTruncate()`
2118: M*/
2119: #if defined(PETSC_USE_64BIT_INDICES)
2120: #define PetscInt64Mult(a, b) ((a) * (b))
2121: #else
2122: #define PetscInt64Mult(a, b) (((PetscInt64)(a)) * ((PetscInt64)(b)))
2123: #endif
2125: /*@C
2126: PetscRealIntMultTruncate - Computes the product of a positive `PetscReal` and a positive
2127: `PetscInt` and truncates the value to slightly less than the maximal possible value.
2129: Not Collective; No Fortran Support
2131: Input Parameters:
2132: + a - The `PetscReal` value
2133: - b - The `PetscInt` value
2135: Level: advanced
2137: Notes:
2138: Returns the result as a `PetscInt` value.
2140: Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`.
2142: Use `PetscIntMultTruncate()` to compute the product of two positive `PetscInt` and truncate
2143: to fit a `PetscInt`.
2145: Use `PetscIntMultError()` to compute the product of two `PetscInt` if you wish to generate an
2146: error if the result will not fit in a `PetscInt`.
2148: Developer Notes:
2149: We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but
2150: requires many more checks.
2152: This is used where we compute approximate sizes for workspace and need to insure the
2153: workspace is index-able.
2155: .seealso: `PetscReal`, `PetscInt`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()`
2156: @*/
2157: static inline PetscInt PetscRealIntMultTruncate(PetscReal a, PetscInt b)
2158: {
2159: PetscInt64 r = (PetscInt64)(a * (PetscReal)b);
2160: if (r > PETSC_INT_MAX - 100) r = PETSC_INT_MAX - 100;
2161: #if defined(PETSC_USE_64BIT_INDICES)
2162: return r;
2163: #else
2164: return (PetscInt)r;
2165: #endif
2166: }
2168: /*@C
2169: PetscIntMultTruncate - Computes the product of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value
2171: Not Collective; No Fortran Support
2173: Input Parameters:
2174: + a - the `PetscInt` value
2175: - b - the second value
2177: Returns:
2178: The result as a `PetscInt` value
2180: Level: advanced
2182: Notes:
2183: Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`
2185: Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt`
2187: 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`
2189: Developer Notes:
2190: We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but requires many more checks.
2192: This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2194: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()`,
2195: `PetscIntSumTruncate()`
2196: @*/
2197: static inline PetscInt PetscIntMultTruncate(PetscInt a, PetscInt b)
2198: {
2199: PetscInt64 r = PetscInt64Mult(a, b);
2200: if (r > PETSC_INT_MAX - 100) r = PETSC_INT_MAX - 100;
2201: #if defined(PETSC_USE_64BIT_INDICES)
2202: return r;
2203: #else
2204: return (PetscInt)r;
2205: #endif
2206: }
2208: /*@C
2209: PetscIntSumTruncate - Computes the sum of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value
2211: Not Collective; No Fortran Support
2213: Input Parameters:
2214: + a - the `PetscInt` value
2215: - b - the second value
2217: Returns:
2218: The result as a `PetscInt` value
2220: Level: advanced
2222: Notes:
2223: Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`
2225: Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt`
2227: 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`
2229: Developer Note:
2230: This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2232: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`
2233: @*/
2234: static inline PetscInt PetscIntSumTruncate(PetscInt a, PetscInt b)
2235: {
2236: PetscInt64 r = a;
2238: r += b;
2239: if (r > PETSC_INT_MAX - 100) r = PETSC_INT_MAX - 100;
2240: #if defined(PETSC_USE_64BIT_INDICES)
2241: return r;
2242: #else
2243: return (PetscInt)r;
2244: #endif
2245: }
2247: /*@C
2248: PetscIntMultError - Computes the product of two positive `PetscInt` and generates an error with overflow.
2250: Not Collective; No Fortran Support
2252: Input Parameters:
2253: + a - the `PetscInt` value
2254: - b - the second value
2256: Output Parameter:
2257: . result - the result as a `PetscInt` value, or `NULL` if you do not want the result, you just want to check if it overflows
2259: Level: advanced
2261: Notes:
2262: Use `PetscInt64Mult()` to compute the product of two `PetscInt` and store in a `PetscInt64`
2264: Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt`
2266: Developer Note:
2267: 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.
2268: `PetscIntSumError()` can be used to check for this situation.
2270: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntSumError()`
2271: @*/
2272: static inline PetscErrorCode PetscIntMultError(PetscInt a, PetscInt b, PetscInt *result)
2273: {
2274: PetscInt64 r = PetscInt64Mult(a, b);
2276: PetscFunctionBegin;
2277: #if defined(PETSC_USE_64BIT_INDICES)
2278: if (result) *result = r;
2279: #else
2280: if (result) *result = (PetscInt)r;
2281: #endif
2282: if (!PetscDefined(USE_64BIT_INDICES)) {
2283: PetscCheck(r <= PETSC_INT_MAX, 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);
2284: }
2285: PetscFunctionReturn(PETSC_SUCCESS);
2286: }
2288: /*@C
2290: PetscIntSumError - Computes the sum of two positive `PetscInt` and generates an error with overflow.
2292: Not Collective; No Fortran Support
2294: Input Parameters:
2295: + a - the `PetscInt` value
2296: - b - the second value
2298: Output Parameter:
2299: . c - the result as a `PetscInt` value, or `NULL` if you do not want the result, you just want to check if it overflows
2301: Level: advanced
2303: Notes:
2304: Use `PetscInt64Mult()` to compute the product of two 32-bit `PetscInt` and store in a `PetscInt64`
2306: Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt`
2308: .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`
2309: @*/
2310: static inline PetscErrorCode PetscIntSumError(PetscInt a, PetscInt b, PetscInt *result)
2311: {
2312: PetscInt64 r = a;
2314: PetscFunctionBegin;
2315: r += b;
2316: #if defined(PETSC_USE_64BIT_INDICES)
2317: if (result) *result = r;
2318: #else
2319: if (result) *result = (PetscInt)r;
2320: #endif
2321: if (!PetscDefined(USE_64BIT_INDICES)) {
2322: PetscCheck(r <= PETSC_INT_MAX, 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);
2323: }
2324: PetscFunctionReturn(PETSC_SUCCESS);
2325: }
2327: /*
2328: The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2329: */
2330: #if defined(hz)
2331: #undef hz
2332: #endif
2334: #if defined(PETSC_HAVE_SYS_TYPES_H)
2335: #include <sys/types.h>
2336: #endif
2338: /*MC
2340: PETSC_VERSION - This manual page provides information about how PETSc documents and uses its version information. This information is available to both C/C++
2341: and Fortran compilers when `petscsys.h` is included.
2343: The current PETSc version and the API for accessing it are defined in <A HREF="../include/petscversion.h.html">include/petscversion.html</A>
2345: The complete version number is given as the triple PETSC_VERSION_MAJOR.PETSC_VERSION_MINOR.PETSC_VERSION_SUBMINOR (in short hand x.y.z)
2347: 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
2348: where only a change in the major version number (x) indicates a change in the API.
2350: 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
2352: 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),
2353: 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
2354: version number (x.y.z).
2356: `PETSC_RELEASE_DATE` is the date the x.y version was released (i.e. the version before any patch releases)
2358: `PETSC_VERSION_DATE` is the date the x.y.z version was released
2360: `PETSC_VERSION_GIT` is the last git commit to the repository given in the form vx.y.z-wwwww
2362: `PETSC_VERSION_DATE_GIT` is the date of the last git commit to the repository
2364: `PETSC_VERSION_()` is deprecated and will eventually be removed.
2366: Level: intermediate
2367: M*/
2369: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[], size_t);
2370: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[], size_t);
2371: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[], size_t);
2372: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[], size_t);
2373: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2374: PETSC_EXTERN PetscErrorCode PetscGetDate(char[], size_t);
2375: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2376: PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt *, PetscInt *, PetscInt *, PetscInt *);
2378: PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscCount, const PetscInt[], PetscBool *);
2379: PETSC_EXTERN PetscErrorCode PetscSortedInt64(PetscCount, const PetscInt64[], PetscBool *);
2380: PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscCount, const PetscMPIInt[], PetscBool *);
2381: PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscCount, const PetscReal[], PetscBool *);
2382: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscCount, PetscInt[]);
2383: PETSC_EXTERN PetscErrorCode PetscSortInt64(PetscCount, PetscInt64[]);
2384: PETSC_EXTERN PetscErrorCode PetscSortCount(PetscCount, PetscCount[]);
2385: PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscCount, PetscInt[]);
2386: PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *, PetscInt[]);
2387: PETSC_EXTERN PetscErrorCode PetscSortedCheckDupsInt(PetscCount, const PetscInt[], PetscBool *);
2388: PETSC_EXTERN PetscErrorCode PetscSortedCheckDupsCount(PetscCount, const PetscCount[], PetscBool *);
2389: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt *, PetscInt[]);
2390: PETSC_EXTERN PetscErrorCode PetscCheckDupsInt(PetscInt, const PetscInt[], PetscBool *);
2391: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscCount, const PetscInt[], PetscInt *);
2392: PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscCount, const PetscMPIInt[], PetscInt *);
2393: PETSC_EXTERN PetscErrorCode PetscFindCount(PetscCount, PetscCount, const PetscCount[], PetscCount *);
2394: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt, const PetscInt[], PetscInt[]);
2395: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt, const char *[], PetscInt[]);
2396: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscCount, PetscInt[], PetscInt[]);
2397: PETSC_EXTERN PetscErrorCode PetscSortIntWithCountArray(PetscCount, PetscInt[], PetscCount[]);
2398: PETSC_EXTERN PetscErrorCode PetscSortIntWithMPIIntArray(PetscCount, PetscInt[], PetscMPIInt[]);
2399: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscCount, PetscInt[], PetscInt[], PetscInt[]);
2400: PETSC_EXTERN PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount, PetscInt[], PetscInt[], PetscCount[]);
2401: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscCount, PetscMPIInt[]);
2402: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *, PetscMPIInt[]);
2403: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscCount, PetscMPIInt[], PetscMPIInt[]);
2404: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscCount, PetscMPIInt[], PetscInt[]);
2405: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscCount, PetscInt[], PetscScalar[]);
2406: PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscCount, PetscInt[], void *, size_t, void *);
2407: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscCount, PetscReal[]);
2408: PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscCount, PetscReal[], PetscInt[]);
2409: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt, const PetscReal[], PetscInt[]);
2410: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt *, PetscReal[]);
2411: PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal, PetscCount, const PetscReal[], PetscReal, PetscInt *);
2412: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt, PetscInt, PetscScalar[], PetscInt[]);
2413: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt, PetscInt, PetscReal[], PetscInt[]);
2414: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt, const PetscBool[], const PetscInt[], PetscInt *, PetscInt *[], PetscInt *[], PetscInt *[], PetscInt *[]);
2415: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt, const PetscInt[], const PetscInt[], PetscInt, const PetscInt[], const PetscInt[], PetscInt *, PetscInt *[], PetscInt *[]);
2416: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscInt *, PetscInt *[]);
2417: PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt, const PetscMPIInt[], PetscInt, const PetscMPIInt[], PetscInt *, PetscMPIInt *[]);
2418: PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *);
2420: PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt, void *, size_t, int (*)(const void *, const void *, void *), void *);
2421: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt, PetscInt[]);
2422: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt, PetscMPIInt[]);
2423: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt, PetscReal[]);
2424: PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt, void *, size_t, void *, size_t, int (*)(const void *, const void *, void *), void *);
2425: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt, PetscInt[], PetscInt[]);
2426: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt, PetscMPIInt[], PetscMPIInt[]);
2427: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt, PetscReal[], PetscInt[]);
2429: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2430: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[], size_t);
2432: /*J
2433: PetscRandomType - String with the name of a PETSc randomizer
2435: Level: beginner
2437: Note:
2438: To use `PETSCSPRNG` or `PETSCRANDOM123` you must have ./configure PETSc
2439: with the option `--download-sprng` or `--download-random123`. We recommend the default provided with PETSc.
2441: .seealso: `PetscRandomSetType()`, `PetscRandom`, `PetscRandomCreate()`
2442: J*/
2443: typedef const char *PetscRandomType;
2444: #define PETSCRAND "rand"
2445: #define PETSCRAND48 "rand48"
2446: #define PETSCSPRNG "sprng"
2447: #define PETSCRANDER48 "rander48"
2448: #define PETSCRANDOM123 "random123"
2449: #define PETSCCURAND "curand"
2451: /* Logging support */
2452: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;
2454: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);
2455: PETSC_EXTERN PetscErrorCode PetscRandomFinalizePackage(void);
2457: /* Dynamic creation and loading functions */
2458: PETSC_EXTERN PetscFunctionList PetscRandomList;
2460: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[], PetscErrorCode (*)(PetscRandom));
2461: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2462: PETSC_EXTERN PetscErrorCode PetscRandomSetOptionsPrefix(PetscRandom, const char[]);
2463: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2464: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType *);
2465: PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom, PetscObject, const char[]);
2466: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom, PetscViewer);
2468: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm, PetscRandom *);
2469: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom, PetscScalar *);
2470: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom, PetscReal *);
2471: PETSC_EXTERN PetscErrorCode PetscRandomGetValues(PetscRandom, PetscInt, PetscScalar *);
2472: PETSC_EXTERN PetscErrorCode PetscRandomGetValuesReal(PetscRandom, PetscInt, PetscReal *);
2473: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom, PetscScalar *, PetscScalar *);
2474: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom, PetscScalar, PetscScalar);
2475: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom, PetscInt64);
2476: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom, PetscInt64 *);
2477: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2478: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom *);
2480: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[], char[], size_t);
2481: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[], char[], size_t);
2482: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[], size_t);
2483: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[], char[]);
2484: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[], size_t);
2485: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[], char, PetscBool *);
2486: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[], char, PetscBool *);
2487: PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2488: PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2489: PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);
2491: /*MC
2492: PetscBinaryBigEndian - indicates if values in memory are stored with big endian format
2494: Synopsis:
2495: #include <petscsys.h>
2496: PetscBool PetscBinaryBigEndian(void);
2498: No Fortran Support
2500: Level: developer
2502: .seealso: `PetscInitialize()`, `PetscFinalize()`, `PetscInitializeCalled`
2503: M*/
2504: static inline PetscBool PetscBinaryBigEndian(void)
2505: {
2506: long _petsc_v = 1;
2507: return ((char *)&_petsc_v)[0] ? PETSC_FALSE : PETSC_TRUE;
2508: }
2510: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int, void *, PetscCount, PetscInt *, PetscDataType);
2511: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm, int, void *, PetscInt, PetscInt *, PetscDataType);
2512: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int, const void *, PetscCount, PetscDataType);
2513: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm, int, const void *, PetscInt, PetscDataType);
2514: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[], PetscFileMode, int *);
2515: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2516: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm, PetscBool *);
2517: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm, PetscBool *);
2518: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm, char[], size_t);
2519: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm, const char[], char[], size_t, PetscBool *);
2520: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm, const char[], char[], size_t, PetscBool *);
2521: #if defined(PETSC_USE_SOCKET_VIEWER)
2522: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[], int, int *);
2523: #endif
2525: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int, off_t, PetscBinarySeekType, off_t *);
2526: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm, int, off_t, PetscBinarySeekType, off_t *);
2527: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *, PetscDataType, PetscCount);
2529: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2530: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[], PetscBool);
2531: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2532: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char *);
2533: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2534: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2535: PETSC_EXTERN PetscErrorCode PetscWaitOnError(void);
2537: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt *);
2538: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], PetscMPIInt *[], PetscMPIInt *[]);
2539: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt *[], PetscMPIInt *[], PetscMPIInt *[]);
2540: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscInt ***, MPI_Request **);
2541: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscScalar ***, MPI_Request **);
2542: 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);
2543: 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 *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(6, 3);
2544: 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 *), PetscCtx ctx) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(6, 3);
2546: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm, PetscBuildTwoSidedType);
2547: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm, PetscBuildTwoSidedType *);
2549: PETSC_DEPRECATED_FUNCTION(3, 24, 0, "PetscSSEIsEnabled()", ) static inline PetscErrorCode PetscSSEIsEnabled(PETSC_UNUSED MPI_Comm comm, PetscBool *lflag, PetscBool *gflag)
2550: {
2551: if (lflag) *lflag = PETSC_FALSE;
2552: if (gflag) *gflag = PETSC_FALSE;
2553: return PETSC_SUCCESS;
2554: }
2556: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);
2558: struct _n_PetscSubcomm {
2559: MPI_Comm parent; /* parent communicator */
2560: MPI_Comm dupparent; /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2561: MPI_Comm child; /* the sub-communicator */
2562: PetscMPIInt n; /* num of subcommunicators under the parent communicator */
2563: PetscMPIInt color; /* color of processors belong to this communicator */
2564: PetscMPIInt *subsize; /* size of subcommunicator[color] */
2565: PetscSubcommType type;
2566: char *subcommprefix;
2567: };
2569: static inline MPI_Comm PetscSubcommParent(PetscSubcomm scomm)
2570: {
2571: return scomm->parent;
2572: }
2573: static inline MPI_Comm PetscSubcommChild(PetscSubcomm scomm)
2574: {
2575: return scomm->child;
2576: }
2577: static inline MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm)
2578: {
2579: return scomm->dupparent;
2580: }
2581: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm, PetscSubcomm *);
2582: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm *);
2583: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm, PetscInt);
2584: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm, PetscSubcommType);
2585: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm, PetscMPIInt, PetscMPIInt);
2586: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm, PetscViewer);
2587: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2588: PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm, const char[]);
2589: PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm, MPI_Comm *);
2590: PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm, MPI_Comm *);
2591: PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm, MPI_Comm *);
2593: PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt, PetscHeap *);
2594: PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap, PetscInt, PetscInt);
2595: PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap, PetscInt *, PetscInt *);
2596: PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap, PetscInt *, PetscInt *);
2597: PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap, PetscInt, PetscInt);
2598: PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2599: PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap *);
2600: PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap, PetscViewer);
2602: PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2603: PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm, PetscShmComm *);
2604: PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm, PetscMPIInt, PetscMPIInt *);
2605: PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm, PetscMPIInt, PetscMPIInt *);
2606: PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm, MPI_Comm *);
2608: /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2609: PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm, PetscInt, PetscOmpCtrl *);
2610: PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl, MPI_Comm *, MPI_Comm *, PetscBool *);
2611: PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl *);
2612: PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2613: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2614: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);
2616: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t, PetscCount, PetscSegBuffer *);
2617: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer *);
2618: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer, PetscCount, void *);
2619: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer, void *);
2620: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer, void *);
2621: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer, void *);
2622: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer, PetscCount *);
2623: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer, PetscCount);
2625: /*MC
2626: PetscSegBufferGetInts - access an array of `PetscInt` from a `PetscSegBuffer`
2628: Synopsis:
2629: #include <petscsys.h>
2630: PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg, size_t count, PetscInt *PETSC_RESTRICT *slot);
2632: No Fortran Support
2634: Input Parameters:
2635: + seg - `PetscSegBuffer` buffer
2636: - count - number of entries needed
2638: Output Parameter:
2639: . buf - address of new buffer for contiguous data
2641: Level: intermediate
2643: Developer Note:
2644: Type-safe wrapper to encourage use of PETSC_RESTRICT. Does not use PetscFunctionBegin because the error handling
2645: prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as
2646: possible.
2648: .seealso: `PetscSegBuffer`, `PetscSegBufferGet()`, `PetscInitialize()`, `PetscFinalize()`, `PetscInitializeCalled`
2649: M*/
2650: static inline PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg, PetscCount count, PetscInt *PETSC_RESTRICT *slot)
2651: {
2652: return PetscSegBufferGet(seg, count, (void **)slot);
2653: }
2655: extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2656: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted *);
2657: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted *);
2658: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted, const char *, const char *, PetscBool *);
2660: #include <stdarg.h>
2661: PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char *, size_t, const char[], size_t *, va_list);
2662: PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE *, const char[], va_list);
2664: PETSC_EXTERN PetscSegBuffer PetscCitationsList;
2666: /*@
2667: PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code.
2669: Not Collective; No Fortran Support
2671: Input Parameters:
2672: + cite - the bibtex item, formatted to displayed on multiple lines nicely
2673: - set - a boolean variable initially set to `PETSC_FALSE`; this is used to insure only a single registration of the citation
2675: Options Database Key:
2676: . -citations [filename] - print out the bibtex entries for the given computation
2678: Level: intermediate
2679: @*/
2680: static inline PetscErrorCode PetscCitationsRegister(const char cit[], PetscBool *set)
2681: {
2682: size_t len;
2683: char *vstring;
2685: PetscFunctionBegin;
2686: if (set && *set) PetscFunctionReturn(PETSC_SUCCESS);
2687: PetscCall(PetscStrlen(cit, &len));
2688: PetscCall(PetscSegBufferGet(PetscCitationsList, (PetscCount)len, &vstring));
2689: PetscCall(PetscArraycpy(vstring, cit, len));
2690: if (set) *set = PETSC_TRUE;
2691: PetscFunctionReturn(PETSC_SUCCESS);
2692: }
2694: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm, char[], char[], size_t);
2695: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm, const char[], char[], size_t);
2696: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm, const char[], const char[]);
2698: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm, char[], char[], size_t);
2699: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm, const char[], char[], char[], size_t);
2700: PETSC_EXTERN PetscErrorCode PetscBoxUpload(MPI_Comm, const char[], const char[]);
2702: PETSC_EXTERN PetscErrorCode PetscGlobusGetTransfers(MPI_Comm, const char[], char[], size_t);
2703: PETSC_EXTERN PetscErrorCode PetscGlobusAuthorize(MPI_Comm, char[], size_t);
2704: PETSC_EXTERN PetscErrorCode PetscGlobusUpload(MPI_Comm, const char[], const char[]);
2706: PETSC_EXTERN PetscErrorCode PetscPullJSONValue(const char[], const char[], char[], size_t, PetscBool *);
2707: PETSC_EXTERN PetscErrorCode PetscPushJSONValue(char[], const char[], const char[], size_t);
2709: #if !defined(PETSC_HAVE_MPI_LARGE_COUNT)
2710: /*
2711: Cast PetscCount <a> to PetscMPIInt <b>, where <a> is likely used for the 'count' argument in MPI routines.
2712: It is similar to PetscMPIIntCast() except that here it returns an MPI error code.
2713: */
2714: #define PetscMPIIntCast_Internal(a, b) \
2715: do { \
2716: *b = 0; \
2717: if (PetscUnlikely(a > (MPIU_Count)PETSC_MPI_INT_MAX)) return MPI_ERR_COUNT; \
2718: *b = (PetscMPIInt)a; \
2719: } while (0)
2721: static inline PetscMPIInt MPIU_Get_count(MPI_Status *status, MPI_Datatype dtype, PetscCount *count)
2722: {
2723: PetscMPIInt count2, err;
2725: *count = 0; /* to prevent incorrect warnings of uninitialized variables */
2726: err = MPI_Get_count(status, dtype, &count2);
2727: *count = count2;
2728: return err;
2729: }
2731: static inline PetscMPIInt MPIU_Send(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt dest, PetscMPIInt tag, MPI_Comm comm)
2732: {
2733: PetscMPIInt count2, err;
2735: PetscMPIIntCast_Internal(count, &count2);
2736: err = MPI_Send((void *)buf, count2, dtype, dest, tag, comm);
2737: return err;
2738: }
2740: static inline PetscMPIInt MPIU_Send_init(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt dest, PetscMPIInt tag, MPI_Comm comm, MPI_Request *request)
2741: {
2742: PetscMPIInt count2, err;
2744: PetscMPIIntCast_Internal(count, &count2);
2745: err = MPI_Send_init((void *)buf, count2, dtype, dest, tag, comm, request);
2746: return err;
2747: }
2749: static inline PetscMPIInt MPIU_Isend(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt dest, PetscMPIInt tag, MPI_Comm comm, MPI_Request *request)
2750: {
2751: PetscMPIInt count2, err;
2753: PetscMPIIntCast_Internal(count, &count2);
2754: err = MPI_Isend((void *)buf, count2, dtype, dest, tag, comm, request);
2755: return err;
2756: }
2758: static inline PetscMPIInt MPIU_Recv(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt source, PetscMPIInt tag, MPI_Comm comm, MPI_Status *status)
2759: {
2760: PetscMPIInt count2, err;
2762: PetscMPIIntCast_Internal(count, &count2);
2763: err = MPI_Recv((void *)buf, count2, dtype, source, tag, comm, status);
2764: return err;
2765: }
2767: static inline PetscMPIInt MPIU_Recv_init(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt source, PetscMPIInt tag, MPI_Comm comm, MPI_Request *request)
2768: {
2769: PetscMPIInt count2, err;
2771: PetscMPIIntCast_Internal(count, &count2);
2772: err = MPI_Recv_init((void *)buf, count2, dtype, source, tag, comm, request);
2773: return err;
2774: }
2776: static inline PetscMPIInt MPIU_Irecv(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt source, PetscMPIInt tag, MPI_Comm comm, MPI_Request *request)
2777: {
2778: PetscMPIInt count2, err;
2780: PetscMPIIntCast_Internal(count, &count2);
2781: err = MPI_Irecv((void *)buf, count2, dtype, source, tag, comm, request);
2782: return err;
2783: }
2785: static inline PetscMPIInt MPIU_Reduce(const void *inbuf, void *outbuf, MPIU_Count count, MPI_Datatype dtype, MPI_Op op, PetscMPIInt root, MPI_Comm comm)
2786: {
2787: PetscMPIInt count2, err;
2789: PetscMPIIntCast_Internal(count, &count2);
2790: err = MPI_Reduce((void *)inbuf, outbuf, count2, dtype, op, root, comm);
2791: return err;
2792: }
2794: #if defined(PETSC_HAVE_MPI_REDUCE_LOCAL)
2795: static inline PetscMPIInt MPIU_Reduce_local(const void *inbuf, void *inoutbuf, MPIU_Count count, MPI_Datatype dtype, MPI_Op op)
2796: {
2797: PetscMPIInt count2, err;
2799: PetscMPIIntCast_Internal(count, &count2);
2800: err = MPI_Reduce_local((void *)inbuf, inoutbuf, count2, dtype, op);
2801: return err;
2802: }
2803: #endif
2805: #if !defined(PETSC_USE_64BIT_INDICES)
2806: #define MPIU_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm) MPI_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm)
2807: #define MPIU_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm) MPI_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm)
2808: #else
2809: #define MPIU_Scatterv(sendbuf, sendcount, displs, sendtype, recvbuf, recvcount, recvtype, root, comm) \
2810: ((void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_SUP, PETSC_ERROR_INITIAL, "Must have MPI 4 support for MPI_Scatterv_c() for this functionality, upgrade your MPI"), MPI_ERR_COUNT)
2811: #define MPIU_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm) \
2812: ((void)PetscError(comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_SUP, PETSC_ERROR_INITIAL, "Must have MPI 4 support for MPI_Scatterv_c() for this functionality, upgrade your MPI"), MPI_ERR_COUNT)
2813: #endif
2815: #else
2817: /* on 32 bit systems MPI_Count maybe 64-bit while PetscCount is 32-bit */
2818: #define PetscCountCast_Internal(a, b) \
2819: do { \
2820: *b = 0; \
2821: if (PetscUnlikely(a > (MPI_Count)PETSC_COUNT_MAX)) return MPI_ERR_COUNT; \
2822: *b = (PetscMPIInt)a; \
2823: } while (0)
2825: static inline PetscMPIInt MPIU_Get_count(MPI_Status *status, MPI_Datatype dtype, PetscCount *count)
2826: {
2827: MPI_Count count2;
2828: PetscMPIInt err;
2830: *count = 0; /* to prevent incorrect warnings of uninitialized variables */
2831: err = MPI_Get_count_c(status, dtype, &count2);
2832: if (err) return err;
2833: PetscCountCast_Internal(count2, count);
2834: return MPI_SUCCESS;
2835: }
2837: #define MPIU_Reduce(inbuf, outbuf, count, dtype, op, root, comm) MPI_Reduce_c(inbuf, outbuf, (MPI_Count)(count), dtype, op, root, comm)
2838: #define MPIU_Send(buf, count, dtype, dest, tag, comm) MPI_Send_c(buf, (MPI_Count)(count), dtype, dest, tag, comm)
2839: #define MPIU_Send_init(buf, count, dtype, dest, tag, comm, request) MPI_Send_init_c(buf, (MPI_Count)(count), dtype, dest, tag, comm, request)
2840: #define MPIU_Isend(buf, count, dtype, dest, tag, comm, request) MPI_Isend_c(buf, (MPI_Count)(count), dtype, dest, tag, comm, request)
2841: #define MPIU_Recv(buf, count, dtype, source, tag, comm, status) MPI_Recv_c(buf, (MPI_Count)(count), dtype, source, tag, comm, status)
2842: #define MPIU_Recv_init(buf, count, dtype, source, tag, comm, request) MPI_Recv_init_c(buf, (MPI_Count)(count), dtype, source, tag, comm, request)
2843: #define MPIU_Irecv(buf, count, dtype, source, tag, comm, request) MPI_Irecv_c(buf, (MPI_Count)(count), dtype, source, tag, comm, request)
2844: #if defined(PETSC_HAVE_MPI_REDUCE_LOCAL)
2845: #define MPIU_Reduce_local(inbuf, inoutbuf, count, dtype, op) MPI_Reduce_local_c(inbuf, inoutbuf, (MPI_Count)(count), dtype, op)
2846: #endif
2848: /*MC
2849: MPIU_Scatterv - A replacement for `MPI_Scatterv()` that can be called with `PetscInt` types when PETSc is built for either 32-bit indices or 64-bit indices.
2851: Synopsis:
2852: #include <petscsys.h>
2853: PetscMPIInt MPIU_Scatterv(const void *sendbuf, const PetscInt sendcounts[], const PetscInt displs[], MPI_Datatype sendtype, void *recvbuf, PetscInt recvcount, MPI_Datatype recvtype, PetscMPIInt root, MPI_Comm comm)
2855: Collective
2857: Input Parameters:
2858: + sendbuf - address of send buffer
2859: . sendcounts - non-negative `PetscInt` array (of length `comm` group size) specifying the number of elements to send to each MPI process
2860: . displs - `PetscInt` array (of length `comm` group size). Entry i specifies the displacement (relative to `sendbuf`) from which to take the outgoing data to process i
2861: . sendtype - data type of `sendbuf` elements
2862: . recvcount - number of elements in `recvbuf` (non-negative integer)
2863: . recvtype - data type of `recvbuf` elements
2864: . root - Rank of the MPI root process, which will dispatch the data to scatter
2865: - comm - `MPI_Comm` communicator
2867: Output Parameter:
2868: . recvbuf - the resulting scattered values on this MPI process
2870: Level: developer
2872: Notes:
2873: Should be wrapped with `PetscCallMPI()` for error checking
2875: This is different than most of the `MPIU_` wrappers in that all the count arguments are in `PetscInt`
2877: .seealso: [](stylePetscCount), `MPI_Allreduce()`, `MPIU_Gatherv()`
2878: M*/
2880: #if !defined(PETSC_USE_64BIT_INDICES)
2881: #define MPIU_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm) MPI_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm)
2882: #define MPIU_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm) MPI_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm)
2883: #else
2884: #define MPIU_Scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm) MPI_Scatterv_c(sendbuf, (const MPI_Count *)(sendcounts), (const MPI_Aint *)(displs), sendtype, recvbuf, recvcount, recvtype, root, comm)
2885: #define MPIU_Gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm) MPI_Gatherv_c(sendbuf, sendcount, sendtype, recvbuf, (const MPI_Count *)(recvcounts), (const MPI_Aint *)(displs), recvtype, root, comm)
2886: #endif
2888: #endif
2890: PETSC_EXTERN PetscMPIInt MPIU_Allreduce_Private(const void *, void *, MPIU_Count, MPI_Datatype, MPI_Op, MPI_Comm);
2891: PETSC_EXTERN PetscErrorCode PetscCheckAllreduceSameLineAndCount_Private(MPI_Comm, const char *, PetscMPIInt, PetscMPIInt);
2893: #if defined(PETSC_USE_DEBUG)
2894: static inline unsigned int PetscStrHash(const char *str)
2895: {
2896: unsigned int c, hash = 5381;
2898: while ((c = (unsigned int)*str++)) hash = ((hash << 5) + hash) + c; /* hash * 33 + c */
2899: return hash;
2900: }
2901: #endif
2903: /*MC
2904: MPIU_Allreduce - A replacement for `MPI_Allreduce()` that (1) performs single-count `MPIU_INT` operations in `PetscInt64` to detect
2905: integer overflows and (2) tries to determine if the call from all the MPI ranks occur in the
2906: same place in the PETSc code. This helps to detect bugs where different MPI ranks follow different code paths
2907: resulting in inconsistent and incorrect calls to `MPI_Allreduce()`.
2909: Synopsis:
2910: #include <petscsys.h>
2911: PetscMPIInt MPIU_Allreduce(void *indata,void *outdata,PetscCount count,MPI_Datatype dtype, MPI_Op op, MPI_Comm comm);
2913: Collective
2915: Input Parameters:
2916: + a - pointer to the input data to be reduced
2917: . count - the number of MPI data items in `a` and `b`
2918: . dtype - the MPI datatype, for example `MPI_INT`
2919: . op - the MPI operation, for example `MPI_SUM`
2920: - comm - the MPI communicator on which the operation occurs
2922: Output Parameter:
2923: . b - the reduced values
2925: Level: developer
2927: Note:
2928: Should be wrapped with `PetscCallMPI()` for error checking
2930: .seealso: [](stylePetscCount), `MPI_Allreduce()`
2931: M*/
2932: #if defined(PETSC_USE_DEBUG)
2933: #define MPIU_Allreduce(a, b, count, dtype, op, comm) \
2934: PetscMacroReturnStandard( \
2935: PetscCall(PetscCheckAllreduceSameLineAndCount_Private((comm), __FILE__, (PetscMPIInt)__LINE__, (PetscMPIInt)(count))); \
2936: PetscCallMPI(MPIU_Allreduce_Private((a), (b), (count), (dtype), (op), (comm)));)
2937: #else
2938: #define MPIU_Allreduce(a, b, count, dtype, op, comm) MPIU_Allreduce_Private((a), (b), (count), (dtype), (op), (comm))
2939: #endif
2941: /* this is a vile hack */
2942: #if defined(PETSC_HAVE_NECMPI)
2943: #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)
2944: #define MPI_Type_free(a) (*(a) = MPI_DATATYPE_NULL, 0);
2945: #endif
2946: #endif
2948: /*
2949: List of external packages and queries on it
2950: */
2951: PETSC_EXTERN PetscErrorCode PetscHasExternalPackage(const char[], PetscBool *);
2953: /* this cannot go here because it may be in a different shared library */
2954: PETSC_EXTERN PetscErrorCode PCMPIServerBegin(void);
2955: PETSC_EXTERN PetscErrorCode PCMPIServerEnd(void);
2956: PETSC_EXTERN PetscBool PCMPIServerActive;
2957: PETSC_EXTERN PetscBool PCMPIServerInSolve;
2958: PETSC_EXTERN PetscBool PCMPIServerUseShmget;
2959: PETSC_EXTERN PetscErrorCode PetscShmgetAllocateArray(size_t, size_t, void **);
2960: PETSC_EXTERN PetscErrorCode PetscShmgetDeallocateArray(void **);
2961: PETSC_EXTERN PetscErrorCode PetscShmgetMapAddresses(MPI_Comm, PetscInt, const void **, void **);
2962: PETSC_EXTERN PetscErrorCode PetscShmgetUnmapAddresses(PetscInt, void **);
2963: PETSC_EXTERN PetscErrorCode PetscShmgetAddressesFinalize(void);
2965: typedef struct {
2966: PetscInt n;
2967: void *addr[3];
2968: } PCMPIServerAddresses;
2969: PETSC_EXTERN PetscErrorCode PCMPIServerAddressesDestroy(PetscCtxRt);
2971: #define PETSC_HAVE_FORTRAN PETSC_DEPRECATED_MACRO(3, 20, 0, "PETSC_USE_FORTRAN_BINDINGS", ) PETSC_USE_FORTRAN_BINDINGS
2973: PETSC_EXTERN PetscErrorCode PetscBLASSetNumThreads(PetscInt);
2974: PETSC_EXTERN PetscErrorCode PetscBLASGetNumThreads(PetscInt *);
2976: /*MC
2977: PetscSafePointerPlusOffset - Checks that a pointer is not `NULL` before applying an offset
2979: Level: beginner
2981: Note:
2982: This is needed to avoid errors with undefined-behavior sanitizers such as
2983: UBSan, assuming PETSc has been configured with `-fsanitize=undefined` as part of the compiler flags
2984: M*/
2985: #define PetscSafePointerPlusOffset(ptr, offset) ((ptr) ? (ptr) + (offset) : NULL)
2987: /* this is required to force PetscDevice to be visible at the system level for the Fortran interface */
2988: #include <petscdevicetypes.h>
2990: #if PetscDefined(USE_DEBUG) && !PetscDefined(HAVE_THREADSAFETY)
2991: PETSC_EXTERN PetscErrorCode PetscStackView(FILE *);
2992: #else
2993: #define PetscStackView(file) PETSC_SUCCESS
2994: #endif