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 8 9 /* ========================================================================== */ 10 /* 11 petscconf.h is contained in ${PETSC_ARCH}/include/petscconf.h it is 12 found automatically by the compiler due to the -I${PETSC_DIR}/${PETSC_ARCH}/include that 13 PETSc's makefiles add to the compiler rules. 14 For --prefix installs the directory ${PETSC_ARCH} does not exist and petscconf.h is in the same 15 directory as the other PETSc include files. 16 */ 17 #include <petscconf.h> 18 #include <petscpkg_version.h> 19 #include <petscconf_poison.h> 20 #include <petscfix.h> 21 #include <petscmacros.h> 22 23 /* SUBMANSEC = Sys */ 24 25 #if defined(PETSC_DESIRE_FEATURE_TEST_MACROS) 26 /* 27 Feature test macros must be included before headers defined by IEEE Std 1003.1-2001 28 We only turn these in PETSc source files that require them by setting PETSC_DESIRE_FEATURE_TEST_MACROS 29 */ 30 #if defined(PETSC__POSIX_C_SOURCE_200112L) && !defined(_POSIX_C_SOURCE) 31 #define _POSIX_C_SOURCE 200112L 32 #endif 33 #if defined(PETSC__BSD_SOURCE) && !defined(_BSD_SOURCE) 34 #define _BSD_SOURCE 35 #endif 36 #if defined(PETSC__DEFAULT_SOURCE) && !defined(_DEFAULT_SOURCE) 37 #define _DEFAULT_SOURCE 38 #endif 39 #if defined(PETSC__GNU_SOURCE) && !defined(_GNU_SOURCE) 40 #define _GNU_SOURCE 41 #endif 42 #endif 43 44 #include <petscsystypes.h> 45 46 /* ========================================================================== */ 47 48 /* 49 Defines the interface to MPI allowing the use of all MPI functions. 50 51 PETSc does not use the C++ binding of MPI at ALL. The following flag 52 makes sure the C++ bindings are not included. The C++ bindings REQUIRE 53 putting mpi.h before ANY C++ include files, we cannot control this 54 with all PETSc users. Users who want to use the MPI C++ bindings can include 55 mpicxx.h directly in their code 56 */ 57 #if !defined(MPICH_SKIP_MPICXX) 58 #define MPICH_SKIP_MPICXX 1 59 #endif 60 #if !defined(OMPI_SKIP_MPICXX) 61 #define OMPI_SKIP_MPICXX 1 62 #endif 63 #if defined(PETSC_HAVE_MPIUNI) 64 #include <petsc/mpiuni/mpi.h> 65 #else 66 #include <mpi.h> 67 #endif 68 69 /* 70 Perform various sanity checks that the correct mpi.h is being included at compile time. 71 This usually happens because 72 * either an unexpected mpi.h is in the default compiler path (i.e. in /usr/include) or 73 * an extra include path -I/something (which contains the unexpected mpi.h) is being passed to the compiler 74 */ 75 #if defined(PETSC_HAVE_MPIUNI) 76 #ifndef MPIUNI_H 77 #error "PETSc was configured with --with-mpi=0 but now appears to be compiling using a different mpi.h" 78 #endif 79 #elif defined(PETSC_HAVE_I_MPI) 80 #if !defined(I_MPI_NUMVERSION) 81 #error "PETSc was configured with I_MPI but now appears to be compiling using a non-I_MPI mpi.h" 82 #elif I_MPI_NUMVERSION != PETSC_PKG_I_MPI_NUMVERSION 83 #error "PETSc was configured with one I_MPI mpi.h version but now appears to be compiling using a different I_MPI mpi.h version" 84 #endif 85 #elif defined(PETSC_HAVE_MVAPICH2) 86 #if !defined(MVAPICH2_NUMVERSION) 87 #error "PETSc was configured with MVAPICH2 but now appears to be compiling using a non-MVAPICH2 mpi.h" 88 #elif MVAPICH2_NUMVERSION != PETSC_PKG_MVAPICH2_NUMVERSION 89 #error "PETSc was configured with one MVAPICH2 mpi.h version but now appears to be compiling using a different MVAPICH2 mpi.h version" 90 #endif 91 #elif defined(PETSC_HAVE_MPICH) 92 #if !defined(MPICH_NUMVERSION) || defined(MVAPICH2_NUMVERSION) || defined(I_MPI_NUMVERSION) 93 #error "PETSc was configured with MPICH but now appears to be compiling using a non-MPICH mpi.h" 94 #elif !PETSC_PKG_MPICH_VERSION_EQ(MPICH_NUMVERSION / 10000000, MPICH_NUMVERSION / 100000 % 100, MPICH_NUMVERSION / 1000 % 100) 95 #error "PETSc was configured with one MPICH mpi.h version but now appears to be compiling using a different MPICH mpi.h version" 96 #endif 97 #elif defined(PETSC_HAVE_OPENMPI) 98 #if !defined(OMPI_MAJOR_VERSION) 99 #error "PETSc was configured with Open MPI but now appears to be compiling using a non-Open MPI mpi.h" 100 #elif !PETSC_PKG_OPENMPI_VERSION_EQ(OMPI_MAJOR_VERSION, OMPI_MINOR_VERSION, OMPI_RELEASE_VERSION) 101 #error "PETSc was configured with one Open MPI mpi.h version but now appears to be compiling using a different Open MPI mpi.h version" 102 #endif 103 #elif defined(PETSC_HAVE_MSMPI_VERSION) 104 #if !defined(MSMPI_VER) 105 #error "PETSc was configured with MSMPI but now appears to be compiling using a non-MSMPI mpi.h" 106 #elif (MSMPI_VER != PETSC_HAVE_MSMPI_VERSION) 107 #error "PETSc was configured with one MSMPI mpi.h version but now appears to be compiling using a different MSMPI mpi.h version" 108 #endif 109 #elif defined(OMPI_MAJOR_VERSION) || defined(MPICH_NUMVERSION) || defined(MSMPI_VER) 110 #error "PETSc was configured with undetermined MPI - but now appears to be compiling using any of Open MPI, MS-MPI or a MPICH variant" 111 #endif 112 113 /* 114 Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler 115 see the top of mpicxx.h in the MPICH2 distribution. 116 */ 117 #include <stdio.h> 118 119 /* MSMPI on 32-bit Microsoft Windows requires this yukky hack - that breaks MPI standard compliance */ 120 #if !defined(MPIAPI) 121 #define MPIAPI 122 #endif 123 124 PETSC_EXTERN MPI_Datatype MPIU_ENUM PETSC_ATTRIBUTE_MPI_TYPE_TAG(PetscEnum); 125 PETSC_EXTERN MPI_Datatype MPIU_BOOL PETSC_ATTRIBUTE_MPI_TYPE_TAG(PetscBool); 126 127 /*MC 128 MPIU_INT - Portable MPI datatype corresponding to `PetscInt` independent of the precision of `PetscInt` 129 130 Level: beginner 131 132 Note: 133 In MPI calls that require an MPI datatype that matches a `PetscInt` or array of `PetscInt` values, pass this value. 134 135 .seealso: `PetscReal`, `PetscScalar`, `PetscComplex`, `PetscInt`, `MPIU_COUNT`, `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX` 136 M*/ 137 138 PETSC_EXTERN MPI_Datatype MPIU_FORTRANADDR; 139 140 #if defined(PETSC_USE_64BIT_INDICES) 141 #define MPIU_INT MPIU_INT64 142 #else 143 #define MPIU_INT MPI_INT 144 #endif 145 146 /*MC 147 MPIU_COUNT - Portable MPI datatype corresponding to `PetscCount` independent of the precision of `PetscCount` 148 149 Level: beginner 150 151 Note: 152 In MPI calls that require an MPI datatype that matches a `PetscCount` or array of `PetscCount` values, pass this value. 153 154 Developer Note: 155 It seems `MPI_AINT` is unsigned so this may be the wrong choice here since `PetscCount` is signed 156 157 .seealso: `PetscReal`, `PetscScalar`, `PetscComplex`, `PetscInt`, `MPIU_INT`, `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX` 158 M*/ 159 #define MPIU_COUNT MPI_AINT 160 161 /* 162 For the rare cases when one needs to send a size_t object with MPI 163 */ 164 PETSC_EXTERN MPI_Datatype MPIU_SIZE_T PETSC_ATTRIBUTE_MPI_TYPE_TAG(size_t); 165 166 /* 167 You can use PETSC_STDOUT as a replacement of stdout. You can also change 168 the value of PETSC_STDOUT to redirect all standard output elsewhere 169 */ 170 PETSC_EXTERN FILE *PETSC_STDOUT; 171 172 /* 173 You can use PETSC_STDERR as a replacement of stderr. You can also change 174 the value of PETSC_STDERR to redirect all standard error elsewhere 175 */ 176 PETSC_EXTERN FILE *PETSC_STDERR; 177 178 /* 179 Handle inclusion when using clang compiler with CUDA support 180 __float128 is not available for the device 181 */ 182 #if defined(__clang__) && (defined(__CUDA_ARCH__) || defined(__HIPCC__)) 183 #define PETSC_SKIP_REAL___FLOAT128 184 #endif 185 186 /* 187 Declare extern C stuff after including external header files 188 */ 189 190 PETSC_EXTERN PetscBool PETSC_RUNNING_ON_VALGRIND; 191 /* 192 Defines elementary mathematics functions and constants. 193 */ 194 #include <petscmath.h> 195 196 /*MC 197 PETSC_IGNORE - same as `NULL`, means PETSc will ignore this argument 198 199 Level: beginner 200 201 Note: 202 Accepted by many PETSc functions to not set a parameter and instead use a default value 203 204 Fortran Note: 205 Use `PETSC_NULL_INTEGER`, `PETSC_NULL_DOUBLE_PRECISION` etc 206 207 .seealso: `PETSC_DECIDE`, `PETSC_DEFAULT`, `PETSC_DETERMINE` 208 M*/ 209 #define PETSC_IGNORE PETSC_NULLPTR 210 #define PETSC_NULL PETSC_DEPRECATED_MACRO(3, 19, 0, "PETSC_NULLPTR", ) PETSC_NULLPTR 211 212 /*MC 213 PETSC_UNLIMITED - standard way of passing an integer or floating point parameter to indicate PETSc there is no bound on the value allowed 214 215 Level: beginner 216 217 Example Usage: 218 .vb 219 KSPSetTolerances(ksp, PETSC_CURRENT, PETSC_CURRENT, PETSC_UNLIMITED, PETSC_UNLIMITED); 220 .ve 221 indicates that the solver is allowed to take any number of iterations and will not stop early no matter how the residual gets. 222 223 Fortran Note: 224 Use `PETSC_UNLIMITED_INTEGER` or `PETSC_UNLIMITED_REAL`. 225 226 .seealso: `PETSC_DEFAULT`, `PETSC_IGNORE`, `PETSC_DETERMINE`, `PETSC_DECIDE` 227 M*/ 228 229 /*MC 230 PETSC_DECIDE - standard way of passing an integer or floating point parameter to indicate PETSc should determine an appropriate value 231 232 Level: beginner 233 234 Example Usage: 235 .vb 236 VecSetSizes(ksp, PETSC_DECIDE, 10); 237 .ve 238 indicates that the global size of the vector is 10 and the local size will be automatically determined so that the sum of the 239 local sizes is the global size, see `PetscSplitOwnership()`. 240 241 Fortran Note: 242 Use `PETSC_DECIDE_INTEGER` or `PETSC_DECIDE_REAL`. 243 244 .seealso: `PETSC_DEFAULT`, `PETSC_IGNORE`, `PETSC_DETERMINE`, `PETSC_UNLIMITED' 245 M*/ 246 247 /*MC 248 PETSC_DETERMINE - standard way of passing an integer or floating point parameter to indicate PETSc should determine an appropriate value 249 250 Level: beginner 251 252 Example Usage: 253 .vb 254 VecSetSizes(ksp, 10, PETSC_DETERMINE); 255 .ve 256 indicates that the local size of the vector is 10 and the global size will be automatically summing up all the local sizes. 257 258 Note: 259 Same as `PETSC_DECIDE` 260 261 Fortran Note: 262 Use `PETSC_DETERMINE_INTEGER` or `PETSC_DETERMINE_REAL`. 263 264 Developer Note: 265 I would like to use const `PetscInt` `PETSC_DETERMINE` = `PETSC_DECIDE`; but for 266 some reason this is not allowed by the standard even though `PETSC_DECIDE` is a constant value. 267 268 .seealso: `PETSC_DECIDE`, `PETSC_DEFAULT`, `PETSC_IGNORE`, `VecSetSizes()`, `PETSC_UNLIMITED' 269 M*/ 270 271 /*MC 272 PETSC_CURRENT - standard way of indicating to an object not to change the current value of the parameter in the object 273 274 Level: beginner 275 276 Note: 277 Use `PETSC_DECIDE` to use the value that was set by PETSc when the object's type was set 278 279 Fortran Note: 280 Use `PETSC_CURRENT_INTEGER` or `PETSC_CURRENT_REAL`. 281 282 .seealso: `PETSC_DECIDE`, `PETSC_IGNORE`, `PETSC_DETERMINE`, `PETSC_DEFAULT`, `PETSC_UNLIMITED' 283 M*/ 284 285 /*MC 286 PETSC_DEFAULT - deprecated, see `PETSC_CURRENT` and `PETSC_DETERMINE` 287 288 Level: beginner 289 290 Note: 291 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. 292 293 Developer Note: 294 Unfortunately this was used for two different purposes in the past, to actually trigger the use of a default value or to continue the 295 use of currently set value (in, for example, `KSPSetTolerances()`. 296 297 .seealso: `PETSC_DECIDE`, `PETSC_IGNORE`, `PETSC_DETERMINE`, `PETSC_CURRENT`, `PETSC_UNLIMITED' 298 M*/ 299 300 /* These MUST be preprocessor defines! see https://gitlab.com/petsc/petsc/-/issues/1370 */ 301 #define PETSC_DECIDE (-1) 302 #define PETSC_DETERMINE PETSC_DECIDE 303 #define PETSC_CURRENT (-2) 304 #define PETSC_UNLIMITED (-3) 305 /* PETSC_DEFAULT is deprecated in favor of PETSC_CURRENT for use in KSPSetTolerances() and similar functions */ 306 #define PETSC_DEFAULT PETSC_CURRENT 307 308 /*MC 309 PETSC_COMM_WORLD - the equivalent of the `MPI_COMM_WORLD` communicator which represents all the processes that PETSc knows about. 310 311 Level: beginner 312 313 Notes: 314 By default `PETSC_COMM_WORLD` and `MPI_COMM_WORLD` are identical unless you wish to 315 run PETSc on ONLY a subset of `MPI_COMM_WORLD`. In that case create your new (smaller) 316 communicator, call it, say comm, and set `PETSC_COMM_WORLD` = comm BEFORE calling 317 `PetscInitialize()`, but after `MPI_Init()` has been called. 318 319 The value of `PETSC_COMM_WORLD` should never be used or accessed before `PetscInitialize()` 320 is called because it may not have a valid value yet. 321 322 .seealso: `PETSC_COMM_SELF` 323 M*/ 324 PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD; 325 326 /*MC 327 PETSC_COMM_SELF - This is always `MPI_COMM_SELF` 328 329 Level: beginner 330 331 Note: 332 Do not USE/access or set this variable before `PetscInitialize()` has been called. 333 334 .seealso: `PETSC_COMM_WORLD` 335 M*/ 336 #define PETSC_COMM_SELF MPI_COMM_SELF 337 338 /*MC 339 PETSC_MPI_THREAD_REQUIRED - the required threading support used if PETSc initializes MPI with `MPI_Init_thread()`. 340 341 No Fortran Support 342 343 Level: beginner 344 345 Note: 346 By default `PETSC_MPI_THREAD_REQUIRED` equals `MPI_THREAD_FUNNELED` when the MPI implementation provides `MPI_Init_thread()`, otherwise it equals `MPI_THREAD_SINGLE` 347 348 .seealso: `PetscInitialize()` 349 M*/ 350 PETSC_EXTERN PetscMPIInt PETSC_MPI_THREAD_REQUIRED; 351 352 /*MC 353 PetscBeganMPI - indicates if PETSc initialized MPI during `PetscInitialize()` or if MPI was already initialized. 354 355 Synopsis: 356 #include <petscsys.h> 357 PetscBool PetscBeganMPI; 358 359 No Fortran Support 360 361 Level: developer 362 363 .seealso: `PetscInitialize()`, `PetscInitializeCalled` 364 M*/ 365 PETSC_EXTERN PetscBool PetscBeganMPI; 366 367 PETSC_EXTERN PetscBool PetscErrorHandlingInitialized; 368 PETSC_EXTERN PetscBool PetscInitializeCalled; 369 PETSC_EXTERN PetscBool PetscFinalizeCalled; 370 PETSC_EXTERN PetscBool PetscViennaCLSynchronize; 371 372 PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm), PetscErrorCode (*)(MPI_Comm)); 373 PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm, MPI_Comm *, int *); 374 PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm *); 375 PETSC_EXTERN PetscErrorCode PetscCommGetComm(MPI_Comm, MPI_Comm *); 376 PETSC_EXTERN PetscErrorCode PetscCommRestoreComm(MPI_Comm, MPI_Comm *); 377 378 #if defined(PETSC_HAVE_KOKKOS) 379 PETSC_EXTERN PetscErrorCode PetscKokkosInitializeCheck(void); /* Initialize Kokkos if not yet. */ 380 #endif 381 382 #if defined(PETSC_HAVE_NVSHMEM) 383 PETSC_EXTERN PetscBool PetscBeganNvshmem; 384 PETSC_EXTERN PetscBool PetscNvshmemInitialized; 385 PETSC_EXTERN PetscErrorCode PetscNvshmemFinalize(void); 386 #endif 387 388 #if defined(PETSC_HAVE_ELEMENTAL) 389 PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void); 390 PETSC_EXTERN PetscErrorCode PetscElementalInitialized(PetscBool *); 391 PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void); 392 #endif 393 394 /*MC 395 PetscMalloc - Allocates memory, One should use `PetscNew()`, `PetscMalloc1()` or `PetscCalloc1()` usually instead of this 396 397 Synopsis: 398 #include <petscsys.h> 399 PetscErrorCode PetscMalloc(size_t m,void **result) 400 401 Not Collective 402 403 Input Parameter: 404 . m - number of bytes to allocate 405 406 Output Parameter: 407 . result - memory allocated 408 409 Level: beginner 410 411 Notes: 412 Memory is always allocated at least double aligned 413 414 It is safe to allocate size 0 and pass the resulting pointer (which may or may not be `NULL`) to `PetscFree()`. 415 416 .seealso: `PetscFree()`, `PetscNew()` 417 M*/ 418 #define PetscMalloc(a, b) ((*PetscTrMalloc)((a), PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, (void **)(b))) 419 420 /*MC 421 PetscRealloc - Reallocates memory 422 423 Synopsis: 424 #include <petscsys.h> 425 PetscErrorCode PetscRealloc(size_t m,void **result) 426 427 Not Collective 428 429 Input Parameters: 430 + m - number of bytes to allocate 431 - result - previous memory 432 433 Output Parameter: 434 . result - new memory allocated 435 436 Level: developer 437 438 Note: 439 Memory is always allocated at least double aligned 440 441 .seealso: `PetscMalloc()`, `PetscFree()`, `PetscNew()` 442 M*/ 443 #define PetscRealloc(a, b) ((*PetscTrRealloc)((a), __LINE__, PETSC_FUNCTION_NAME, __FILE__, (void **)(b))) 444 445 /*MC 446 PetscAddrAlign - Rounds up an address to `PETSC_MEMALIGN` alignment 447 448 Synopsis: 449 #include <petscsys.h> 450 void *PetscAddrAlign(void *addr) 451 452 Not Collective 453 454 Input Parameter: 455 . addr - address to align (any pointer type) 456 457 Level: developer 458 459 .seealso: `PetscMallocAlign()` 460 M*/ 461 #define PetscAddrAlign(a) ((void *)((((PETSC_UINTPTR_T)(a)) + (PETSC_MEMALIGN - 1)) & ~(PETSC_MEMALIGN - 1))) 462 463 /*MC 464 PetscCalloc - Allocates a cleared (zeroed) memory region aligned to `PETSC_MEMALIGN` 465 466 Synopsis: 467 #include <petscsys.h> 468 PetscErrorCode PetscCalloc(size_t m,void **result) 469 470 Not Collective 471 472 Input Parameter: 473 . m - number of bytes to allocate 474 475 Output Parameter: 476 . result - memory allocated 477 478 Level: beginner 479 480 Notes: 481 Memory is always allocated at least double aligned. This macro is useful in allocating memory pointed by void pointers 482 483 It is safe to allocate size 0 and pass the resulting pointer (which may or may not be `NULL`) to `PetscFree()`. 484 485 .seealso: `PetscFree()`, `PetscNew()` 486 M*/ 487 #define PetscCalloc(m, result) PetscMallocA(1, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)m), (result)) 488 489 /*MC 490 PetscMalloc1 - Allocates an array of memory aligned to `PETSC_MEMALIGN` 491 492 Synopsis: 493 #include <petscsys.h> 494 PetscErrorCode PetscMalloc1(size_t m1,type **r1) 495 496 Not Collective 497 498 Input Parameter: 499 . m1 - number of elements to allocate (may be zero) 500 501 Output Parameter: 502 . r1 - memory allocated 503 504 Level: beginner 505 506 Note: 507 This uses `sizeof()` of the memory type requested to determine the total memory to be allocated; therefore, you should not 508 multiply the number of elements requested by the `sizeof()` the type. For example, use 509 .vb 510 PetscInt *id; 511 PetscMalloc1(10,&id); 512 .ve 513 not 514 .vb 515 PetscInt *id; 516 PetscMalloc1(10*sizeof(PetscInt),&id); 517 .ve 518 519 Does not zero the memory allocated, use `PetscCalloc1()` to obtain memory that has been zeroed. 520 521 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 522 32bit `int` computation - while computation in 64bit `size_t` would not overflow! 523 It's best if any arithmetic that is done for size computations is done with `size_t` type - avoiding arithmetic overflow! 524 525 `PetscMalloc[N]()` and `PetscCalloc[N]()` attempt to work-around this by casting the first variable to `size_t`. 526 This works for most expressions, but not all, such as 527 .vb 528 PetscInt *id, a, b; 529 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 530 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 531 .ve 532 533 These expressions should either be avoided, or appropriately cast variables to `size_t`: 534 .vb 535 PetscInt *id, a, b; 536 PetscMalloc1(use_a_squared ? (size_t)a * a * b : (size_t)a * b, &id); // Cast a to size_t before multiplication 537 PetscMalloc1(b * b + a, &id); // b is automatically cast to size_t and order-of-operations ensures size_t precision is maintained 538 .ve 539 540 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc2()` 541 M*/ 542 #define PetscMalloc1(m1, r1) PetscMallocA(1, PETSC_FALSE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1)) 543 544 /*MC 545 PetscCalloc1 - Allocates a cleared (zeroed) array of memory aligned to `PETSC_MEMALIGN` 546 547 Synopsis: 548 #include <petscsys.h> 549 PetscErrorCode PetscCalloc1(size_t m1,type **r1) 550 551 Not Collective 552 553 Input Parameter: 554 . m1 - number of elements to allocate in 1st chunk (may be zero) 555 556 Output Parameter: 557 . r1 - memory allocated 558 559 Level: beginner 560 561 Note: 562 See `PetsMalloc1()` for more details on usage. 563 564 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc2()` 565 M*/ 566 #define PetscCalloc1(m1, r1) PetscMallocA(1, PETSC_TRUE, __LINE__, PETSC_FUNCTION_NAME, __FILE__, ((size_t)((size_t)m1) * sizeof(**(r1))), (r1)) 567 568 /*MC 569 PetscMalloc2 - Allocates 2 arrays of memory both aligned to `PETSC_MEMALIGN` 570 571 Synopsis: 572 #include <petscsys.h> 573 PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2) 574 575 Not Collective 576 577 Input Parameters: 578 + m1 - number of elements to allocate in 1st chunk (may be zero) 579 - m2 - number of elements to allocate in 2nd chunk (may be zero) 580 581 Output Parameters: 582 + r1 - memory allocated in first chunk 583 - r2 - memory allocated in second chunk 584 585 Level: developer 586 587 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc2()` 588 M*/ 589 #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)) 590 591 /*MC 592 PetscCalloc2 - Allocates 2 cleared (zeroed) arrays of memory both aligned to `PETSC_MEMALIGN` 593 594 Synopsis: 595 #include <petscsys.h> 596 PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2) 597 598 Not Collective 599 600 Input Parameters: 601 + m1 - number of elements to allocate in 1st chunk (may be zero) 602 - m2 - number of elements to allocate in 2nd chunk (may be zero) 603 604 Output Parameters: 605 + r1 - memory allocated in first chunk 606 - r2 - memory allocated in second chunk 607 608 Level: developer 609 610 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc2()` 611 M*/ 612 #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)) 613 614 /*MC 615 PetscMalloc3 - Allocates 3 arrays of memory, all aligned to `PETSC_MEMALIGN` 616 617 Synopsis: 618 #include <petscsys.h> 619 PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3) 620 621 Not Collective 622 623 Input Parameters: 624 + m1 - number of elements to allocate in 1st chunk (may be zero) 625 . m2 - number of elements to allocate in 2nd chunk (may be zero) 626 - m3 - number of elements to allocate in 3rd chunk (may be zero) 627 628 Output Parameters: 629 + r1 - memory allocated in first chunk 630 . r2 - memory allocated in second chunk 631 - r3 - memory allocated in third chunk 632 633 Level: developer 634 635 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc3()`, `PetscFree3()` 636 M*/ 637 #define PetscMalloc3(m1, r1, m2, r2, m3, r3) \ 638 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)) 639 640 /*MC 641 PetscCalloc3 - Allocates 3 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN` 642 643 Synopsis: 644 #include <petscsys.h> 645 PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3) 646 647 Not Collective 648 649 Input Parameters: 650 + m1 - number of elements to allocate in 1st chunk (may be zero) 651 . m2 - number of elements to allocate in 2nd chunk (may be zero) 652 - m3 - number of elements to allocate in 3rd chunk (may be zero) 653 654 Output Parameters: 655 + r1 - memory allocated in first chunk 656 . r2 - memory allocated in second chunk 657 - r3 - memory allocated in third chunk 658 659 Level: developer 660 661 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscCalloc2()`, `PetscMalloc3()`, `PetscFree3()` 662 M*/ 663 #define PetscCalloc3(m1, r1, m2, r2, m3, r3) \ 664 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)) 665 666 /*MC 667 PetscMalloc4 - Allocates 4 arrays of memory, all aligned to `PETSC_MEMALIGN` 668 669 Synopsis: 670 #include <petscsys.h> 671 PetscErrorCode PetscMalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4) 672 673 Not Collective 674 675 Input Parameters: 676 + m1 - number of elements to allocate in 1st chunk (may be zero) 677 . m2 - number of elements to allocate in 2nd chunk (may be zero) 678 . m3 - number of elements to allocate in 3rd chunk (may be zero) 679 - m4 - number of elements to allocate in 4th chunk (may be zero) 680 681 Output Parameters: 682 + r1 - memory allocated in first chunk 683 . r2 - memory allocated in second chunk 684 . r3 - memory allocated in third chunk 685 - r4 - memory allocated in fourth chunk 686 687 Level: developer 688 689 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc4()`, `PetscFree4()` 690 M*/ 691 #define PetscMalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \ 692 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)) 693 694 /*MC 695 PetscCalloc4 - Allocates 4 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN` 696 697 Synopsis: 698 #include <petscsys.h> 699 PetscErrorCode PetscCalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4) 700 701 Not Collective 702 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) 707 - m4 - number of elements to allocate in 4th chunk (may be zero) 708 709 Output Parameters: 710 + r1 - memory allocated in first chunk 711 . r2 - memory allocated in second chunk 712 . r3 - memory allocated in third chunk 713 - r4 - memory allocated in fourth chunk 714 715 Level: developer 716 717 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc4()`, `PetscFree4()` 718 M*/ 719 #define PetscCalloc4(m1, r1, m2, r2, m3, r3, m4, r4) \ 720 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)) 721 722 /*MC 723 PetscMalloc5 - Allocates 5 arrays of memory, all aligned to `PETSC_MEMALIGN` 724 725 Synopsis: 726 #include <petscsys.h> 727 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) 728 729 Not Collective 730 731 Input Parameters: 732 + m1 - number of elements to allocate in 1st chunk (may be zero) 733 . m2 - number of elements to allocate in 2nd chunk (may be zero) 734 . m3 - number of elements to allocate in 3rd chunk (may be zero) 735 . m4 - number of elements to allocate in 4th chunk (may be zero) 736 - m5 - number of elements to allocate in 5th chunk (may be zero) 737 738 Output Parameters: 739 + r1 - memory allocated in first chunk 740 . r2 - memory allocated in second chunk 741 . r3 - memory allocated in third chunk 742 . r4 - memory allocated in fourth chunk 743 - r5 - memory allocated in fifth chunk 744 745 Level: developer 746 747 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc5()`, `PetscFree5()` 748 M*/ 749 #define PetscMalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \ 750 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)) 751 752 /*MC 753 PetscCalloc5 - Allocates 5 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN` 754 755 Synopsis: 756 #include <petscsys.h> 757 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) 758 759 Not Collective 760 761 Input Parameters: 762 + m1 - number of elements to allocate in 1st chunk (may be zero) 763 . m2 - number of elements to allocate in 2nd chunk (may be zero) 764 . m3 - number of elements to allocate in 3rd chunk (may be zero) 765 . m4 - number of elements to allocate in 4th chunk (may be zero) 766 - m5 - number of elements to allocate in 5th chunk (may be zero) 767 768 Output Parameters: 769 + r1 - memory allocated in first chunk 770 . r2 - memory allocated in second chunk 771 . r3 - memory allocated in third chunk 772 . r4 - memory allocated in fourth chunk 773 - r5 - memory allocated in fifth chunk 774 775 Level: developer 776 777 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc5()`, `PetscFree5()` 778 M*/ 779 #define PetscCalloc5(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5) \ 780 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)) 781 782 /*MC 783 PetscMalloc6 - Allocates 6 arrays of memory, all aligned to `PETSC_MEMALIGN` 784 785 Synopsis: 786 #include <petscsys.h> 787 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) 788 789 Not Collective 790 791 Input Parameters: 792 + m1 - number of elements to allocate in 1st chunk (may be zero) 793 . m2 - number of elements to allocate in 2nd chunk (may be zero) 794 . m3 - number of elements to allocate in 3rd chunk (may be zero) 795 . m4 - number of elements to allocate in 4th chunk (may be zero) 796 . m5 - number of elements to allocate in 5th chunk (may be zero) 797 - m6 - number of elements to allocate in 6th chunk (may be zero) 798 799 Output Parameteasr: 800 + r1 - memory allocated in first chunk 801 . r2 - memory allocated in second chunk 802 . r3 - memory allocated in third chunk 803 . r4 - memory allocated in fourth chunk 804 . r5 - memory allocated in fifth chunk 805 - r6 - memory allocated in sixth chunk 806 807 Level: developer 808 809 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc6()`, `PetscFree3()`, `PetscFree4()`, `PetscFree5()`, `PetscFree6()` 810 M*/ 811 #define PetscMalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \ 812 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)) 813 814 /*MC 815 PetscCalloc6 - Allocates 6 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN` 816 817 Synopsis: 818 #include <petscsys.h> 819 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) 820 821 Not Collective 822 823 Input Parameters: 824 + m1 - number of elements to allocate in 1st chunk (may be zero) 825 . m2 - number of elements to allocate in 2nd chunk (may be zero) 826 . m3 - number of elements to allocate in 3rd chunk (may be zero) 827 . m4 - number of elements to allocate in 4th chunk (may be zero) 828 . m5 - number of elements to allocate in 5th chunk (may be zero) 829 - m6 - number of elements to allocate in 6th chunk (may be zero) 830 831 Output Parameters: 832 + r1 - memory allocated in first chunk 833 . r2 - memory allocated in second chunk 834 . r3 - memory allocated in third chunk 835 . r4 - memory allocated in fourth chunk 836 . r5 - memory allocated in fifth chunk 837 - r6 - memory allocated in sixth chunk 838 839 Level: developer 840 841 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscMalloc6()`, `PetscFree6()` 842 M*/ 843 #define PetscCalloc6(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6) \ 844 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)) 845 846 /*MC 847 PetscMalloc7 - Allocates 7 arrays of memory, all aligned to `PETSC_MEMALIGN` 848 849 Synopsis: 850 #include <petscsys.h> 851 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) 852 853 Not Collective 854 855 Input Parameters: 856 + m1 - number of elements to allocate in 1st chunk (may be zero) 857 . m2 - number of elements to allocate in 2nd chunk (may be zero) 858 . m3 - number of elements to allocate in 3rd chunk (may be zero) 859 . m4 - number of elements to allocate in 4th chunk (may be zero) 860 . m5 - number of elements to allocate in 5th chunk (may be zero) 861 . m6 - number of elements to allocate in 6th chunk (may be zero) 862 - m7 - number of elements to allocate in 7th chunk (may be zero) 863 864 Output Parameters: 865 + r1 - memory allocated in first chunk 866 . r2 - memory allocated in second chunk 867 . r3 - memory allocated in third chunk 868 . r4 - memory allocated in fourth chunk 869 . r5 - memory allocated in fifth chunk 870 . r6 - memory allocated in sixth chunk 871 - r7 - memory allocated in seventh chunk 872 873 Level: developer 874 875 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscCalloc7()`, `PetscFree7()` 876 M*/ 877 #define PetscMalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \ 878 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)) 879 880 /*MC 881 PetscCalloc7 - Allocates 7 cleared (zeroed) arrays of memory, all aligned to `PETSC_MEMALIGN` 882 883 Synopsis: 884 #include <petscsys.h> 885 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) 886 887 Not Collective 888 889 Input Parameters: 890 + m1 - number of elements to allocate in 1st chunk (may be zero) 891 . m2 - number of elements to allocate in 2nd chunk (may be zero) 892 . m3 - number of elements to allocate in 3rd chunk (may be zero) 893 . m4 - number of elements to allocate in 4th chunk (may be zero) 894 . m5 - number of elements to allocate in 5th chunk (may be zero) 895 . m6 - number of elements to allocate in 6th chunk (may be zero) 896 - m7 - number of elements to allocate in 7th chunk (may be zero) 897 898 Output Parameters: 899 + r1 - memory allocated in first chunk 900 . r2 - memory allocated in second chunk 901 . r3 - memory allocated in third chunk 902 . r4 - memory allocated in fourth chunk 903 . r5 - memory allocated in fifth chunk 904 . r6 - memory allocated in sixth chunk 905 - r7 - memory allocated in seventh chunk 906 907 Level: developer 908 909 .seealso: `PetscFree()`, `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscMalloc7()`, `PetscFree7()` 910 M*/ 911 #define PetscCalloc7(m1, r1, m2, r2, m3, r3, m4, r4, m5, r5, m6, r6, m7, r7) \ 912 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)) 913 914 /*MC 915 PetscNew - Allocates memory of a particular type, zeros the memory! Aligned to `PETSC_MEMALIGN` 916 917 Synopsis: 918 #include <petscsys.h> 919 PetscErrorCode PetscNew(type **result) 920 921 Not Collective 922 923 Output Parameter: 924 . result - memory allocated, sized to match pointer type 925 926 Level: beginner 927 928 .seealso: `PetscFree()`, `PetscMalloc()`, `PetscCalloc1()`, `PetscMalloc1()` 929 M*/ 930 #define PetscNew(b) PetscCalloc1(1, (b)) 931 932 #define PetscNewLog(o, b) PETSC_DEPRECATED_MACRO(3, 18, 0, "PetscNew()", ) PetscNew(b) 933 934 /*MC 935 PetscFree - Frees memory 936 937 Synopsis: 938 #include <petscsys.h> 939 PetscErrorCode PetscFree(void *memory) 940 941 Not Collective 942 943 Input Parameter: 944 . memory - memory to free (the pointer is ALWAYS set to `NULL` upon success) 945 946 Level: beginner 947 948 Notes: 949 Do not free memory obtained with `PetscMalloc2()`, `PetscCalloc2()` etc, they must be freed with `PetscFree2()` etc. 950 951 It is safe to call `PetscFree()` on a `NULL` pointer. 952 953 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc1()`, `PetscCalloc1()` 954 M*/ 955 #define PetscFree(a) ((PetscErrorCode)((*PetscTrFree)((void *)(a), __LINE__, PETSC_FUNCTION_NAME, __FILE__) || ((a) = PETSC_NULLPTR, PETSC_SUCCESS))) 956 957 /*MC 958 PetscFree2 - Frees 2 chunks of memory obtained with `PetscMalloc2()` 959 960 Synopsis: 961 #include <petscsys.h> 962 PetscErrorCode PetscFree2(void *memory1,void *memory2) 963 964 Not Collective 965 966 Input Parameters: 967 + memory1 - memory to free 968 - memory2 - 2nd memory to free 969 970 Level: developer 971 972 Notes: 973 Memory must have been obtained with `PetscMalloc2()` 974 975 The arguments need to be in the same order as they were in the call to `PetscMalloc2()` 976 977 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()` 978 M*/ 979 #define PetscFree2(m1, m2) PetscFreeA(2, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2)) 980 981 /*MC 982 PetscFree3 - Frees 3 chunks of memory obtained with `PetscMalloc3()` 983 984 Synopsis: 985 #include <petscsys.h> 986 PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3) 987 988 Not Collective 989 990 Input Parameters: 991 + memory1 - memory to free 992 . memory2 - 2nd memory to free 993 - memory3 - 3rd memory to free 994 995 Level: developer 996 997 Notes: 998 Memory must have been obtained with `PetscMalloc3()` 999 1000 The arguments need to be in the same order as they were in the call to `PetscMalloc3()` 1001 1002 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()` 1003 M*/ 1004 #define PetscFree3(m1, m2, m3) PetscFreeA(3, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3)) 1005 1006 /*MC 1007 PetscFree4 - Frees 4 chunks of memory obtained with `PetscMalloc4()` 1008 1009 Synopsis: 1010 #include <petscsys.h> 1011 PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4) 1012 1013 Not Collective 1014 1015 Input Parameters: 1016 + m1 - memory to free 1017 . m2 - 2nd memory to free 1018 . m3 - 3rd memory to free 1019 - m4 - 4th memory to free 1020 1021 Level: developer 1022 1023 Notes: 1024 Memory must have been obtained with `PetscMalloc4()` 1025 1026 The arguments need to be in the same order as they were in the call to `PetscMalloc4()` 1027 1028 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()` 1029 M*/ 1030 #define PetscFree4(m1, m2, m3, m4) PetscFreeA(4, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4)) 1031 1032 /*MC 1033 PetscFree5 - Frees 5 chunks of memory obtained with `PetscMalloc5()` 1034 1035 Synopsis: 1036 #include <petscsys.h> 1037 PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5) 1038 1039 Not Collective 1040 1041 Input Parameters: 1042 + m1 - memory to free 1043 . m2 - 2nd memory to free 1044 . m3 - 3rd memory to free 1045 . m4 - 4th memory to free 1046 - m5 - 5th memory to free 1047 1048 Level: developer 1049 1050 Notes: 1051 Memory must have been obtained with `PetscMalloc5()` 1052 1053 The arguments need to be in the same order as they were in the call to `PetscMalloc5()` 1054 1055 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()` 1056 M*/ 1057 #define PetscFree5(m1, m2, m3, m4, m5) PetscFreeA(5, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5)) 1058 1059 /*MC 1060 PetscFree6 - Frees 6 chunks of memory obtained with `PetscMalloc6()` 1061 1062 Synopsis: 1063 #include <petscsys.h> 1064 PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6) 1065 1066 Not Collective 1067 1068 Input Parameters: 1069 + m1 - memory to free 1070 . m2 - 2nd memory to free 1071 . m3 - 3rd memory to free 1072 . m4 - 4th memory to free 1073 . m5 - 5th memory to free 1074 - m6 - 6th memory to free 1075 1076 Level: developer 1077 1078 Notes: 1079 Memory must have been obtained with `PetscMalloc6()` 1080 1081 The arguments need to be in the same order as they were in the call to `PetscMalloc6()` 1082 1083 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()` 1084 M*/ 1085 #define PetscFree6(m1, m2, m3, m4, m5, m6) PetscFreeA(6, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6)) 1086 1087 /*MC 1088 PetscFree7 - Frees 7 chunks of memory obtained with `PetscMalloc7()` 1089 1090 Synopsis: 1091 #include <petscsys.h> 1092 PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7) 1093 1094 Not Collective 1095 1096 Input Parameters: 1097 + m1 - memory to free 1098 . m2 - 2nd memory to free 1099 . m3 - 3rd memory to free 1100 . m4 - 4th memory to free 1101 . m5 - 5th memory to free 1102 . m6 - 6th memory to free 1103 - m7 - 7th memory to free 1104 1105 Level: developer 1106 1107 Notes: 1108 Memory must have been obtained with `PetscMalloc7()` 1109 1110 The arguments need to be in the same order as they were in the call to `PetscMalloc7()` 1111 1112 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()`, 1113 `PetscMalloc7()` 1114 M*/ 1115 #define PetscFree7(m1, m2, m3, m4, m5, m6, m7) PetscFreeA(7, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6), &(m7)) 1116 1117 PETSC_EXTERN PetscErrorCode PetscMallocA(int, PetscBool, int, const char *, const char *, size_t, void *, ...); 1118 PETSC_EXTERN PetscErrorCode PetscFreeA(int, int, const char *, const char *, void *, ...); 1119 PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t, PetscBool, int, const char[], const char[], void **); 1120 PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void *, int, const char[], const char[]); 1121 PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t, int, const char[], const char[], void **); 1122 PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool); 1123 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 **)); 1124 PETSC_EXTERN PetscErrorCode PetscMallocClear(void); 1125 1126 /* 1127 Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair. 1128 */ 1129 PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void); 1130 PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void); 1131 #if defined(PETSC_HAVE_CUDA) 1132 PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void); 1133 PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void); 1134 #endif 1135 #if defined(PETSC_HAVE_HIP) 1136 PETSC_EXTERN PetscErrorCode PetscMallocSetHIPHost(void); 1137 PETSC_EXTERN PetscErrorCode PetscMallocResetHIPHost(void); 1138 #endif 1139 1140 #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE 1141 #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION 1142 1143 /* 1144 Routines for tracing memory corruption/bleeding with default PETSc memory allocation 1145 */ 1146 PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *); 1147 PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *); 1148 PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *); 1149 PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *); 1150 PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int); 1151 PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int, PetscLogDouble *); 1152 PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool, PetscBool); 1153 PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool *, PetscBool *, PetscBool *); 1154 PETSC_EXTERN PetscErrorCode PetscMallocValidate(int, const char[], const char[]); 1155 PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble); 1156 PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool *); 1157 PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool); 1158 PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool *); 1159 1160 PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType, MPI_Datatype *); 1161 PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype, PetscDataType *); 1162 PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType, size_t *); 1163 PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char *, PetscDataType *, PetscBool *); 1164 1165 /* 1166 These are MPI operations for MPI_Allreduce() etc 1167 */ 1168 PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP; 1169 #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 1170 PETSC_EXTERN MPI_Op MPIU_SUM; 1171 PETSC_EXTERN MPI_Op MPIU_MAX; 1172 PETSC_EXTERN MPI_Op MPIU_MIN; 1173 #else 1174 #define MPIU_SUM MPI_SUM 1175 #define MPIU_MAX MPI_MAX 1176 #define MPIU_MIN MPI_MIN 1177 #endif 1178 PETSC_EXTERN MPI_Op Petsc_Garbage_SetIntersectOp; 1179 PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm, const PetscInt[], PetscInt *, PetscInt *); 1180 1181 #if (defined(PETSC_HAVE_REAL___FLOAT128) && !defined(PETSC_SKIP_REAL___FLOAT128)) || (defined(PETSC_HAVE_REAL___FP16) && !defined(PETSC_SKIP_REAL___FP16)) 1182 /*MC 1183 MPIU_SUM___FP16___FLOAT128 - MPI_Op that acts as a replacement for `MPI_SUM` with 1184 custom `MPI_Datatype` `MPIU___FLOAT128`, `MPIU___COMPLEX128`, and `MPIU___FP16`. 1185 1186 Level: advanced 1187 1188 Developer Note: 1189 This should be unified with `MPIU_SUM` 1190 1191 .seealso: `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX` 1192 M*/ 1193 PETSC_EXTERN MPI_Op MPIU_SUM___FP16___FLOAT128; 1194 #endif 1195 PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm, const PetscInt[], PetscInt *, PetscInt *); 1196 1197 PETSC_EXTERN PetscErrorCode MPIULong_Send(void *, PetscInt, MPI_Datatype, PetscMPIInt, PetscMPIInt, MPI_Comm) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(1, 3); 1198 PETSC_EXTERN PetscErrorCode MPIULong_Recv(void *, PetscInt, MPI_Datatype, PetscMPIInt, PetscMPIInt, MPI_Comm) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(1, 3); 1199 1200 /* 1201 Defines PETSc error handling. 1202 */ 1203 #include <petscerror.h> 1204 1205 PETSC_EXTERN PetscBool PetscCIEnabled; /* code is running in the PETSc test harness CI */ 1206 PETSC_EXTERN PetscBool PetscCIEnabledPortableErrorOutput; /* error output is stripped to ensure portability of error messages across systems */ 1207 PETSC_EXTERN const char *PetscCIFilename(const char *); 1208 PETSC_EXTERN int PetscCILinenumber(int); 1209 1210 #define PETSC_SMALLEST_CLASSID ((PetscClassId)1211211) 1211 PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID; 1212 PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID; 1213 PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[], PetscClassId *); 1214 PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject, PetscObjectId *); 1215 PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject, PetscObjectId, PetscBool *); 1216 1217 /* 1218 Routines that get memory usage information from the OS 1219 */ 1220 PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *); 1221 PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *); 1222 PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void); 1223 PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]); 1224 1225 PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal); 1226 1227 /* 1228 Initialization of PETSc 1229 */ 1230 PETSC_EXTERN PetscErrorCode PetscInitialize(int *, char ***, const char[], const char[]); 1231 PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int, char **, const char[], const char[]); 1232 PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void); 1233 PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *); 1234 PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *); 1235 PETSC_EXTERN PetscErrorCode PetscFinalize(void); 1236 PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void); 1237 PETSC_EXTERN PetscErrorCode PetscGetArgs(int *, char ***); 1238 PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***); 1239 PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **); 1240 1241 PETSC_EXTERN PetscErrorCode PetscEnd(void); 1242 PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void); 1243 PETSC_EXTERN PetscErrorCode PetscSysFinalizePackage(void); 1244 1245 PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[], const char[]); 1246 PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void); 1247 PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void); 1248 PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject, const char[]); 1249 1250 PETSC_EXTERN PetscErrorCode PetscMonitorCompare(PetscErrorCode (*)(void), void *, PetscErrorCode (*)(void **), PetscErrorCode (*)(void), void *, PetscErrorCode (*)(void **), PetscBool *); 1251 1252 /* 1253 These are so that in extern C code we can cast function pointers to non-extern C 1254 function pointers. Since the regular C++ code expects its function pointers to be C++ 1255 */ 1256 1257 /*S 1258 PetscVoidFn - A prototype of a void (fn)(void) function 1259 1260 Level: developer 1261 1262 Notes: 1263 The deprecated `PetscVoidFunction` works as a replacement for `PetscVoidFn` *. 1264 1265 The deprecated `PetscVoidStarFunction` works as a replacement for `PetscVoidFn` **. 1266 1267 .seealso: `PetscObject`, `PetscObjectDestroy()` 1268 S*/ 1269 PETSC_EXTERN_TYPEDEF typedef void(PetscVoidFn)(void); 1270 1271 PETSC_EXTERN_TYPEDEF typedef PetscVoidFn *PetscVoidFunction; 1272 PETSC_EXTERN_TYPEDEF typedef PetscVoidFn **PetscVoidStarFunction; 1273 1274 /*S 1275 PetscErrorCodeFn - A prototype of a PetscErrorCode (fn)(void) function 1276 1277 Level: developer 1278 1279 Notes: 1280 The deprecated `PetscErrorCodeFunction` works as a replacement for `PetscErrorCodeFn` *. 1281 1282 .seealso: `PetscObject`, `PetscObjectDestroy()` 1283 S*/ 1284 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(PetscErrorCodeFn)(void); 1285 1286 PETSC_EXTERN_TYPEDEF typedef PetscErrorCodeFn *PetscErrorCodeFunction; 1287 1288 /* 1289 Functions that can act on any PETSc object. 1290 */ 1291 PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject *); 1292 PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject, MPI_Comm *); 1293 PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject, PetscClassId *); 1294 PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject, const char *[]); 1295 PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject, const char *[]); 1296 PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject, const char[]); 1297 PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject, const char *[]); 1298 PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject, PetscInt); 1299 PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject, PetscInt *); 1300 PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject, PetscObject, PetscInt); 1301 PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject); 1302 PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject, PetscInt *); 1303 PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject); 1304 PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject, PetscMPIInt *); 1305 PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject, const char[], PetscObject); 1306 PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject, const char[]); 1307 PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject, const char[], PetscObject *); 1308 PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject, const char[], void (*)(void)); 1309 #define PetscObjectComposeFunction(a, b, ...) PetscObjectComposeFunction_Private((a), (b), (PetscVoidFn *)(__VA_ARGS__)) 1310 PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject); 1311 PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject); 1312 PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject); 1313 PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject, PetscObject); 1314 PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm, PetscMPIInt *); 1315 1316 /*MC 1317 PetscObjectParameterSetDefault - sets a parameter default value in a `PetscObject` to a new default value. 1318 If the current value matches the old default value, then the current value is also set to the new value. 1319 1320 No Fortran Support 1321 1322 Synopsis: 1323 #include <petscsys.h> 1324 PetscBool PetscObjectParameterSetDefault(PetscObject obj, char* NAME, PetscReal value); 1325 1326 Input Parameters: 1327 + obj - the `PetscObject` 1328 . NAME - the name of the parameter, unquoted 1329 - value - the new value 1330 1331 Level: developer 1332 1333 Notes: 1334 The defaults for an object are the values set when the object's type is set. 1335 1336 This should only be used in object constructors, such as, `SNESCreate_NGS()`. 1337 1338 This only works for parameters that are declared in the struct with `PetscObjectParameterDeclare()` 1339 1340 .seealso: `PetscObjectParameterDeclare()`, `PetscInitialize()`, `PetscFinalize()`, `PetscObject`, `SNESParametersInitialize()` 1341 M*/ 1342 #define PetscObjectParameterSetDefault(obj, NAME, value) \ 1343 do { \ 1344 if (obj->NAME == obj->default_##NAME) obj->NAME = value; \ 1345 obj->default_##NAME = value; \ 1346 } while (0) 1347 1348 /*MC 1349 PetscObjectParameterDeclare - declares a parameter in a `PetscObject` and a location to store its default 1350 1351 No Fortran Support 1352 1353 Synopsis: 1354 #include <petscsys.h> 1355 PetscBool PetscObjectParameterDeclare(type, char* NAME) 1356 1357 Input Parameters: 1358 + type - the type of the parameter, for example `PetscInt` 1359 - NAME - the name of the parameter, unquoted 1360 1361 Level: developer. 1362 1363 .seealso: `PetscObjectParameterSetDefault()`, `PetscInitialize()`, `PetscFinalize()`, `PetscObject`, `SNESParametersInitialize()` 1364 M*/ 1365 #define PetscObjectParameterDeclare(type, NAME) type NAME, default_##NAME 1366 1367 #include <petscviewertypes.h> 1368 #include <petscoptions.h> 1369 1370 PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer, PetscBool, PetscLogDouble); 1371 PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool *); 1372 1373 PETSC_EXTERN PetscErrorCode PetscObjectsListGetGlobalNumbering(MPI_Comm, PetscInt, PetscObject *, PetscInt *, PetscInt *); 1374 1375 PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer, const char[]); 1376 PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject, PetscViewer); 1377 PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject, PetscViewer); 1378 #define PetscObjectQueryFunction(obj, name, fptr) PetscObjectQueryFunction_Private((obj), (name), (PetscVoidFn **)(fptr)) 1379 PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject, const char[], void (**)(void)); 1380 PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject, const char[]); 1381 PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject, const char[]); 1382 PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject, const char[]); 1383 PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject, const char *[]); 1384 PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject, const char[]); 1385 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject); 1386 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void); 1387 PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject, PetscObject, const char[]); 1388 PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject); 1389 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject, const char[], PetscBool *); 1390 PETSC_EXTERN PetscErrorCode PetscObjectObjectTypeCompare(PetscObject, PetscObject, PetscBool *); 1391 PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject, const char[], PetscBool *); 1392 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject, PetscBool *, const char[], ...); 1393 PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject, PetscBool *, const char[], ...); 1394 PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void)); 1395 PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void); 1396 1397 #if defined(PETSC_HAVE_SAWS) 1398 PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void); 1399 PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject); 1400 PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject, PetscBool); 1401 PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject); 1402 PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject); 1403 PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject); 1404 PETSC_EXTERN void PetscStackSAWsGrantAccess(void); 1405 PETSC_EXTERN void PetscStackSAWsTakeAccess(void); 1406 PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void); 1407 PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void); 1408 1409 #else 1410 #define PetscSAWsBlock() PETSC_SUCCESS 1411 #define PetscObjectSAWsViewOff(obj) PETSC_SUCCESS 1412 #define PetscObjectSAWsSetBlock(obj, flg) PETSC_SUCCESS 1413 #define PetscObjectSAWsBlock(obj) PETSC_SUCCESS 1414 #define PetscObjectSAWsGrantAccess(obj) PETSC_SUCCESS 1415 #define PetscObjectSAWsTakeAccess(obj) PETSC_SUCCESS 1416 #define PetscStackViewSAWs() PETSC_SUCCESS 1417 #define PetscStackSAWsViewOff() PETSC_SUCCESS 1418 #define PetscStackSAWsTakeAccess() 1419 #define PetscStackSAWsGrantAccess() 1420 1421 #endif 1422 1423 PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[], PetscDLMode, PetscDLHandle *); 1424 PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *); 1425 PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle, const char[], void **); 1426 PETSC_EXTERN PetscErrorCode PetscDLAddr(void (*)(void), char **); 1427 #ifdef PETSC_HAVE_CXX 1428 PETSC_EXTERN PetscErrorCode PetscDemangleSymbol(const char *, char **); 1429 #endif 1430 1431 PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void *, PetscStack **); 1432 1433 PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE *, PetscBool); 1434 PETSC_EXTERN PetscErrorCode PetscObjectsView(PetscViewer); 1435 PETSC_EXTERN PetscErrorCode PetscObjectsGetObject(const char *, PetscObject *, const char **); 1436 PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList *); 1437 PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList, const char[], PetscObject *); 1438 PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList, PetscObject, char **, PetscBool *); 1439 PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *, const char[], PetscObject); 1440 PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *, const char[]); 1441 PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList, PetscObjectList *); 1442 1443 /* 1444 Dynamic library lists. Lists of names of routines in objects or in dynamic 1445 link libraries that will be loaded as needed. 1446 */ 1447 1448 #define PetscFunctionListAdd(list, name, fptr) PetscFunctionListAdd_Private((list), (name), (PetscVoidFn *)(fptr)) 1449 PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList *, const char[], PetscVoidFn *); 1450 PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList *); 1451 PETSC_EXTERN PetscErrorCode PetscFunctionListClear(PetscFunctionList); 1452 #define PetscFunctionListFind(list, name, fptr) PetscFunctionListFind_Private((list), (name), (PetscVoidFn **)(fptr)) 1453 PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList, const char[], PetscVoidFn **); 1454 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm, FILE *, const char[], const char[], const char[], const char[], PetscFunctionList, const char[], const char[]); 1455 PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList, PetscFunctionList *); 1456 PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList, PetscViewer); 1457 PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList, const char ***, int *); 1458 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintNonEmpty(PetscFunctionList); 1459 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintAll(void); 1460 1461 PETSC_EXTERN PetscDLLibrary PetscDLLibrariesLoaded; 1462 PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm, PetscDLLibrary *, const char[]); 1463 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm, PetscDLLibrary *, const char[]); 1464 PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm, PetscDLLibrary *, const char[], const char[], void **); 1465 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary); 1466 PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm, const char[], char *, size_t, PetscBool *); 1467 PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm, const char[], PetscDLLibrary *); 1468 PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary); 1469 1470 /* 1471 Useful utility routines 1472 */ 1473 PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm, PetscInt *, PetscInt *); 1474 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm, PetscInt, PetscInt *, PetscInt *); 1475 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm, PetscInt *, PetscInt *); 1476 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm, PetscMPIInt); 1477 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm, PetscMPIInt); 1478 PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject); 1479 PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE *); 1480 PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm, const PetscInt[2], PetscInt[2]); 1481 PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm, const PetscReal[2], PetscReal[2]); 1482 1483 /*MC 1484 PetscNot - negates a logical type value and returns result as a `PetscBool` 1485 1486 Level: beginner 1487 1488 Note: 1489 This is useful in cases like 1490 .vb 1491 int *a; 1492 PetscBool flag = PetscNot(a) 1493 .ve 1494 where !a would not return a `PetscBool` because we cannot provide a cast from int to `PetscBool` in C. 1495 1496 .seealso: `PetscBool`, `PETSC_TRUE`, `PETSC_FALSE` 1497 M*/ 1498 #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE) 1499 1500 /*MC 1501 PetscHelpPrintf - Prints help messages. 1502 1503 Synopsis: 1504 #include <petscsys.h> 1505 PetscErrorCode (*PetscHelpPrintf)(MPI_Comm comm, const char format[],args); 1506 1507 Not Collective, only applies on MPI rank 0; No Fortran Support 1508 1509 Input Parameters: 1510 + comm - the MPI communicator over which the help message is printed 1511 . format - the usual printf() format string 1512 - args - arguments to be printed 1513 1514 Level: developer 1515 1516 Notes: 1517 You can change how help messages are printed by replacing the function pointer with a function that does not simply write to stdout. 1518 1519 To use, write your own function, for example, 1520 .vb 1521 PetscErrorCode mypetschelpprintf(MPI_Comm comm,const char format[],....) 1522 { 1523 PetscFunctionReturn(PETSC_SUCCESS); 1524 } 1525 .ve 1526 then do the assignment 1527 .vb 1528 PetscHelpPrintf = mypetschelpprintf; 1529 .ve 1530 1531 You can do the assignment before `PetscInitialize()`. 1532 1533 The default routine used is called `PetscHelpPrintfDefault()`. 1534 1535 .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscErrorPrintf()`, `PetscHelpPrintfDefault()` 1536 M*/ 1537 PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3); 1538 1539 /* 1540 Defines PETSc profiling. 1541 */ 1542 #include <petsclog.h> 1543 1544 /* 1545 Simple PETSc parallel IO for ASCII printing 1546 */ 1547 PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[], char[]); 1548 PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm, const char[], const char[], FILE **); 1549 PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm, FILE *); 1550 PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4); 1551 PETSC_EXTERN PetscErrorCode PetscFFlush(FILE *); 1552 PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3); 1553 PETSC_EXTERN PetscErrorCode PetscSNPrintf(char *, size_t, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4); 1554 PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char *, size_t, const char[], size_t *, ...) PETSC_ATTRIBUTE_FORMAT(3, 5); 1555 PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[], size_t, const char *, PetscInt, const PetscReal[]); 1556 1557 PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2); 1558 PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2); 1559 PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3); 1560 1561 PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char *, size_t *); 1562 PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char *, char *); 1563 1564 PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm, const char[], const char[], const char[], FILE **); 1565 PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm, FILE *); 1566 PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]); 1567 1568 PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3); 1569 PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4); 1570 PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm, FILE *); 1571 PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm, FILE *, size_t, char[]); 1572 PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm, const char[], const char[], FILE **); 1573 PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char *[]); 1574 1575 PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID; 1576 PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer, void **); 1577 PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer, void *); 1578 PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer *); 1579 PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm, PetscContainer *); 1580 PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void *)); 1581 PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void *); 1582 PETSC_EXTERN PetscErrorCode PetscObjectContainerCompose(PetscObject, const char *name, void *, PetscErrorCode (*)(void *)); 1583 PETSC_EXTERN PetscErrorCode PetscObjectContainerQuery(PetscObject, const char *, void **); 1584 1585 /* 1586 For use in debuggers 1587 */ 1588 PETSC_EXTERN PetscMPIInt PetscGlobalRank; 1589 PETSC_EXTERN PetscMPIInt PetscGlobalSize; 1590 PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt, const PetscInt[], PetscViewer); 1591 PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt, const PetscReal[], PetscViewer); 1592 PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt, const PetscScalar[], PetscViewer); 1593 1594 /* 1595 Basic memory and string operations. These are usually simple wrappers 1596 around the basic Unix system calls, but a few of them have additional 1597 functionality and/or error checking. 1598 */ 1599 #include <petscstring.h> 1600 1601 #include <stddef.h> 1602 #include <stdlib.h> 1603 1604 #if defined(PETSC_CLANG_STATIC_ANALYZER) 1605 #define PetscPrefetchBlock(a, b, c, d) 1606 #else 1607 /*MC 1608 PetscPrefetchBlock - Prefetches a block of memory 1609 1610 Synopsis: 1611 #include <petscsys.h> 1612 void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t) 1613 1614 Not Collective 1615 1616 Input Parameters: 1617 + a - pointer to first element to fetch (any type but usually `PetscInt` or `PetscScalar`) 1618 . n - number of elements to fetch 1619 . rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors) 1620 - t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note 1621 1622 Level: developer 1623 1624 Notes: 1625 The last two arguments (`rw` and `t`) must be compile-time constants. 1626 1627 Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality. Not all architectures offer 1628 equivalent locality hints, but the following macros are always defined to their closest analogue. 1629 + `PETSC_PREFETCH_HINT_NTA` - Non-temporal. Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched). 1630 . `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. 1631 . `PETSC_PREFETCH_HINT_T1` - Fetch to level 2 and higher (not L1). 1632 - `PETSC_PREFETCH_HINT_T2` - Fetch to high-level cache only. (On many systems, T0 and T1 are equivalent.) 1633 1634 This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid 1635 address). 1636 1637 M*/ 1638 #define PetscPrefetchBlock(a, n, rw, t) \ 1639 do { \ 1640 const char *_p = (const char *)(a), *_end = (const char *)((a) + (n)); \ 1641 for (; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p, (rw), (t)); \ 1642 } while (0) 1643 #endif 1644 /* 1645 Determine if some of the kernel computation routines use 1646 Fortran (rather than C) for the numerical calculations. On some machines 1647 and compilers (like complex numbers) the Fortran version of the routines 1648 is faster than the C/C++ versions. The flag --with-fortran-kernels 1649 should be used with ./configure to turn these on. 1650 */ 1651 #if defined(PETSC_USE_FORTRAN_KERNELS) 1652 1653 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL) 1654 #define PETSC_USE_FORTRAN_KERNEL_MULTCRL 1655 #endif 1656 1657 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ) 1658 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ 1659 #endif 1660 1661 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ) 1662 #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ 1663 #endif 1664 1665 #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY) 1666 #define PETSC_USE_FORTRAN_KERNEL_MAXPY 1667 #endif 1668 1669 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 1670 #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ 1671 #endif 1672 1673 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ) 1674 #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ 1675 #endif 1676 1677 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ) 1678 #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ 1679 #endif 1680 1681 #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT) 1682 #define PETSC_USE_FORTRAN_KERNEL_MDOT 1683 #endif 1684 1685 #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY) 1686 #define PETSC_USE_FORTRAN_KERNEL_XTIMESY 1687 #endif 1688 1689 #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX) 1690 #define PETSC_USE_FORTRAN_KERNEL_AYPX 1691 #endif 1692 1693 #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY) 1694 #define PETSC_USE_FORTRAN_KERNEL_WAXPY 1695 #endif 1696 1697 #endif 1698 1699 /* 1700 Macros for indicating code that should be compiled with a C interface, 1701 rather than a C++ interface. Any routines that are dynamically loaded 1702 (such as the PCCreate_XXX() routines) must be wrapped so that the name 1703 mangler does not change the functions symbol name. This just hides the 1704 ugly extern "C" {} wrappers. 1705 */ 1706 #if defined(__cplusplus) 1707 #define EXTERN_C_BEGIN extern "C" { 1708 #define EXTERN_C_END } 1709 #else 1710 #define EXTERN_C_BEGIN 1711 #define EXTERN_C_END 1712 #endif 1713 1714 /*MC 1715 MPI_Comm - the basic object used by MPI to determine which processes are involved in a 1716 communication 1717 1718 Level: beginner 1719 1720 Note: 1721 This manual page is a place-holder because MPICH does not have a manual page for `MPI_Comm` 1722 1723 .seealso: `PETSC_COMM_WORLD`, `PETSC_COMM_SELF` 1724 M*/ 1725 1726 #if defined(PETSC_HAVE_MPIIO) 1727 PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4); 1728 PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4); 1729 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); 1730 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); 1731 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); 1732 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); 1733 #endif 1734 1735 #if defined(PETSC_HAVE_MPI_COUNT) 1736 typedef MPI_Count MPIU_Count; 1737 #else 1738 typedef PetscInt64 MPIU_Count; 1739 #endif 1740 1741 /*@C 1742 PetscIntCast - casts a `MPI_Count`, `PetscInt64`, `PetscCount`, or `size_t` to a `PetscInt` (which may be 32-bits in size), generates an 1743 error if the `PetscInt` is not large enough to hold the number. 1744 1745 Not Collective; No Fortran Support 1746 1747 Input Parameter: 1748 . a - the `PetscInt64` value 1749 1750 Output Parameter: 1751 . b - the resulting `PetscInt` value, or `NULL` if the result is not needed 1752 1753 Level: advanced 1754 1755 Note: 1756 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 1757 1758 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscIntMultError()`, `PetscIntSumError()` 1759 @*/ 1760 static inline PetscErrorCode PetscIntCast(MPIU_Count a, PetscInt *b) 1761 { 1762 PetscFunctionBegin; 1763 if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */ 1764 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); 1765 if (b) *b = (PetscInt)a; 1766 PetscFunctionReturn(PETSC_SUCCESS); 1767 } 1768 1769 /*@C 1770 PetscBLASIntCast - casts a `MPI_Count`, `PetscInt`, `PetscCount` or `PetscInt64` to a `PetscBLASInt` (which may be 32-bits in size), generates an 1771 error if the `PetscBLASInt` is not large enough to hold the number. 1772 1773 Not Collective; No Fortran Support 1774 1775 Input Parameter: 1776 . a - the `PetscInt` value 1777 1778 Output Parameter: 1779 . b - the resulting `PetscBLASInt` value, or `NULL` if the result is not needed 1780 1781 Level: advanced 1782 1783 Note: 1784 Errors if the integer is negative since PETSc calls to BLAS/LAPACK never need to cast negative integer inputs 1785 1786 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscIntCast()` 1787 @*/ 1788 static inline PetscErrorCode PetscBLASIntCast(MPIU_Count a, PetscBLASInt *b) 1789 { 1790 PetscFunctionBegin; 1791 if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */ 1792 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); 1793 PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer to BLAS/LAPACK routine"); 1794 if (b) *b = (PetscBLASInt)a; 1795 PetscFunctionReturn(PETSC_SUCCESS); 1796 } 1797 1798 /*@C 1799 PetscCuBLASIntCast - like `PetscBLASIntCast()`, but for `PetscCuBLASInt`. 1800 1801 Not Collective; No Fortran Support 1802 1803 Input Parameter: 1804 . a - the `PetscInt` value 1805 1806 Output Parameter: 1807 . b - the resulting `PetscCuBLASInt` value, or `NULL` if the result is not needed 1808 1809 Level: advanced 1810 1811 Note: 1812 Errors if the integer is negative since PETSc calls to cuBLAS and friends never need to cast negative integer inputs 1813 1814 .seealso: `PetscCuBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscIntCast()` 1815 @*/ 1816 static inline PetscErrorCode PetscCuBLASIntCast(MPIU_Count a, PetscCuBLASInt *b) 1817 { 1818 PetscFunctionBegin; 1819 if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */ 1820 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); 1821 PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer %" PetscInt64_FMT "to cuBLAS routine", (PetscInt64)a); 1822 if (b) *b = (PetscCuBLASInt)a; 1823 PetscFunctionReturn(PETSC_SUCCESS); 1824 } 1825 1826 /*@C 1827 PetscHipBLASIntCast - like `PetscBLASIntCast()`, but for `PetscHipBLASInt`. 1828 1829 Not Collective; No Fortran Support 1830 1831 Input Parameter: 1832 . a - the `PetscInt` value 1833 1834 Output Parameter: 1835 . b - the resulting `PetscHipBLASInt` value, or `NULL` if the result is not needed 1836 1837 Level: advanced 1838 1839 Note: 1840 Errors if the integer is negative since PETSc calls to hipBLAS and friends never need to cast negative integer inputs 1841 1842 .seealso: `PetscHipBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscIntCast()` 1843 @*/ 1844 static inline PetscErrorCode PetscHipBLASIntCast(MPIU_Count a, PetscHipBLASInt *b) 1845 { 1846 PetscFunctionBegin; 1847 if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */ 1848 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); 1849 PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer %" PetscInt64_FMT "to hipBLAS routine", (PetscInt64)a); 1850 if (b) *b = (PetscHipBLASInt)a; 1851 PetscFunctionReturn(PETSC_SUCCESS); 1852 } 1853 1854 /*@C 1855 PetscMPIIntCast - casts a `MPI_Count`, `PetscInt`, `PetscCount`, or `PetscInt64` to a `PetscMPIInt` (which is always 32-bits in size), generates an 1856 error if the `PetscMPIInt` is not large enough to hold the number. 1857 1858 Not Collective; No Fortran Support 1859 1860 Input Parameter: 1861 . a - the `PetscInt` value 1862 1863 Output Parameter: 1864 . b - the resulting `PetscMPIInt` value, or `NULL` if the result is not needed 1865 1866 Level: advanced 1867 1868 .seealso: [](stylePetscCount), `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscIntCast()` 1869 @*/ 1870 static inline PetscErrorCode PetscMPIIntCast(MPIU_Count a, PetscMPIInt *b) 1871 { 1872 PetscFunctionBegin; 1873 if (b) *b = 0; /* to prevent compilers erroneously suggesting uninitialized variable */ 1874 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); 1875 if (b) *b = (PetscMPIInt)a; 1876 PetscFunctionReturn(PETSC_SUCCESS); 1877 } 1878 1879 #define PetscInt64Mult(a, b) (((PetscInt64)(a)) * ((PetscInt64)(b))) 1880 1881 /*@C 1882 PetscRealIntMultTruncate - Computes the product of a positive `PetscReal` and a positive 1883 `PetscInt` and truncates the value to slightly less than the maximal possible value. 1884 1885 Not Collective; No Fortran Support 1886 1887 Input Parameters: 1888 + a - The `PetscReal` value 1889 - b - The `PetscInt` value 1890 1891 Level: advanced 1892 1893 Notes: 1894 Returns the result as a `PetscInt` value. 1895 1896 Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64`. 1897 1898 Use `PetscIntMultTruncate()` to compute the product of two positive `PetscInt` and truncate 1899 to fit a `PetscInt`. 1900 1901 Use `PetscIntMultError()` to compute the product of two `PetscInt` if you wish to generate an 1902 error if the result will not fit in a `PetscInt`. 1903 1904 Developer Notes: 1905 We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but 1906 requires many more checks. 1907 1908 This is used where we compute approximate sizes for workspace and need to insure the 1909 workspace is index-able. 1910 1911 .seealso: `PetscReal`, `PetscInt`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()` 1912 @*/ 1913 static inline PetscInt PetscRealIntMultTruncate(PetscReal a, PetscInt b) 1914 { 1915 PetscInt64 r = (PetscInt64)(a * (PetscReal)b); 1916 if (r > PETSC_INT_MAX - 100) r = PETSC_INT_MAX - 100; 1917 return (PetscInt)r; 1918 } 1919 1920 /*@C 1921 PetscIntMultTruncate - Computes the product of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value 1922 1923 Not Collective; No Fortran Support 1924 1925 Input Parameters: 1926 + a - the `PetscInt` value 1927 - b - the second value 1928 1929 Returns: 1930 The result as a `PetscInt` value 1931 1932 Level: advanced 1933 1934 Notes: 1935 Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64` 1936 1937 Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt` 1938 1939 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` 1940 1941 Developer Notes: 1942 We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but requires many more checks. 1943 1944 This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able. 1945 1946 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()`, 1947 `PetscIntSumTruncate()` 1948 @*/ 1949 static inline PetscInt PetscIntMultTruncate(PetscInt a, PetscInt b) 1950 { 1951 PetscInt64 r = PetscInt64Mult(a, b); 1952 if (r > PETSC_INT_MAX - 100) r = PETSC_INT_MAX - 100; 1953 return (PetscInt)r; 1954 } 1955 1956 /*@C 1957 PetscIntSumTruncate - Computes the sum of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value 1958 1959 Not Collective; No Fortran Support 1960 1961 Input Parameters: 1962 + a - the `PetscInt` value 1963 - b - the second value 1964 1965 Returns: 1966 The result as a `PetscInt` value 1967 1968 Level: advanced 1969 1970 Notes: 1971 Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64` 1972 1973 Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt` 1974 1975 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` 1976 1977 Developer Note: 1978 This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able. 1979 1980 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()` 1981 @*/ 1982 static inline PetscInt PetscIntSumTruncate(PetscInt a, PetscInt b) 1983 { 1984 PetscInt64 r = ((PetscInt64)a) + ((PetscInt64)b); 1985 if (r > PETSC_INT_MAX - 100) r = PETSC_INT_MAX - 100; 1986 return (PetscInt)r; 1987 } 1988 1989 /*@C 1990 PetscIntMultError - Computes the product of two positive `PetscInt` and generates an error with overflow. 1991 1992 Not Collective; No Fortran Support 1993 1994 Input Parameters: 1995 + a - the `PetscInt` value 1996 - b - the second value 1997 1998 Output Parameter: 1999 . result - the result as a `PetscInt` value, or `NULL` if you do not want the result, you just want to check if it overflows 2000 2001 Level: advanced 2002 2003 Notes: 2004 Use `PetscInt64Mult()` to compute the product of two `PetscInt` and store in a `PetscInt64` 2005 2006 Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt` 2007 2008 Developer Note: 2009 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. 2010 `PetscIntSumError()` can be used to check for this situation. 2011 2012 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscIntMult64()`, `PetscIntSumError()` 2013 @*/ 2014 static inline PetscErrorCode PetscIntMultError(PetscInt a, PetscInt b, PetscInt *result) 2015 { 2016 PetscInt64 r = PetscInt64Mult(a, b); 2017 2018 PetscFunctionBegin; 2019 if (result) *result = (PetscInt)r; 2020 if (!PetscDefined(USE_64BIT_INDICES)) { 2021 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); 2022 } 2023 PetscFunctionReturn(PETSC_SUCCESS); 2024 } 2025 2026 /*@C 2027 2028 PetscIntSumError - Computes the sum of two positive `PetscInt` and generates an error with overflow. 2029 2030 Not Collective; No Fortran Support 2031 2032 Input Parameters: 2033 + a - the `PetscInt` value 2034 - b - the second value 2035 2036 Output Parameter: 2037 . c - the result as a `PetscInt` value, or `NULL` if you do not want the result, you just want to check if it overflows 2038 2039 Level: advanced 2040 2041 Notes: 2042 Use `PetscInt64Mult()` to compute the product of two 32-bit `PetscInt` and store in a `PetscInt64` 2043 2044 Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt` 2045 2046 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()` 2047 @*/ 2048 static inline PetscErrorCode PetscIntSumError(PetscInt a, PetscInt b, PetscInt *result) 2049 { 2050 PetscInt64 r = ((PetscInt64)a) + ((PetscInt64)b); 2051 2052 PetscFunctionBegin; 2053 if (result) *result = (PetscInt)r; 2054 if (!PetscDefined(USE_64BIT_INDICES)) { 2055 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); 2056 } 2057 PetscFunctionReturn(PETSC_SUCCESS); 2058 } 2059 2060 /* 2061 The IBM include files define hz, here we hide it so that it may be used as a regular user variable. 2062 */ 2063 #if defined(hz) 2064 #undef hz 2065 #endif 2066 2067 #if defined(PETSC_HAVE_SYS_TYPES_H) 2068 #include <sys/types.h> 2069 #endif 2070 2071 /*MC 2072 2073 PETSC_VERSION - This manual page provides information about how PETSc documents and uses its version information. This information is available to both C/C++ 2074 and Fortran compilers when `petscsys.h` is included. 2075 2076 The current PETSc version and the API for accessing it are defined in <A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/include/petscversion.h.html">include/petscversion.html</A> 2077 2078 The complete version number is given as the triple PETSC_VERSION_MAJOR.PETSC_VERSION_MINOR.PETSC_VERSION_SUBMINOR (in short hand x.y.z) 2079 2080 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 2081 where only a change in the major version number (x) indicates a change in the API. 2082 2083 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 2084 2085 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), 2086 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 2087 version number (x.y.z). 2088 2089 `PETSC_RELEASE_DATE` is the date the x.y version was released (i.e. the version before any patch releases) 2090 2091 `PETSC_VERSION_DATE` is the date the x.y.z version was released 2092 2093 `PETSC_VERSION_GIT` is the last git commit to the repository given in the form vx.y.z-wwwww 2094 2095 `PETSC_VERSION_DATE_GIT` is the date of the last git commit to the repository 2096 2097 `PETSC_VERSION_()` is deprecated and will eventually be removed. 2098 2099 Level: intermediate 2100 M*/ 2101 2102 PETSC_EXTERN PetscErrorCode PetscGetArchType(char[], size_t); 2103 PETSC_EXTERN PetscErrorCode PetscGetHostName(char[], size_t); 2104 PETSC_EXTERN PetscErrorCode PetscGetUserName(char[], size_t); 2105 PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[], size_t); 2106 PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]); 2107 PETSC_EXTERN PetscErrorCode PetscGetDate(char[], size_t); 2108 PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t); 2109 PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt *, PetscInt *, PetscInt *, PetscInt *); 2110 2111 PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscCount, const PetscInt[], PetscBool *); 2112 PETSC_EXTERN PetscErrorCode PetscSortedInt64(PetscCount, const PetscInt64[], PetscBool *); 2113 PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscCount, const PetscMPIInt[], PetscBool *); 2114 PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscCount, const PetscReal[], PetscBool *); 2115 PETSC_EXTERN PetscErrorCode PetscSortInt(PetscCount, PetscInt[]); 2116 PETSC_EXTERN PetscErrorCode PetscSortInt64(PetscCount, PetscInt64[]); 2117 PETSC_EXTERN PetscErrorCode PetscSortCount(PetscCount, PetscCount[]); 2118 PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscCount, PetscInt[]); 2119 PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *, PetscInt[]); 2120 PETSC_EXTERN PetscErrorCode PetscSortedCheckDupsInt(PetscCount, const PetscInt[], PetscBool *); 2121 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt *, PetscInt[]); 2122 PETSC_EXTERN PetscErrorCode PetscCheckDupsInt(PetscInt, const PetscInt[], PetscBool *); 2123 PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscCount, const PetscInt[], PetscInt *); 2124 PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscCount, const PetscMPIInt[], PetscInt *); 2125 PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt, const PetscInt[], PetscInt[]); 2126 PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt, const char *[], PetscInt[]); 2127 PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscCount, PetscInt[], PetscInt[]); 2128 PETSC_EXTERN PetscErrorCode PetscSortIntWithCountArray(PetscCount, PetscInt[], PetscCount[]); 2129 PETSC_EXTERN PetscErrorCode PetscSortIntWithMPIIntArray(PetscCount, PetscInt[], PetscMPIInt[]); 2130 PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscCount, PetscInt[], PetscInt[], PetscInt[]); 2131 PETSC_EXTERN PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount, PetscInt[], PetscInt[], PetscCount[]); 2132 PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscCount, PetscMPIInt[]); 2133 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *, PetscMPIInt[]); 2134 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscCount, PetscMPIInt[], PetscMPIInt[]); 2135 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscCount, PetscMPIInt[], PetscInt[]); 2136 PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscCount, PetscInt[], PetscScalar[]); 2137 PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscCount, PetscInt[], void *, size_t, void *); 2138 PETSC_EXTERN PetscErrorCode PetscSortReal(PetscCount, PetscReal[]); 2139 PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscCount, PetscReal[], PetscInt[]); 2140 PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt, const PetscReal[], PetscInt[]); 2141 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt *, PetscReal[]); 2142 PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal, PetscCount, const PetscReal[], PetscReal, PetscInt *); 2143 PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt, PetscInt, PetscScalar[], PetscInt[]); 2144 PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt, PetscInt, PetscReal[], PetscInt[]); 2145 PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt, const PetscBool[], const PetscInt[], PetscInt *, PetscInt **, PetscInt **, PetscInt **, PetscInt **); 2146 PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt, const PetscInt[], const PetscInt[], PetscInt, const PetscInt[], const PetscInt[], PetscInt *, PetscInt **, PetscInt **); 2147 PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscInt *, PetscInt **); 2148 PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt, const PetscMPIInt[], PetscInt, const PetscMPIInt[], PetscInt *, PetscMPIInt **); 2149 PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *); 2150 2151 PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt, void *, size_t, int (*)(const void *, const void *, void *), void *); 2152 PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt, PetscInt[]); 2153 PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt, PetscMPIInt[]); 2154 PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt, PetscReal[]); 2155 PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt, void *, size_t, void *, size_t, int (*)(const void *, const void *, void *), void *); 2156 PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt, PetscInt[], PetscInt[]); 2157 PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt, PetscMPIInt[], PetscMPIInt[]); 2158 PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt, PetscReal[], PetscInt[]); 2159 2160 PETSC_EXTERN PetscErrorCode PetscSetDisplay(void); 2161 PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[], size_t); 2162 2163 /*J 2164 PetscRandomType - String with the name of a PETSc randomizer 2165 2166 Level: beginner 2167 2168 Note: 2169 To use `PETSCSPRNG` or `PETSCRANDOM123` you must have ./configure PETSc 2170 with the option `--download-sprng` or `--download-random123`. We recommend the default provided with PETSc. 2171 2172 .seealso: `PetscRandomSetType()`, `PetscRandom`, `PetscRandomCreate()` 2173 J*/ 2174 typedef const char *PetscRandomType; 2175 #define PETSCRAND "rand" 2176 #define PETSCRAND48 "rand48" 2177 #define PETSCSPRNG "sprng" 2178 #define PETSCRANDER48 "rander48" 2179 #define PETSCRANDOM123 "random123" 2180 #define PETSCCURAND "curand" 2181 2182 /* Logging support */ 2183 PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID; 2184 2185 PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void); 2186 PETSC_EXTERN PetscErrorCode PetscRandomFinalizePackage(void); 2187 2188 /* Dynamic creation and loading functions */ 2189 PETSC_EXTERN PetscFunctionList PetscRandomList; 2190 2191 PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[], PetscErrorCode (*)(PetscRandom)); 2192 PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType); 2193 PETSC_EXTERN PetscErrorCode PetscRandomSetOptionsPrefix(PetscRandom, const char[]); 2194 PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom); 2195 PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType *); 2196 PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom, PetscObject, const char[]); 2197 PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom, PetscViewer); 2198 2199 PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm, PetscRandom *); 2200 PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom, PetscScalar *); 2201 PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom, PetscReal *); 2202 PETSC_EXTERN PetscErrorCode PetscRandomGetValues(PetscRandom, PetscInt, PetscScalar *); 2203 PETSC_EXTERN PetscErrorCode PetscRandomGetValuesReal(PetscRandom, PetscInt, PetscReal *); 2204 PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom, PetscScalar *, PetscScalar *); 2205 PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom, PetscScalar, PetscScalar); 2206 PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom, unsigned long); 2207 PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom, unsigned long *); 2208 PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom); 2209 PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom *); 2210 2211 PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[], char[], size_t); 2212 PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[], char[], size_t); 2213 PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[], size_t); 2214 PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[], char[]); 2215 PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[], size_t); 2216 PETSC_EXTERN PetscErrorCode PetscTestFile(const char[], char, PetscBool *); 2217 PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[], char, PetscBool *); 2218 PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]); 2219 PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]); 2220 PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]); 2221 2222 /*MC 2223 PetscBinaryBigEndian - indicates if values in memory are stored with big endian format 2224 2225 Synopsis: 2226 #include <petscsys.h> 2227 PetscBool PetscBinaryBigEndian(void); 2228 2229 No Fortran Support 2230 2231 Level: developer 2232 2233 .seealso: `PetscInitialize()`, `PetscFinalize()`, `PetscInitializeCalled` 2234 M*/ 2235 static inline PetscBool PetscBinaryBigEndian(void) 2236 { 2237 long _petsc_v = 1; 2238 return ((char *)&_petsc_v)[0] ? PETSC_FALSE : PETSC_TRUE; 2239 } 2240 2241 PETSC_EXTERN PetscErrorCode PetscBinaryRead(int, void *, PetscCount, PetscInt *, PetscDataType); 2242 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm, int, void *, PetscInt, PetscInt *, PetscDataType); 2243 PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int, const void *, PetscCount, PetscDataType); 2244 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm, int, const void *, PetscInt, PetscDataType); 2245 PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[], PetscFileMode, int *); 2246 PETSC_EXTERN PetscErrorCode PetscBinaryClose(int); 2247 PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm, PetscBool *); 2248 PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm, PetscBool *); 2249 PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm, char[], size_t); 2250 PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm, const char[], char[], size_t, PetscBool *); 2251 PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm, const char[], char[], size_t, PetscBool *); 2252 #if defined(PETSC_USE_SOCKET_VIEWER) 2253 PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[], int, int *); 2254 #endif 2255 2256 PETSC_EXTERN PetscErrorCode PetscBinarySeek(int, off_t, PetscBinarySeekType, off_t *); 2257 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm, int, off_t, PetscBinarySeekType, off_t *); 2258 PETSC_EXTERN PetscErrorCode PetscByteSwap(void *, PetscDataType, PetscCount); 2259 2260 PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]); 2261 PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[], PetscBool); 2262 PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void); 2263 PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char *); 2264 PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void); 2265 PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void); 2266 PETSC_EXTERN PetscErrorCode PetscWaitOnError(void); 2267 2268 PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt *); 2269 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], PetscMPIInt **, PetscMPIInt **); 2270 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt **, PetscMPIInt **, PetscMPIInt **); 2271 PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscInt ***, MPI_Request **); 2272 PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscScalar ***, MPI_Request **); 2273 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); 2274 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm, PetscMPIInt, MPI_Datatype, PetscMPIInt, const PetscMPIInt[], const void *, PetscMPIInt *, PetscMPIInt **, void *, PetscMPIInt, PetscErrorCode (*send)(MPI_Comm, const PetscMPIInt[], PetscMPIInt, PetscMPIInt, void *, MPI_Request[], void *), PetscErrorCode (*recv)(MPI_Comm, const PetscMPIInt[], PetscMPIInt, void *, MPI_Request[], void *), void *ctx) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(6, 3); 2275 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm, PetscMPIInt, MPI_Datatype, PetscMPIInt, const PetscMPIInt[], const void *, PetscMPIInt *, PetscMPIInt **, void *, PetscMPIInt, MPI_Request **, MPI_Request **, PetscErrorCode (*send)(MPI_Comm, const PetscMPIInt[], PetscMPIInt, PetscMPIInt, void *, MPI_Request[], void *), PetscErrorCode (*recv)(MPI_Comm, const PetscMPIInt[], PetscMPIInt, void *, MPI_Request[], void *), void *ctx) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(6, 3); 2276 2277 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm, PetscBuildTwoSidedType); 2278 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm, PetscBuildTwoSidedType *); 2279 2280 PETSC_EXTERN PetscErrorCode PetscSSEIsEnabled(MPI_Comm, PetscBool *, PetscBool *); 2281 2282 PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject); 2283 2284 struct _n_PetscSubcomm { 2285 MPI_Comm parent; /* parent communicator */ 2286 MPI_Comm dupparent; /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */ 2287 MPI_Comm child; /* the sub-communicator */ 2288 PetscMPIInt n; /* num of subcommunicators under the parent communicator */ 2289 PetscMPIInt color; /* color of processors belong to this communicator */ 2290 PetscMPIInt *subsize; /* size of subcommunicator[color] */ 2291 PetscSubcommType type; 2292 char *subcommprefix; 2293 }; 2294 2295 static inline MPI_Comm PetscSubcommParent(PetscSubcomm scomm) 2296 { 2297 return scomm->parent; 2298 } 2299 static inline MPI_Comm PetscSubcommChild(PetscSubcomm scomm) 2300 { 2301 return scomm->child; 2302 } 2303 static inline MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm) 2304 { 2305 return scomm->dupparent; 2306 } 2307 PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm, PetscSubcomm *); 2308 PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm *); 2309 PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm, PetscInt); 2310 PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm, PetscSubcommType); 2311 PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm, PetscMPIInt, PetscMPIInt); 2312 PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm, PetscViewer); 2313 PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm); 2314 PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm, const char[]); 2315 PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm, MPI_Comm *); 2316 PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm, MPI_Comm *); 2317 PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm, MPI_Comm *); 2318 2319 PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt, PetscHeap *); 2320 PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap, PetscInt, PetscInt); 2321 PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap, PetscInt *, PetscInt *); 2322 PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap, PetscInt *, PetscInt *); 2323 PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap, PetscInt, PetscInt); 2324 PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap); 2325 PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap *); 2326 PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap, PetscViewer); 2327 2328 PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer); 2329 PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm, PetscShmComm *); 2330 PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm, PetscMPIInt, PetscMPIInt *); 2331 PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm, PetscMPIInt, PetscMPIInt *); 2332 PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm, MPI_Comm *); 2333 2334 /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */ 2335 PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm, PetscInt, PetscOmpCtrl *); 2336 PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl, MPI_Comm *, MPI_Comm *, PetscBool *); 2337 PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl *); 2338 PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl); 2339 PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl); 2340 PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl); 2341 2342 PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t, PetscCount, PetscSegBuffer *); 2343 PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer *); 2344 PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer, PetscCount, void *); 2345 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer, void *); 2346 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer, void *); 2347 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer, void *); 2348 PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer, PetscCount *); 2349 PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer, PetscCount); 2350 2351 /*MC 2352 PetscSegBufferGetInts - access an array of `PetscInt` from a `PetscSegBuffer` 2353 2354 Synopsis: 2355 #include <petscsys.h> 2356 PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg, size_t count, PetscInt *PETSC_RESTRICT *slot); 2357 2358 No Fortran Support 2359 2360 Input Parameters: 2361 + seg - `PetscSegBuffer` buffer 2362 - count - number of entries needed 2363 2364 Output Parameter: 2365 . buf - address of new buffer for contiguous data 2366 2367 Level: intermediate 2368 2369 Developer Note: 2370 Type-safe wrapper to encourage use of PETSC_RESTRICT. Does not use PetscFunctionBegin because the error handling 2371 prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as 2372 possible. 2373 2374 .seealso: `PetscSegBuffer`, `PetscSegBufferGet()`, `PetscInitialize()`, `PetscFinalize()`, `PetscInitializeCalled` 2375 M*/ 2376 /* */ 2377 static inline PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg, PetscCount count, PetscInt *PETSC_RESTRICT *slot) 2378 { 2379 return PetscSegBufferGet(seg, count, (void **)slot); 2380 } 2381 2382 extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton; 2383 PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted *); 2384 PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted *); 2385 PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted, const char *, const char *, PetscBool *); 2386 2387 #include <stdarg.h> 2388 PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char *, size_t, const char[], size_t *, va_list); 2389 PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE *, const char[], va_list); 2390 2391 PETSC_EXTERN PetscSegBuffer PetscCitationsList; 2392 2393 /*@ 2394 PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code. 2395 2396 Not Collective; No Fortran Support 2397 2398 Input Parameters: 2399 + cite - the bibtex item, formatted to displayed on multiple lines nicely 2400 - set - a boolean variable initially set to `PETSC_FALSE`; this is used to insure only a single registration of the citation 2401 2402 Options Database Key: 2403 . -citations [filename] - print out the bibtex entries for the given computation 2404 2405 Level: intermediate 2406 @*/ 2407 static inline PetscErrorCode PetscCitationsRegister(const char cit[], PetscBool *set) 2408 { 2409 size_t len; 2410 char *vstring; 2411 2412 PetscFunctionBegin; 2413 if (set && *set) PetscFunctionReturn(PETSC_SUCCESS); 2414 PetscCall(PetscStrlen(cit, &len)); 2415 PetscCall(PetscSegBufferGet(PetscCitationsList, (PetscCount)len, &vstring)); 2416 PetscCall(PetscArraycpy(vstring, cit, len)); 2417 if (set) *set = PETSC_TRUE; 2418 PetscFunctionReturn(PETSC_SUCCESS); 2419 } 2420 2421 PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm, char[], char[], size_t); 2422 PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm, const char[], char[], size_t); 2423 PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm, const char[], const char[]); 2424 2425 PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm, char[], char[], size_t); 2426 PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm, const char[], char[], char[], size_t); 2427 PETSC_EXTERN PetscErrorCode PetscBoxUpload(MPI_Comm, const char[], const char[]); 2428 2429 PETSC_EXTERN PetscErrorCode PetscGlobusGetTransfers(MPI_Comm, const char[], char[], size_t); 2430 PETSC_EXTERN PetscErrorCode PetscGlobusAuthorize(MPI_Comm, char[], size_t); 2431 PETSC_EXTERN PetscErrorCode PetscGlobusUpload(MPI_Comm, const char[], const char[]); 2432 2433 PETSC_EXTERN PetscErrorCode PetscPullJSONValue(const char[], const char[], char[], size_t, PetscBool *); 2434 PETSC_EXTERN PetscErrorCode PetscPushJSONValue(char[], const char[], const char[], size_t); 2435 2436 #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) 2437 PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint, PetscMPIInt, MPI_Info, MPI_Comm, void *, MPI_Win *); 2438 PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win, PetscMPIInt, MPI_Aint *, PetscMPIInt *, void *); 2439 #endif 2440 2441 #if !defined(PETSC_HAVE_MPI_LARGE_COUNT) 2442 /* 2443 Cast PetscCount <a> to PetscMPIInt <b>, where <a> is likely used for the 'count' argument in MPI routines. 2444 It is similar to PetscMPIIntCast() except that here it returns an MPI error code. 2445 */ 2446 #define PetscMPIIntCast_Internal(a, b) \ 2447 do { \ 2448 *b = 0; \ 2449 if (PetscUnlikely(a > (MPIU_Count)PETSC_MPI_INT_MAX)) return MPI_ERR_COUNT; \ 2450 *b = (PetscMPIInt)a; \ 2451 } while (0) 2452 2453 static inline PetscMPIInt MPIU_Get_count(MPI_Status *status, MPI_Datatype dtype, PetscCount *count) 2454 { 2455 PetscMPIInt count2, err; 2456 2457 *count = 0; /* to prevent incorrect warnings of uninitialized variables */ 2458 err = MPI_Get_count(status, dtype, &count2); 2459 *count = count2; 2460 return err; 2461 } 2462 2463 static inline PetscMPIInt MPIU_Send(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt dest, PetscMPIInt tag, MPI_Comm comm) 2464 { 2465 PetscMPIInt count2, err; 2466 2467 PetscMPIIntCast_Internal(count, &count2); 2468 err = MPI_Send((void *)buf, count2, dtype, dest, tag, comm); 2469 return err; 2470 } 2471 2472 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) 2473 { 2474 PetscMPIInt count2, err; 2475 2476 PetscMPIIntCast_Internal(count, &count2); 2477 err = MPI_Send_init((void *)buf, count2, dtype, dest, tag, comm, request); 2478 return err; 2479 } 2480 2481 static inline PetscMPIInt MPIU_Isend(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt dest, PetscMPIInt tag, MPI_Comm comm, MPI_Request *request) 2482 { 2483 PetscMPIInt count2, err; 2484 2485 PetscMPIIntCast_Internal(count, &count2); 2486 err = MPI_Isend((void *)buf, count2, dtype, dest, tag, comm, request); 2487 return err; 2488 } 2489 2490 static inline PetscMPIInt MPIU_Recv(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt source, PetscMPIInt tag, MPI_Comm comm, MPI_Status *status) 2491 { 2492 PetscMPIInt count2, err; 2493 2494 PetscMPIIntCast_Internal(count, &count2); 2495 err = MPI_Recv((void *)buf, count2, dtype, source, tag, comm, status); 2496 return err; 2497 } 2498 2499 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) 2500 { 2501 PetscMPIInt count2, err; 2502 2503 PetscMPIIntCast_Internal(count, &count2); 2504 err = MPI_Recv_init((void *)buf, count2, dtype, source, tag, comm, request); 2505 return err; 2506 } 2507 2508 static inline PetscMPIInt MPIU_Irecv(const void *buf, MPIU_Count count, MPI_Datatype dtype, PetscMPIInt source, PetscMPIInt tag, MPI_Comm comm, MPI_Request *request) 2509 { 2510 PetscMPIInt count2, err; 2511 2512 PetscMPIIntCast_Internal(count, &count2); 2513 err = MPI_Irecv((void *)buf, count2, dtype, source, tag, comm, request); 2514 return err; 2515 } 2516 2517 static inline PetscMPIInt MPIU_Reduce(const void *inbuf, void *outbuf, MPIU_Count count, MPI_Datatype dtype, MPI_Op op, PetscMPIInt root, MPI_Comm comm) 2518 { 2519 PetscMPIInt count2, err; 2520 2521 PetscMPIIntCast_Internal(count, &count2); 2522 err = MPI_Reduce((void *)inbuf, outbuf, count2, dtype, op, root, comm); 2523 return err; 2524 } 2525 2526 #if defined(PETSC_HAVE_MPI_REDUCE_LOCAL) 2527 static inline PetscMPIInt MPIU_Reduce_local(const void *inbuf, void *inoutbuf, MPIU_Count count, MPI_Datatype dtype, MPI_Op op) 2528 { 2529 PetscMPIInt count2, err; 2530 2531 PetscMPIIntCast_Internal(count, &count2); 2532 err = MPI_Reduce_local((void *)inbuf, inoutbuf, count2, dtype, op); 2533 return err; 2534 } 2535 #endif 2536 2537 #else 2538 2539 /* on 32 bit systems MPI_Count maybe 64-bit while PetscCount is 32-bit */ 2540 #define PetscCountCast_Internal(a, b) \ 2541 do { \ 2542 *b = 0; \ 2543 if (PetscUnlikely(a > (MPI_Count)PETSC_COUNT_MAX)) return MPI_ERR_COUNT; \ 2544 *b = (PetscMPIInt)a; \ 2545 } while (0) 2546 2547 static inline PetscMPIInt MPIU_Get_count(MPI_Status *status, MPI_Datatype dtype, PetscCount *count) 2548 { 2549 MPI_Count count2; 2550 PetscMPIInt err; 2551 2552 *count = 0; /* to prevent incorrect warnings of uninitialized variables */ 2553 err = MPI_Get_count_c(status, dtype, &count2); 2554 if (err) return err; 2555 PetscCountCast_Internal(count2, count); 2556 return MPI_SUCCESS; 2557 } 2558 2559 #define MPIU_Reduce(inbuf, outbuf, count, dtype, op, root, comm) MPI_Reduce_c(inbuf, outbuf, (MPI_Count)(count), dtype, op, root, comm) 2560 #define MPIU_Send(buf, count, dtype, dest, tag, comm) MPI_Send_c(buf, (MPI_Count)(count), dtype, dest, tag, comm) 2561 #define MPIU_Send_init(buf, count, dtype, dest, tag, comm, request) MPI_Send_init_c(buf, (MPI_Count)(count), dtype, dest, tag, comm, request) 2562 #define MPIU_Isend(buf, count, dtype, dest, tag, comm, request) MPI_Isend_c(buf, (MPI_Count)(count), dtype, dest, tag, comm, request) 2563 #define MPIU_Recv(buf, count, dtype, source, tag, comm, status) MPI_Recv_c(buf, (MPI_Count)(count), dtype, source, tag, comm, status) 2564 #define MPIU_Recv_init(buf, count, dtype, source, tag, comm, request) MPI_Recv_init_c(buf, (MPI_Count)(count), dtype, source, tag, comm, request) 2565 #define MPIU_Irecv(buf, count, dtype, source, tag, comm, request) MPI_Irecv_c(buf, (MPI_Count)(count), dtype, source, tag, comm, request) 2566 #if defined(PETSC_HAVE_MPI_REDUCE_LOCAL) 2567 #define MPIU_Reduce_local(inbuf, inoutbuf, count, dtype, op) MPI_Reduce_local_c(inbuf, inoutbuf, (MPI_Count)(count), dtype, op) 2568 #endif 2569 #endif 2570 2571 PETSC_EXTERN PetscMPIInt MPIU_Allreduce_Private(const void *, void *, MPIU_Count, MPI_Datatype, MPI_Op, MPI_Comm); 2572 2573 #if defined(PETSC_USE_DEBUG) 2574 static inline unsigned int PetscStrHash(const char *str) 2575 { 2576 unsigned int c, hash = 5381; 2577 2578 while ((c = (unsigned int)*str++)) hash = ((hash << 5) + hash) + c; /* hash * 33 + c */ 2579 return hash; 2580 } 2581 #endif 2582 2583 /*MC 2584 MPIU_Allreduce - A replacement for `MPI_Allreduce()` that (1) performs single-count `MPIU_INT` operations in `PetscInt64` to detect 2585 integer overflows and (2) tries to determine if the call from all the MPI ranks occur in the 2586 same place in the PETSc code. This helps to detect bugs where different MPI ranks follow different code paths 2587 resulting in inconsistent and incorrect calls to `MPI_Allreduce()`. 2588 2589 Synopsis: 2590 #include <petscsys.h> 2591 PetscMPIInt MPIU_Allreduce(void *indata,void *outdata,PetscCount count,MPI_Datatype dtype, MPI_Op op, MPI_Comm comm); 2592 2593 Collective 2594 2595 Input Parameters: 2596 + a - pointer to the input data to be reduced 2597 . count - the number of MPI data items in `a` and `b` 2598 . dtype - the MPI datatype, for example `MPI_INT` 2599 . op - the MPI operation, for example `MPI_SUM` 2600 - comm - the MPI communicator on which the operation occurs 2601 2602 Output Parameter: 2603 . b - the reduced values 2604 2605 Level: developer 2606 2607 Note: 2608 Should be wrapped with `PetscCallMPI()` for error checking 2609 2610 .seealso: [](stylePetscCount), `MPI_Allreduce()` 2611 M*/ 2612 #if defined(PETSC_USE_DEBUG) 2613 #define MPIU_Allreduce(a, b, c, d, e, fcomm) \ 2614 PetscMacroReturnStandard( \ 2615 PetscMPIInt a_b1[6], a_b2[6]; \ 2616 int _mpiu_allreduce_c_int = (int)(c); \ 2617 a_b1[0] = -(PetscMPIInt)__LINE__; \ 2618 a_b1[1] = -a_b1[0]; \ 2619 a_b1[2] = -(PetscMPIInt)PetscStrHash(PETSC_FUNCTION_NAME); \ 2620 a_b1[3] = -a_b1[2]; \ 2621 a_b1[4] = -(PetscMPIInt)(c); \ 2622 a_b1[5] = -a_b1[4]; \ 2623 \ 2624 PetscCallMPI(MPI_Allreduce(a_b1, a_b2, 6, MPI_INT, MPI_MAX, fcomm)); \ 2625 PetscCheck(-a_b2[0] == a_b2[1], (fcomm), PETSC_ERR_PLIB, "MPIU_Allreduce() called in different locations (code lines) on different processors"); \ 2626 PetscCheck(-a_b2[2] == a_b2[3], (fcomm), PETSC_ERR_PLIB, "MPIU_Allreduce() called in different locations (functions) on different processors"); \ 2627 PetscCheck(-a_b2[4] == a_b2[5], (fcomm), PETSC_ERR_PLIB, "MPIU_Allreduce() called with different counts %d on different processors", _mpiu_allreduce_c_int); \ 2628 PetscCallMPI(MPIU_Allreduce_Private((a), (b), (c), (d), (e), (fcomm)));) 2629 #else 2630 #define MPIU_Allreduce(a, b, c, d, e, fcomm) MPIU_Allreduce_Private((a), (b), (c), (d), (e), (fcomm)) 2631 #endif 2632 2633 #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) 2634 PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint, PetscMPIInt, MPI_Info, MPI_Comm, void *, MPI_Win *); 2635 PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win, PetscMPIInt, MPI_Aint *, PetscMPIInt *, void *); 2636 #endif 2637 2638 /* this is a vile hack */ 2639 #if defined(PETSC_HAVE_NECMPI) 2640 #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) 2641 #define MPI_Type_free(a) (*(a) = MPI_DATATYPE_NULL, 0); 2642 #endif 2643 #endif 2644 2645 /* 2646 List of external packages and queries on it 2647 */ 2648 PETSC_EXTERN PetscErrorCode PetscHasExternalPackage(const char[], PetscBool *); 2649 2650 /* this cannot go here because it may be in a different shared library */ 2651 PETSC_EXTERN PetscErrorCode PCMPIServerBegin(void); 2652 PETSC_EXTERN PetscErrorCode PCMPIServerEnd(void); 2653 PETSC_EXTERN PetscBool PCMPIServerActive; 2654 PETSC_EXTERN PetscBool PCMPIServerInSolve; 2655 PETSC_EXTERN PetscBool PCMPIServerUseShmget; 2656 PETSC_EXTERN PetscErrorCode PetscShmgetAllocateArray(size_t, size_t, void **); 2657 PETSC_EXTERN PetscErrorCode PetscShmgetDeallocateArray(void **); 2658 PETSC_EXTERN PetscErrorCode PetscShmgetMapAddresses(MPI_Comm, PetscInt, const void **, void **); 2659 PETSC_EXTERN PetscErrorCode PetscShmgetUnmapAddresses(PetscInt, void **); 2660 PETSC_EXTERN PetscErrorCode PetscShmgetAddressesFinalize(void); 2661 2662 typedef struct { 2663 PetscInt n; 2664 void *addr[3]; 2665 } PCMPIServerAddresses; 2666 PETSC_EXTERN PetscErrorCode PCMPIServerAddressesDestroy(void *); 2667 2668 #define PETSC_HAVE_FORTRAN PETSC_DEPRECATED_MACRO(3, 20, 0, "PETSC_USE_FORTRAN_BINDINGS", ) PETSC_USE_FORTRAN_BINDINGS 2669 2670 PETSC_EXTERN PetscErrorCode PetscBLASSetNumThreads(PetscInt); 2671 PETSC_EXTERN PetscErrorCode PetscBLASGetNumThreads(PetscInt *); 2672 2673 /*MC 2674 PetscSafePointerPlusOffset - Checks that a pointer is not `NULL` before applying an offset 2675 2676 Level: beginner 2677 2678 Note: 2679 This is needed to avoid errors with undefined-behavior sanitizers such as 2680 UBSan, assuming PETSc has been configured with `-fsanitize=undefined` as part of the compiler flags 2681 M*/ 2682 #define PetscSafePointerPlusOffset(ptr, offset) ((ptr) ? (ptr) + (offset) : NULL) 2683