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