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