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