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