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