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 /*MC 875 PetscNewLog - Allocates memory of a type matching pointer, zeros the memory! Aligned to `PETSC_MEMALIGN`. Associates the memory allocated 876 with the given object using `PetscLogObjectMemory()`. 877 878 Synopsis: 879 #include <petscsys.h> 880 PetscErrorCode PetscNewLog(PetscObject obj,type **result) 881 882 Not Collective 883 884 Input Parameter: 885 . obj - object memory is logged to 886 887 Output Parameter: 888 . result - memory allocated, sized to match pointer type 889 890 Level: developer 891 892 .seealso: `PetscFree()`, `PetscMalloc()`, `PetscNew()`, `PetscLogObjectMemory()`, `PetscCalloc1()`, `PetscMalloc1()` 893 894 M*/ 895 #define PetscNewLog(o, b) (PetscNew((b)) || PetscLogObjectMemory((PetscObject)o, sizeof(**(b)))) 896 897 /*MC 898 PetscFree - Frees memory 899 900 Synopsis: 901 #include <petscsys.h> 902 PetscErrorCode PetscFree(void *memory) 903 904 Not Collective 905 906 Input Parameter: 907 . memory - memory to free (the pointer is ALWAYS set to NULL upon success) 908 909 Level: beginner 910 911 Notes: 912 Do not free memory obtained with `PetscMalloc2()`, `PetscCalloc2()` etc, they must be freed with `PetscFree2()` etc. 913 914 It is safe to call `PetscFree()` on a NULL pointer. 915 916 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscNewLog()`, `PetscMalloc1()`, `PetscCalloc1()` 917 918 M*/ 919 #define PetscFree(a) ((*PetscTrFree)((void *)(a), __LINE__, PETSC_FUNCTION_NAME, __FILE__) || ((a) = NULL, 0)) 920 921 /*MC 922 PetscFree2 - Frees 2 chunks of memory obtained with `PetscMalloc2()` 923 924 Synopsis: 925 #include <petscsys.h> 926 PetscErrorCode PetscFree2(void *memory1,void *memory2) 927 928 Not Collective 929 930 Input Parameters: 931 + memory1 - memory to free 932 - memory2 - 2nd memory to free 933 934 Level: developer 935 936 Notes: 937 Memory must have been obtained with `PetscMalloc2()` 938 939 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()` 940 941 M*/ 942 #define PetscFree2(m1, m2) PetscFreeA(2, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2)) 943 944 /*MC 945 PetscFree3 - Frees 3 chunks of memory obtained with `PetscMalloc3()` 946 947 Synopsis: 948 #include <petscsys.h> 949 PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3) 950 951 Not Collective 952 953 Input Parameters: 954 + memory1 - memory to free 955 . memory2 - 2nd memory to free 956 - memory3 - 3rd memory to free 957 958 Level: developer 959 960 Notes: 961 Memory must have been obtained with `PetscMalloc3()` 962 963 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()` 964 965 M*/ 966 #define PetscFree3(m1, m2, m3) PetscFreeA(3, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3)) 967 968 /*MC 969 PetscFree4 - Frees 4 chunks of memory obtained with `PetscMalloc4()` 970 971 Synopsis: 972 #include <petscsys.h> 973 PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4) 974 975 Not Collective 976 977 Input Parameters: 978 + m1 - memory to free 979 . m2 - 2nd memory to free 980 . m3 - 3rd memory to free 981 - m4 - 4th memory to free 982 983 Level: developer 984 985 Notes: 986 Memory must have been obtained with `PetscMalloc4()` 987 988 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()` 989 990 M*/ 991 #define PetscFree4(m1, m2, m3, m4) PetscFreeA(4, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4)) 992 993 /*MC 994 PetscFree5 - Frees 5 chunks of memory obtained with `PetscMalloc5()` 995 996 Synopsis: 997 #include <petscsys.h> 998 PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5) 999 1000 Not Collective 1001 1002 Input Parameters: 1003 + m1 - memory to free 1004 . m2 - 2nd memory to free 1005 . m3 - 3rd memory to free 1006 . m4 - 4th memory to free 1007 - m5 - 5th memory to free 1008 1009 Level: developer 1010 1011 Notes: 1012 Memory must have been obtained with `PetscMalloc5()` 1013 1014 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()` 1015 1016 M*/ 1017 #define PetscFree5(m1, m2, m3, m4, m5) PetscFreeA(5, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5)) 1018 1019 /*MC 1020 PetscFree6 - Frees 6 chunks of memory obtained with `PetscMalloc6()` 1021 1022 Synopsis: 1023 #include <petscsys.h> 1024 PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6) 1025 1026 Not Collective 1027 1028 Input Parameters: 1029 + m1 - memory to free 1030 . m2 - 2nd memory to free 1031 . m3 - 3rd memory to free 1032 . m4 - 4th memory to free 1033 . m5 - 5th memory to free 1034 - m6 - 6th memory to free 1035 1036 Level: developer 1037 1038 Notes: 1039 Memory must have been obtained with `PetscMalloc6()` 1040 1041 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()` 1042 1043 M*/ 1044 #define PetscFree6(m1, m2, m3, m4, m5, m6) PetscFreeA(6, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6)) 1045 1046 /*MC 1047 PetscFree7 - Frees 7 chunks of memory obtained with `PetscMalloc7()` 1048 1049 Synopsis: 1050 #include <petscsys.h> 1051 PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7) 1052 1053 Not Collective 1054 1055 Input Parameters: 1056 + m1 - memory to free 1057 . m2 - 2nd memory to free 1058 . m3 - 3rd memory to free 1059 . m4 - 4th memory to free 1060 . m5 - 5th memory to free 1061 . m6 - 6th memory to free 1062 - m7 - 7th memory to free 1063 1064 Level: developer 1065 1066 Notes: 1067 Memory must have been obtained with `PetscMalloc7()` 1068 1069 .seealso: `PetscNew()`, `PetscMalloc()`, `PetscMalloc2()`, `PetscFree()`, `PetscMalloc3()`, `PetscMalloc4()`, `PetscMalloc5()`, `PetscMalloc6()`, 1070 `PetscMalloc7()` 1071 1072 M*/ 1073 #define PetscFree7(m1, m2, m3, m4, m5, m6, m7) PetscFreeA(7, __LINE__, PETSC_FUNCTION_NAME, __FILE__, &(m1), &(m2), &(m3), &(m4), &(m5), &(m6), &(m7)) 1074 1075 PETSC_EXTERN PetscErrorCode PetscMallocA(int, PetscBool, int, const char *, const char *, size_t, void *, ...); 1076 PETSC_EXTERN PetscErrorCode PetscFreeA(int, int, const char *, const char *, void *, ...); 1077 PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t, PetscBool, int, const char[], const char[], void **); 1078 PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void *, int, const char[], const char[]); 1079 PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t, int, const char[], const char[], void **); 1080 PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool); 1081 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 **)); 1082 PETSC_EXTERN PetscErrorCode PetscMallocClear(void); 1083 1084 /* 1085 Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair. 1086 */ 1087 PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void); 1088 PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void); 1089 #if defined(PETSC_HAVE_CUDA) 1090 PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void); 1091 PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void); 1092 #endif 1093 #if defined(PETSC_HAVE_HIP) 1094 PETSC_EXTERN PetscErrorCode PetscMallocSetHIPHost(void); 1095 PETSC_EXTERN PetscErrorCode PetscMallocResetHIPHost(void); 1096 #endif 1097 1098 #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE 1099 #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION 1100 1101 /* 1102 Routines for tracing memory corruption/bleeding with default PETSc memory allocation 1103 */ 1104 PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *); 1105 PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *); 1106 PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *); 1107 PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *); 1108 PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int); 1109 PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int, PetscLogDouble *); 1110 PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool, PetscBool); 1111 PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool *, PetscBool *, PetscBool *); 1112 PETSC_EXTERN PetscErrorCode PetscMallocValidate(int, const char[], const char[]); 1113 PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble); 1114 PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool *); 1115 PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool); 1116 PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool *); 1117 1118 PETSC_EXTERN const char *const PetscDataTypes[]; 1119 PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType, MPI_Datatype *); 1120 PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype, PetscDataType *); 1121 PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType, size_t *); 1122 PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char *, PetscDataType *, PetscBool *); 1123 1124 /* 1125 Basic memory and string operations. These are usually simple wrappers 1126 around the basic Unix system calls, but a few of them have additional 1127 functionality and/or error checking. 1128 */ 1129 PETSC_EXTERN PetscErrorCode PetscMemcmp(const void *, const void *, size_t, PetscBool *); 1130 PETSC_EXTERN PetscErrorCode PetscStrlen(const char[], size_t *); 1131 PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[], char, int *, char ***); 1132 PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int, char **); 1133 PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[], const char[], PetscBool *); 1134 PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[], const char[], PetscBool *); 1135 PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[], const char[], PetscBool *); 1136 PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[], const char[], size_t, PetscBool *); 1137 PETSC_EXTERN PetscErrorCode PetscStrcpy(char[], const char[]); 1138 PETSC_EXTERN PetscErrorCode PetscStrcat(char[], const char[]); 1139 PETSC_EXTERN PetscErrorCode PetscStrlcat(char[], const char[], size_t); 1140 PETSC_EXTERN PetscErrorCode PetscStrncpy(char[], const char[], size_t); 1141 PETSC_EXTERN PetscErrorCode PetscStrchr(const char[], char, char *[]); 1142 PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]); 1143 PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]); 1144 PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[], char, char *[]); 1145 PETSC_EXTERN PetscErrorCode PetscStrstr(const char[], const char[], char *[]); 1146 PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[], const char[], char *[]); 1147 PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[], const char[], PetscBool *); 1148 PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[], const char[], PetscBool *); 1149 PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[], const char *const *, PetscInt *); 1150 PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[], char *[]); 1151 PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const *, char ***); 1152 PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char ***); 1153 PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt, const char *const *, char ***); 1154 PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt, char ***); 1155 PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm, const char[], char[], size_t); 1156 1157 PETSC_EXTERN void PetscStrcmpNoError(const char[], const char[], PetscBool *); 1158 1159 PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[], const char, PetscToken *); 1160 PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken, char *[]); 1161 PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken *); 1162 1163 PETSC_EXTERN PetscErrorCode PetscStrInList(const char[], const char[], char, PetscBool *); 1164 PETSC_EXTERN const char *PetscBasename(const char[]); 1165 PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt, const char *const *, const char *, PetscInt *, PetscBool *); 1166 PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const *, const char *, PetscEnum *, PetscBool *); 1167 1168 /* 1169 These are MPI operations for MPI_Allreduce() etc 1170 */ 1171 PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP; 1172 #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 1173 PETSC_EXTERN MPI_Op MPIU_SUM; 1174 PETSC_EXTERN MPI_Op MPIU_MAX; 1175 PETSC_EXTERN MPI_Op MPIU_MIN; 1176 #else 1177 #define MPIU_SUM MPI_SUM 1178 #define MPIU_MAX MPI_MAX 1179 #define MPIU_MIN MPI_MIN 1180 #endif 1181 #if defined(PETSC_HAVE_REAL___FLOAT128) || defined(PETSC_HAVE_REAL___FP16) 1182 /*MC 1183 MPIU_SUM___FP16___FLOAT128 - MPI_Op that acts as a replacement for MPI_SUM with 1184 custom MPI_Datatype `MPIU___FLOAT128`, `MPIU___COMPLEX128`, and `MPIU___FP16`. 1185 1186 Level: advanced 1187 1188 Developer Note: 1189 This should be unified with `MPIU_SUM` 1190 1191 .seealso: `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX` 1192 M*/ 1193 PETSC_EXTERN MPI_Op MPIU_SUM___FP16___FLOAT128; 1194 #endif 1195 PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm, const PetscInt[], PetscInt *, PetscInt *); 1196 1197 PETSC_EXTERN PetscErrorCode MPIULong_Send(void *, PetscInt, MPI_Datatype, PetscMPIInt, PetscMPIInt, MPI_Comm) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(1, 3); 1198 PETSC_EXTERN PetscErrorCode MPIULong_Recv(void *, PetscInt, MPI_Datatype, PetscMPIInt, PetscMPIInt, MPI_Comm) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(1, 3); 1199 1200 PETSC_EXTERN const char *const PetscFileModes[]; 1201 1202 /* 1203 Defines PETSc error handling. 1204 */ 1205 #include <petscerror.h> 1206 1207 PETSC_EXTERN PetscBool PetscCIEnabled; /* code is running in the PETSc test harness CI */ 1208 PETSC_EXTERN PetscBool PetscCIEnabledPortableErrorOutput; /* error output is stripped to ensure portability of error messages across systems */ 1209 PETSC_EXTERN const char *PetscCIFilename(const char *); 1210 PETSC_EXTERN int PetscCILinenumber(int); 1211 1212 #define PETSC_SMALLEST_CLASSID 1211211 1213 PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID; 1214 PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID; 1215 PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[], PetscClassId *); 1216 PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject, PetscObjectId *); 1217 PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject, PetscObjectId, PetscBool *); 1218 1219 /* 1220 Routines that get memory usage information from the OS 1221 */ 1222 PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *); 1223 PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *); 1224 PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void); 1225 PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]); 1226 1227 PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal); 1228 1229 /* 1230 Initialization of PETSc 1231 */ 1232 PETSC_EXTERN PetscErrorCode PetscInitialize(int *, char ***, const char[], const char[]); 1233 PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int, char **, const char[], const char[]); 1234 PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void); 1235 PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *); 1236 PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *); 1237 PETSC_EXTERN PetscErrorCode PetscFinalize(void); 1238 PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void); 1239 PETSC_EXTERN PetscErrorCode PetscGetArgs(int *, char ***); 1240 PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***); 1241 PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **); 1242 1243 PETSC_EXTERN PetscErrorCode PetscEnd(void); 1244 PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void); 1245 1246 PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[], const char[]); 1247 PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void); 1248 PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void); 1249 PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject, const char[]); 1250 1251 PETSC_EXTERN PetscErrorCode PetscMonitorCompare(PetscErrorCode (*)(void), void *, PetscErrorCode (*)(void **), PetscErrorCode (*)(void), void *, PetscErrorCode (*)(void **), PetscBool *); 1252 1253 /* 1254 These are so that in extern C code we can caste function pointers to non-extern C 1255 function pointers. Since the regular C++ code expects its function pointers to be C++ 1256 */ 1257 PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void); 1258 PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void); 1259 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void); 1260 1261 /* 1262 Functions that can act on any PETSc object. 1263 */ 1264 PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject *); 1265 PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject, MPI_Comm *); 1266 PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject, PetscClassId *); 1267 PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject, const char *[]); 1268 PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject, const char *[]); 1269 PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject, const char[]); 1270 PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject, const char *[]); 1271 PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject, PetscInt); 1272 PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject, PetscInt *); 1273 PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject, PetscObject, PetscInt); 1274 PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject); 1275 PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject, PetscInt *); 1276 PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject); 1277 PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject, PetscMPIInt *); 1278 PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject, const char[], PetscObject); 1279 PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject, const char[]); 1280 PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject, const char[], PetscObject *); 1281 PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject, const char[], void (*)(void)); 1282 #define PetscObjectComposeFunction(a, b, d) PetscObjectComposeFunction_Private(a, b, (PetscVoidFunction)(d)) 1283 PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject); 1284 PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject); 1285 PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject); 1286 PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject, PetscObject); 1287 PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm, PetscMPIInt *); 1288 1289 #include <petscviewertypes.h> 1290 #include <petscoptions.h> 1291 1292 PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer, PetscBool, PetscLogDouble); 1293 PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool *); 1294 1295 PETSC_EXTERN PetscErrorCode PetscObjectsListGetGlobalNumbering(MPI_Comm, PetscInt, PetscObject *, PetscInt *, PetscInt *); 1296 1297 PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer, const char[]); 1298 PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject, PetscViewer); 1299 PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject, PetscViewer); 1300 #define PetscObjectQueryFunction(obj, name, fptr) PetscObjectQueryFunction_Private((obj), (name), (PetscVoidFunction *)(fptr)) 1301 PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject, const char[], void (**)(void)); 1302 PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject, const char[]); 1303 PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject, const char[]); 1304 PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject, const char[]); 1305 PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject, const char *[]); 1306 PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject, const char[]); 1307 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject); 1308 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void); 1309 PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject, PetscObject, const char[]); 1310 PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject); 1311 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject, const char[], PetscBool *); 1312 PETSC_EXTERN PetscErrorCode PetscObjectObjectTypeCompare(PetscObject, PetscObject, PetscBool *); 1313 PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject, const char[], PetscBool *); 1314 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject, PetscBool *, const char[], ...); 1315 PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject, PetscBool *, const char[], ...); 1316 PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void)); 1317 PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void); 1318 1319 #if defined(PETSC_HAVE_SAWS) 1320 PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void); 1321 PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject); 1322 PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject, PetscBool); 1323 PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject); 1324 PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject); 1325 PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject); 1326 PETSC_EXTERN void PetscStackSAWsGrantAccess(void); 1327 PETSC_EXTERN void PetscStackSAWsTakeAccess(void); 1328 PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void); 1329 PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void); 1330 1331 #else 1332 #define PetscSAWsBlock() 0 1333 #define PetscObjectSAWsViewOff(obj) 0 1334 #define PetscObjectSAWsSetBlock(obj, flg) 0 1335 #define PetscObjectSAWsBlock(obj) 0 1336 #define PetscObjectSAWsGrantAccess(obj) 0 1337 #define PetscObjectSAWsTakeAccess(obj) 0 1338 #define PetscStackViewSAWs() 0 1339 #define PetscStackSAWsViewOff() 0 1340 #define PetscStackSAWsTakeAccess() 1341 #define PetscStackSAWsGrantAccess() 1342 1343 #endif 1344 1345 PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[], PetscDLMode, PetscDLHandle *); 1346 PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *); 1347 PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle, const char[], void **); 1348 PETSC_EXTERN PetscErrorCode PetscDLAddr(void (*)(void), char **); 1349 #ifdef PETSC_HAVE_CXX 1350 PETSC_EXTERN PetscErrorCode PetscDemangleSymbol(const char *, char **); 1351 #endif 1352 1353 #if defined(PETSC_USE_DEBUG) 1354 PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void *, PetscStack **); 1355 #endif 1356 PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE *, PetscBool); 1357 1358 PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList *); 1359 PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList, const char[], PetscObject *); 1360 PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList, PetscObject, char **, PetscBool *); 1361 PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *, const char[], PetscObject); 1362 PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *, const char[]); 1363 PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList, PetscObjectList *); 1364 1365 /* 1366 Dynamic library lists. Lists of names of routines in objects or in dynamic 1367 link libraries that will be loaded as needed. 1368 */ 1369 1370 #define PetscFunctionListAdd(list, name, fptr) PetscFunctionListAdd_Private((list), (name), (PetscVoidFunction)(fptr)) 1371 PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList *, const char[], void (*)(void)); 1372 PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList *); 1373 PETSC_EXTERN PetscErrorCode PetscFunctionListClear(PetscFunctionList); 1374 #define PetscFunctionListFind(list, name, fptr) PetscFunctionListFind_Private((list), (name), (PetscVoidFunction *)(fptr)) 1375 PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList, const char[], void (**)(void)); 1376 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm, FILE *, const char[], const char[], const char[], const char[], PetscFunctionList, const char[], const char[]); 1377 PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList, PetscFunctionList *); 1378 PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList, PetscViewer); 1379 PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList, const char ***, int *); 1380 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintNonEmpty(PetscFunctionList); 1381 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintAll(void); 1382 1383 PETSC_EXTERN PetscDLLibrary PetscDLLibrariesLoaded; 1384 PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm, PetscDLLibrary *, const char[]); 1385 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm, PetscDLLibrary *, const char[]); 1386 PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm, PetscDLLibrary *, const char[], const char[], void **); 1387 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary); 1388 PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm, const char[], char *, size_t, PetscBool *); 1389 PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm, const char[], PetscDLLibrary *); 1390 PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary); 1391 1392 /* 1393 Useful utility routines 1394 */ 1395 PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm, PetscInt *, PetscInt *); 1396 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm, PetscInt, PetscInt *, PetscInt *); 1397 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm, PetscInt *, PetscInt *); 1398 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm, PetscMPIInt); 1399 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm, PetscMPIInt); 1400 PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject); 1401 PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE *); 1402 PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm, const PetscInt[2], PetscInt[2]); 1403 PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm, const PetscReal[2], PetscReal[2]); 1404 1405 /*MC 1406 PetscNot - negates a logical type value and returns result as a `PetscBool` 1407 1408 Notes: 1409 This is useful in cases like 1410 $ int *a; 1411 $ PetscBool flag = PetscNot(a) 1412 where !a would not return a PetscBool because we cannot provide a cast from int to PetscBool in C. 1413 1414 Level: beginner 1415 1416 .seealso `:` `PetscBool`, `PETSC_TRUE`, `PETSC_FALSE` 1417 M*/ 1418 #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE) 1419 1420 /*MC 1421 PetscHelpPrintf - Prints help messages. 1422 1423 Synopsis: 1424 #include <petscsys.h> 1425 PetscErrorCode (*PetscHelpPrintf)(MPI_Comm comm, const char format[],args); 1426 1427 Collective on comm 1428 1429 Input Parameters: 1430 + comm - the MPI communicator over which the help message is printed 1431 . format - the usual printf() format string 1432 - args - arguments to be printed 1433 1434 Level: developer 1435 1436 Fortran Note: 1437 This routine is not supported in Fortran. 1438 1439 Note: 1440 You can change how help messages are printed by replacing the function pointer with a function that does not simply write to stdout. 1441 1442 To use, write your own function, for example, 1443 $PetscErrorCode mypetschelpprintf(MPI_Comm comm,const char format[],....) 1444 ${ 1445 $ PetscFunctionReturn(0); 1446 $} 1447 then do the assigment 1448 $ PetscHelpPrintf = mypetschelpprintf; 1449 You can do the assignment before `PetscInitialize()`. 1450 1451 The default routine used is called `PetscHelpPrintfDefault()`. 1452 1453 .seealso: `PetscFPrintf()`, `PetscSynchronizedPrintf()`, `PetscErrorPrintf()` 1454 M*/ 1455 PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3); 1456 1457 /* 1458 Defines PETSc profiling. 1459 */ 1460 #include <petsclog.h> 1461 1462 /* 1463 Simple PETSc parallel IO for ASCII printing 1464 */ 1465 PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[], char[]); 1466 PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm, const char[], const char[], FILE **); 1467 PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm, FILE *); 1468 PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4); 1469 PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3); 1470 PETSC_EXTERN PetscErrorCode PetscSNPrintf(char *, size_t, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4); 1471 PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char *, size_t, const char[], size_t *, ...) PETSC_ATTRIBUTE_FORMAT(3, 5); 1472 PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[], size_t, const char *, PetscInt, const PetscReal[]); 1473 1474 PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2); 1475 PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char[], ...) PETSC_ATTRIBUTE_FORMAT(1, 2); 1476 PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3); 1477 1478 PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char *, size_t *); 1479 PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char *, char *); 1480 1481 #if defined(PETSC_HAVE_POPEN) 1482 PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm, const char[], const char[], const char[], FILE **); 1483 PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm, FILE *); 1484 PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]); 1485 #endif 1486 1487 PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm, const char[], ...) PETSC_ATTRIBUTE_FORMAT(2, 3); 1488 PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm, FILE *, const char[], ...) PETSC_ATTRIBUTE_FORMAT(3, 4); 1489 PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm, FILE *); 1490 PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm, FILE *, size_t, char[]); 1491 PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm, const char[], const char[], FILE **); 1492 PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm, const char[], const char[], FILE **); 1493 PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char *[]); 1494 1495 PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID; 1496 PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer, void **); 1497 PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer, void *); 1498 PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer *); 1499 PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm, PetscContainer *); 1500 PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void *)); 1501 PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void *); 1502 1503 /* 1504 For use in debuggers 1505 */ 1506 PETSC_EXTERN PetscMPIInt PetscGlobalRank; 1507 PETSC_EXTERN PetscMPIInt PetscGlobalSize; 1508 PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt, const PetscInt[], PetscViewer); 1509 PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt, const PetscReal[], PetscViewer); 1510 PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt, const PetscScalar[], PetscViewer); 1511 1512 #include <stddef.h> 1513 #include <string.h> /* for memcpy, memset */ 1514 #include <stdlib.h> 1515 1516 #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__) 1517 #include <xmmintrin.h> 1518 #endif 1519 1520 /*@C 1521 PetscMemmove - Copies n bytes, beginning at location b, to the space 1522 beginning at location a. Copying between regions that overlap will 1523 take place correctly. Use `PetscMemcpy()` if the locations do not overlap 1524 1525 Not Collective 1526 1527 Input Parameters: 1528 + b - pointer to initial memory space 1529 . a - pointer to copy space 1530 - n - length (in bytes) of space to copy 1531 1532 Level: intermediate 1533 1534 Notes: 1535 `PetscArraymove()` is preferred 1536 1537 This routine is analogous to `memmove()`. 1538 1539 Developers Note: 1540 This is inlined for performance 1541 1542 .seealso: `PetscMemcpy()`, `PetscMemcmp()`, `PetscArrayzero()`, `PetscMemzero()`, `PetscArraycmp()`, `PetscArraycpy()`, `PetscStrallocpy()`, 1543 `PetscArraymove()` 1544 @*/ 1545 static inline PetscErrorCode PetscMemmove(void *a, const void *b, size_t n) { 1546 PetscFunctionBegin; 1547 if (n > 0) { 1548 PetscCheck(a, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Trying to copy to null pointer"); 1549 PetscCheck(b, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Trying to copy from a null pointer"); 1550 } 1551 #if !defined(PETSC_HAVE_MEMMOVE) 1552 if (a < b) { 1553 if (a <= b - n) memcpy(a, b, n); 1554 else { 1555 memcpy(a, b, (int)(b - a)); 1556 PetscMemmove(b, b + (int)(b - a), n - (int)(b - a)); 1557 } 1558 } else { 1559 if (b <= a - n) memcpy(a, b, n); 1560 else { 1561 memcpy(b + n, b + (n - (int)(a - b)), (int)(a - b)); 1562 PetscMemmove(a, b, n - (int)(a - b)); 1563 } 1564 } 1565 #else 1566 memmove((char *)(a), (char *)(b), n); 1567 #endif 1568 PetscFunctionReturn(0); 1569 } 1570 1571 /*@C 1572 PetscMemcpy - Copies n bytes, beginning at location b, to the space 1573 beginning at location a. The two memory regions CANNOT overlap, use 1574 `PetscMemmove()` in that case. 1575 1576 Not Collective 1577 1578 Input Parameters: 1579 + b - pointer to initial memory space 1580 - n - length (in bytes) of space to copy 1581 1582 Output Parameter: 1583 . a - pointer to copy space 1584 1585 Level: intermediate 1586 1587 Compile Option: 1588 `PETSC_PREFER_DCOPY_FOR_MEMCPY` will cause the BLAS dcopy() routine to be used 1589 for memory copies on double precision values. 1590 `PETSC_PREFER_COPY_FOR_MEMCPY` will cause C code to be used 1591 for memory copies on double precision values. 1592 `PETSC_PREFER_FORTRAN_FORMEMCPY` will cause Fortran code to be used 1593 for memory copies on double precision values. 1594 1595 Notes: 1596 Prefer `PetscArraycpy()` 1597 1598 This routine is analogous to `memcpy()`. 1599 1600 Not available from Fortran 1601 1602 Developer Note: 1603 This is inlined for fastest performance 1604 1605 .seealso: `PetscMemzero()`, `PetscMemcmp()`, `PetscArrayzero()`, `PetscArraycmp()`, `PetscArraycpy()`, `PetscMemmove()`, `PetscStrallocpy()` 1606 1607 @*/ 1608 static inline PetscErrorCode PetscMemcpy(void *a, const void *b, size_t n) { 1609 #if defined(PETSC_USE_DEBUG) 1610 size_t al = (size_t)a, bl = (size_t)b; 1611 size_t nl = (size_t)n; 1612 PetscFunctionBegin; 1613 if (n > 0) { 1614 PetscCheck(b, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Trying to copy from a null pointer"); 1615 PetscCheck(a, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Trying to copy to a null pointer"); 1616 } 1617 #else 1618 PetscFunctionBegin; 1619 #endif 1620 if (a != b && n > 0) { 1621 #if defined(PETSC_USE_DEBUG) 1622 PetscCheck(!((al > bl && (al - bl) < nl) || (bl - al) < nl), PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Memory regions overlap: either use PetscMemmov()\n\ 1623 or make sure your copy regions and lengths are correct. \n\ 1624 Length (bytes) %zu first address %zu second address %zu", 1625 nl, al, bl); 1626 #endif 1627 #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY)) 1628 if (!(a % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) { 1629 size_t len = n / sizeof(PetscScalar); 1630 #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) 1631 PetscBLASInt one = 1, blen; 1632 PetscCall(PetscBLASIntCast(len, &blen)); 1633 PetscCallBLAS("BLAScopy", BLAScopy_(&blen, (PetscScalar *)b, &one, (PetscScalar *)a, &one)); 1634 #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY) 1635 fortrancopy_(&len, (PetscScalar *)b, (PetscScalar *)a); 1636 #else 1637 PetscScalar *x = (PetscScalar *)b, *y = (PetscScalar *)a; 1638 for (size_t i = 0; i < len; i++) y[i] = x[i]; 1639 #endif 1640 } else { 1641 memcpy((char *)(a), (char *)(b), n); 1642 } 1643 #else 1644 memcpy((char *)(a), (char *)(b), n); 1645 #endif 1646 } 1647 PetscFunctionReturn(0); 1648 } 1649 1650 /*@C 1651 PetscMemzero - Zeros the specified memory. 1652 1653 Not Collective 1654 1655 Input Parameters: 1656 + a - pointer to beginning memory location 1657 - n - length (in bytes) of memory to initialize 1658 1659 Level: intermediate 1660 1661 Compile Option: 1662 `PETSC_PREFER_BZERO` - on certain machines (the IBM RS6000) the bzero() routine happens 1663 to be faster than the memset() routine. This flag causes the bzero() routine to be used. 1664 1665 Notes: 1666 Not available from Fortran 1667 1668 Prefer `PetscArrayzero()` 1669 1670 Developer Note: 1671 This is inlined for fastest performance 1672 1673 .seealso: `PetscMemcpy()`, `PetscMemcmp()`, `PetscArrayzero()`, `PetscArraycmp()`, `PetscArraycpy()`, `PetscMemmove()`, `PetscStrallocpy()` 1674 @*/ 1675 static inline PetscErrorCode PetscMemzero(void *a, size_t n) { 1676 if (n > 0) { 1677 #if defined(PETSC_USE_DEBUG) 1678 PetscCheck(a, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Trying to zero at a null pointer with %zu bytes", n); 1679 #endif 1680 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) 1681 if (!(((long)a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) { 1682 size_t i, len = n / sizeof(PetscScalar); 1683 PetscScalar *x = (PetscScalar *)a; 1684 for (i = 0; i < len; i++) x[i] = 0.0; 1685 } else { 1686 #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO) 1687 if (!(((long)a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) { 1688 PetscInt len = n / sizeof(PetscScalar); 1689 fortranzero_(&len, (PetscScalar *)a); 1690 } else { 1691 #endif 1692 #if defined(PETSC_PREFER_BZERO) 1693 bzero((char *)a, n); 1694 #else 1695 memset((char *)a, 0, n); 1696 #endif 1697 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO) 1698 } 1699 #endif 1700 } 1701 return 0; 1702 } 1703 1704 /*MC 1705 PetscArraycmp - Compares two arrays in memory. 1706 1707 Synopsis: 1708 #include <petscsys.h> 1709 PetscErrorCode PetscArraycmp(const anytype *str1,const anytype *str2,size_t cnt,PetscBool *e) 1710 1711 Not Collective 1712 1713 Input Parameters: 1714 + str1 - First array 1715 . str2 - Second array 1716 - cnt - Count of the array, not in bytes, but number of entries in the arrays 1717 1718 Output Parameters: 1719 . e - `PETSC_TRUE` if equal else `PETSC_FALSE`. 1720 1721 Level: intermediate 1722 1723 Notes: 1724 This routine is a preferred replacement to `PetscMemcmp()` 1725 1726 The arrays must be of the same type 1727 1728 .seealso: `PetscMemcpy()`, `PetscMemcmp()`, `PetscArrayzero()`, `PetscMemzero()`, `PetscArraycpy()`, `PetscMemmove()`, `PetscStrallocpy()`, 1729 `PetscArraymove()` 1730 M*/ 1731 #define PetscArraycmp(str1, str2, cnt, e) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcmp(str1, str2, (size_t)(cnt) * sizeof(*(str1)), e)) 1732 1733 /*MC 1734 PetscArraymove - Copies from one array in memory to another, the arrays may overlap. Use `PetscArraycpy()` when the arrays 1735 do not overlap 1736 1737 Synopsis: 1738 #include <petscsys.h> 1739 PetscErrorCode PetscArraymove(anytype *str1,const anytype *str2,size_t cnt) 1740 1741 Not Collective 1742 1743 Input Parameters: 1744 + str1 - First array 1745 . str2 - Second array 1746 - cnt - Count of the array, not in bytes, but number of entries in the arrays 1747 1748 Level: intermediate 1749 1750 Notes: 1751 This routine is a preferred replacement to `PetscMemmove()` 1752 1753 The arrays must be of the same type 1754 1755 .seealso: `PetscMemcpy()`, `PetscMemcmp()`, `PetscArrayzero()`, `PetscMemzero()`, `PetscArraycpy()`, `PetscMemmove()`, `PetscArraycmp()`, `PetscStrallocpy()` 1756 M*/ 1757 #define PetscArraymove(str1, str2, cnt) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemmove(str1, str2, (size_t)(cnt) * sizeof(*(str1)))) 1758 1759 /*MC 1760 PetscArraycpy - Copies from one array in memory to another 1761 1762 Synopsis: 1763 #include <petscsys.h> 1764 PetscErrorCode PetscArraycpy(anytype *str1,const anytype *str2,size_t cnt) 1765 1766 Not Collective 1767 1768 Input Parameters: 1769 + str1 - First array (destination) 1770 . str2 - Second array (source) 1771 - cnt - Count of the array, not in bytes, but number of entries in the arrays 1772 1773 Level: intermediate 1774 1775 Notes: 1776 This routine is a preferred replacement to `PetscMemcpy()` 1777 1778 The arrays must be of the same type 1779 1780 .seealso: `PetscMemcpy()`, `PetscMemcmp()`, `PetscArrayzero()`, `PetscMemzero()`, `PetscArraymove()`, `PetscMemmove()`, `PetscArraycmp()`, `PetscStrallocpy()` 1781 M*/ 1782 #define PetscArraycpy(str1, str2, cnt) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcpy(str1, str2, (size_t)(cnt) * sizeof(*(str1)))) 1783 1784 /*MC 1785 PetscArrayzero - Zeros an array in memory. 1786 1787 Synopsis: 1788 #include <petscsys.h> 1789 PetscErrorCode PetscArrayzero(anytype *str1,size_t cnt) 1790 1791 Not Collective 1792 1793 Input Parameters: 1794 + str1 - array 1795 - cnt - Count of the array, not in bytes, but number of entries in the array 1796 1797 Level: intermediate 1798 1799 Note: 1800 This routine is a preferred replacement to `PetscMemzero()` 1801 1802 .seealso: `PetscMemcpy()`, `PetscMemcmp()`, `PetscMemzero()`, `PetscArraycmp()`, `PetscArraycpy()`, `PetscMemmove()`, `PetscStrallocpy()`, `PetscArraymove()` 1803 M*/ 1804 #define PetscArrayzero(str1, cnt) PetscMemzero(str1, (size_t)(cnt) * sizeof(*(str1))) 1805 1806 /*MC 1807 PetscPrefetchBlock - Prefetches a block of memory 1808 1809 Synopsis: 1810 #include <petscsys.h> 1811 void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t) 1812 1813 Not Collective 1814 1815 Input Parameters: 1816 + a - pointer to first element to fetch (any type but usually PetscInt or PetscScalar) 1817 . n - number of elements to fetch 1818 . rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors) 1819 - t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note 1820 1821 Level: developer 1822 1823 Notes: 1824 The last two arguments (rw and t) must be compile-time constants. 1825 1826 Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality. Not all architectures offer 1827 equivalent locality hints, but the following macros are always defined to their closest analogue. 1828 + `PETSC_PREFETCH_HINT_NTA` - Non-temporal. Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched). 1829 . `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. 1830 . `PETSC_PREFETCH_HINT_T1` - Fetch to level 2 and higher (not L1). 1831 - `PETSC_PREFETCH_HINT_T2` - Fetch to high-level cache only. (On many systems, T0 and T1 are equivalent.) 1832 1833 This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid 1834 address). 1835 1836 M*/ 1837 #define PetscPrefetchBlock(a, n, rw, t) \ 1838 do { \ 1839 const char *_p = (const char *)(a), *_end = (const char *)((a) + (n)); \ 1840 for (; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p, (rw), (t)); \ 1841 } while (0) 1842 1843 /* 1844 Determine if some of the kernel computation routines use 1845 Fortran (rather than C) for the numerical calculations. On some machines 1846 and compilers (like complex numbers) the Fortran version of the routines 1847 is faster than the C/C++ versions. The flag --with-fortran-kernels 1848 should be used with ./configure to turn these on. 1849 */ 1850 #if defined(PETSC_USE_FORTRAN_KERNELS) 1851 1852 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL) 1853 #define PETSC_USE_FORTRAN_KERNEL_MULTCRL 1854 #endif 1855 1856 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM) 1857 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM 1858 #endif 1859 1860 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ) 1861 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ 1862 #endif 1863 1864 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ) 1865 #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ 1866 #endif 1867 1868 #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM) 1869 #define PETSC_USE_FORTRAN_KERNEL_NORM 1870 #endif 1871 1872 #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY) 1873 #define PETSC_USE_FORTRAN_KERNEL_MAXPY 1874 #endif 1875 1876 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 1877 #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ 1878 #endif 1879 1880 #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ) 1881 #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ 1882 #endif 1883 1884 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ) 1885 #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ 1886 #endif 1887 1888 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ) 1889 #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ 1890 #endif 1891 1892 #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT) 1893 #define PETSC_USE_FORTRAN_KERNEL_MDOT 1894 #endif 1895 1896 #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY) 1897 #define PETSC_USE_FORTRAN_KERNEL_XTIMESY 1898 #endif 1899 1900 #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX) 1901 #define PETSC_USE_FORTRAN_KERNEL_AYPX 1902 #endif 1903 1904 #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY) 1905 #define PETSC_USE_FORTRAN_KERNEL_WAXPY 1906 #endif 1907 1908 #endif 1909 1910 /* 1911 Macros for indicating code that should be compiled with a C interface, 1912 rather than a C++ interface. Any routines that are dynamically loaded 1913 (such as the PCCreate_XXX() routines) must be wrapped so that the name 1914 mangler does not change the functions symbol name. This just hides the 1915 ugly extern "C" {} wrappers. 1916 */ 1917 #if defined(__cplusplus) 1918 #define EXTERN_C_BEGIN extern "C" { 1919 #define EXTERN_C_END } 1920 #else 1921 #define EXTERN_C_BEGIN 1922 #define EXTERN_C_END 1923 #endif 1924 1925 /* --------------------------------------------------------------------*/ 1926 1927 /*MC 1928 MPI_Comm - the basic object used by MPI to determine which processes are involved in a 1929 communication 1930 1931 Level: beginner 1932 1933 Note: This manual page is a place-holder because MPICH does not have a manual page for MPI_Comm 1934 1935 .seealso: `PETSC_COMM_WORLD`, `PETSC_COMM_SELF` 1936 M*/ 1937 1938 #if defined(PETSC_HAVE_MPIIO) 1939 PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4); 1940 PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File, void *, PetscMPIInt, MPI_Datatype, MPI_Status *) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(2, 4); 1941 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); 1942 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); 1943 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); 1944 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); 1945 #endif 1946 1947 /* the following petsc_static_inline require petscerror.h */ 1948 1949 /* Limit MPI to 32-bits */ 1950 #define PETSC_MPI_INT_MAX 2147483647 1951 #define PETSC_MPI_INT_MIN -2147483647 1952 /* Limit BLAS to 32-bits */ 1953 #define PETSC_BLAS_INT_MAX 2147483647 1954 #define PETSC_BLAS_INT_MIN -2147483647 1955 #define PETSC_CUBLAS_INT_MAX 2147483647 1956 1957 /*@C 1958 PetscIntCast - casts a `PetscInt64` (which is 64 bits in size) to a `PetscInt` (which may be 32 bits in size), generates an 1959 error if the `PetscInt` is not large enough to hold the number. 1960 1961 Not Collective 1962 1963 Input Parameter: 1964 . a - the `PetscInt64` value 1965 1966 Output Parameter: 1967 . b - the resulting `PetscInt` value 1968 1969 Level: advanced 1970 1971 Notes: 1972 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 1973 1974 Not available from Fortran 1975 1976 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscIntMultError()`, `PetscIntSumError()` 1977 @*/ 1978 static inline PetscErrorCode PetscIntCast(PetscInt64 a, PetscInt *b) { 1979 PetscFunctionBegin; 1980 #if !defined(PETSC_USE_64BIT_INDICES) 1981 if (a > PETSC_MAX_INT) { 1982 *b = 0; 1983 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); 1984 } 1985 #endif 1986 *b = (PetscInt)(a); 1987 PetscFunctionReturn(0); 1988 } 1989 1990 /*@C 1991 PetscCountCast - casts a `PetscCount` to a `PetscInt` (which may be 32 bits in size), generates an 1992 error if the `PetscInt` is not large enough to hold the number. 1993 1994 Not Collective 1995 1996 Input Parameter: 1997 . a - the `PetscCount` value 1998 1999 Output Parameter: 2000 . b - the resulting `PetscInt` value 2001 2002 Level: advanced 2003 2004 Notes: 2005 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 2006 2007 Not available from Fortran 2008 2009 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscBLASIntCast()`, `PetscIntMultError()`, `PetscIntSumError()`, `PetscIntCast()` 2010 @*/ 2011 static inline PetscErrorCode PetscCountCast(PetscCount a, PetscInt *b) { 2012 PetscFunctionBegin; 2013 if (sizeof(PetscCount) > sizeof(PetscInt) && a > PETSC_MAX_INT) { 2014 *b = 0; 2015 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); 2016 } 2017 *b = (PetscInt)(a); 2018 PetscFunctionReturn(0); 2019 } 2020 2021 /*@C 2022 PetscBLASIntCast - casts a `PetscInt` (which may be 64 bits in size) to a `PetscBLASInt` (which may be 32 bits in size), generates an 2023 error if the `PetscBLASInt` is not large enough to hold the number. 2024 2025 Not Collective 2026 2027 Input Parameter: 2028 . a - the `PetscInt` value 2029 2030 Output Parameter: 2031 . b - the resulting `PetscBLASInt` value 2032 2033 Level: advanced 2034 2035 Notes: 2036 Not available from Fortran 2037 2038 Errors if the integer is negative since PETSc calls to BLAS/LAPACK never need to cast negative integer inputs 2039 2040 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscMPIIntCast()`, `PetscIntCast()`, `PetscCountCast()` 2041 @*/ 2042 static inline PetscErrorCode PetscBLASIntCast(PetscInt a, PetscBLASInt *b) { 2043 PetscFunctionBegin; 2044 *b = 0; 2045 if (PetscDefined(USE_64BIT_INDICES) && !PetscDefined(HAVE_64BIT_BLAS_INDICES)) { 2046 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); 2047 } 2048 PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer to BLAS/LAPACK routine"); 2049 *b = (PetscBLASInt)(a); 2050 PetscFunctionReturn(0); 2051 } 2052 2053 /*@C 2054 PetscCuBLASIntCast - like `PetscBLASIntCast()`, but for `PetscCuBLASInt`. 2055 2056 Not Collective 2057 2058 Input Parameter: 2059 . a - the `PetscInt` value 2060 2061 Output Parameter: 2062 . b - the resulting `PetscCuBLASInt` value 2063 2064 Level: advanced 2065 2066 Notes: 2067 Errors if the integer is negative since PETSc calls to cuBLAS and friends never need to cast negative integer inputs 2068 2069 .seealso: `PetscCuBLASInt`, `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscMPIIntCast()`, `PetscIntCast()` 2070 @*/ 2071 static inline PetscErrorCode PetscCuBLASIntCast(PetscInt a, PetscCuBLASInt *b) { 2072 PetscFunctionBegin; 2073 *b = 0; 2074 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); 2075 PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Passing negative integer to cuBLAS routine"); 2076 *b = (PetscCuBLASInt)(a); 2077 PetscFunctionReturn(0); 2078 } 2079 2080 /*@C 2081 PetscMPIIntCast - casts a `PetscInt` (which may be 64 bits in size) to a `PetscMPIInt` (which may be 32 bits in size), generates an 2082 error if the `PetscMPIInt` is not large enough to hold the number. 2083 2084 Not Collective 2085 2086 Input Parameter: 2087 . a - the `PetscInt` value 2088 2089 Output Parameter: 2090 . b - the resulting `PetscMPIInt` value 2091 2092 Level: advanced 2093 2094 Not available from Fortran 2095 2096 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscIntCast()` 2097 @*/ 2098 static inline PetscErrorCode PetscMPIIntCast(PetscInt a, PetscMPIInt *b) { 2099 PetscFunctionBegin; 2100 *b = 0; 2101 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); 2102 *b = (PetscMPIInt)(a); 2103 PetscFunctionReturn(0); 2104 } 2105 2106 #define PetscInt64Mult(a, b) ((PetscInt64)(a)) * ((PetscInt64)(b)) 2107 2108 /*@C 2109 2110 PetscRealIntMultTruncate - Computes the product of a positive `PetscReal` and a positive `PetscInt` and truncates the value to slightly less than the maximal possible value 2111 2112 Not Collective 2113 2114 Input Parameters: 2115 + a - the `PetscReal` value 2116 - b - the second value 2117 2118 Returns: 2119 the result as a `PetscInt` value 2120 2121 Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64` 2122 Use `PetscIntMultTruncate()` to compute the product of two positive `PetscInt` and truncate to fit a `PetscInt` 2123 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` 2124 2125 Developers Note: 2126 We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but requires many more checks. 2127 2128 This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able. 2129 2130 Not available from Fortran 2131 2132 Level: advanced 2133 2134 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()` 2135 @*/ 2136 static inline PetscInt PetscRealIntMultTruncate(PetscReal a, PetscInt b) { 2137 PetscInt64 r = (PetscInt64)(a * (PetscReal)b); 2138 if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100; 2139 return (PetscInt)r; 2140 } 2141 2142 /*@C 2143 2144 PetscIntMultTruncate - Computes the product of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value 2145 2146 Not Collective 2147 2148 Input Parameters: 2149 + a - the PetscInt value 2150 - b - the second value 2151 2152 Returns: 2153 the result as a `PetscInt` value 2154 2155 Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64` 2156 Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt` 2157 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` 2158 2159 Not available from Fortran 2160 2161 Developers Note: 2162 We currently assume that `PetscInt` addition can never overflow, this is obviously wrong but requires many more checks. 2163 2164 This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able. 2165 2166 Level: advanced 2167 2168 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()`, `PetscIntSumError()` 2169 @*/ 2170 static inline PetscInt PetscIntMultTruncate(PetscInt a, PetscInt b) { 2171 PetscInt64 r = PetscInt64Mult(a, b); 2172 if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100; 2173 return (PetscInt)r; 2174 } 2175 2176 /*@C 2177 2178 PetscIntSumTruncate - Computes the sum of two positive `PetscInt` and truncates the value to slightly less than the maximal possible value 2179 2180 Not Collective 2181 2182 Input Parameters: 2183 + a - the `PetscInt` value 2184 - b - the second value 2185 2186 Returns: 2187 the result as a `PetscInt` value 2188 2189 Use `PetscInt64Mult()` to compute the product of two `PetscInt` as a `PetscInt64` 2190 Use `PetscRealIntMultTruncate()` to compute the product of a `PetscReal` and a `PetscInt` and truncate to fit a `PetscInt` 2191 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` 2192 2193 This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able. 2194 2195 Not available from Fortran 2196 2197 Level: advanced 2198 2199 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()` 2200 @*/ 2201 static inline PetscInt PetscIntSumTruncate(PetscInt a, PetscInt b) { 2202 PetscInt64 r = ((PetscInt64)a) + ((PetscInt64)b); 2203 if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100; 2204 return (PetscInt)r; 2205 } 2206 2207 /*@C 2208 2209 PetscIntMultError - Computes the product of two positive `PetscInt` and generates an error with overflow. 2210 2211 Not Collective 2212 2213 Input Parameters: 2214 + a - the `PetscInt` value 2215 - b - the second value 2216 2217 Output Parameter: 2218 . result - the result as a `PetscInt` value, or NULL if you do not want the result, you just want to check if it overflows 2219 2220 Use `PetscInt64Mult()` to compute the product of two `PetscInt` and store in a `PetscInt64` 2221 Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt` 2222 2223 Not available from Fortran 2224 2225 Developers Note: 2226 We currently assume that `PetscInt` addition does not overflow, this is obviously wrong but requires many more checks. 2227 2228 Level: advanced 2229 2230 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscIntMult64()`, `PetscIntSumError()` 2231 @*/ 2232 static inline PetscErrorCode PetscIntMultError(PetscInt a, PetscInt b, PetscInt *result) { 2233 PetscInt64 r; 2234 2235 PetscFunctionBegin; 2236 r = PetscInt64Mult(a, b); 2237 #if !defined(PETSC_USE_64BIT_INDICES) 2238 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); 2239 #endif 2240 if (result) *result = (PetscInt)r; 2241 PetscFunctionReturn(0); 2242 } 2243 2244 /*@C 2245 2246 PetscIntSumError - Computes the sum of two positive `PetscInt` and generates an error with overflow. 2247 2248 Not Collective 2249 2250 Input Parameters: 2251 + a - the `PetscInt` value 2252 - b - the second value 2253 2254 Output Parameter: 2255 . c - the result as a `PetscInt` value, or NULL if you do not want the result, you just want to check if it overflows 2256 2257 Use `PetscInt64Mult()` to compute the product of two 32 bit PetscInt and store in a `PetscInt64` 2258 Use `PetscIntMultTruncate()` to compute the product of two `PetscInt` and truncate it to fit in a `PetscInt` 2259 2260 Not available from Fortran 2261 2262 Level: advanced 2263 2264 .seealso: `PetscBLASInt`, `PetscMPIInt`, `PetscInt`, `PetscBLASIntCast()`, `PetscInt64Mult()`, `PetscIntMultError()` 2265 @*/ 2266 static inline PetscErrorCode PetscIntSumError(PetscInt a, PetscInt b, PetscInt *result) { 2267 PetscInt64 r; 2268 2269 PetscFunctionBegin; 2270 r = ((PetscInt64)a) + ((PetscInt64)b); 2271 #if !defined(PETSC_USE_64BIT_INDICES) 2272 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); 2273 #endif 2274 if (result) *result = (PetscInt)r; 2275 PetscFunctionReturn(0); 2276 } 2277 2278 /* 2279 The IBM include files define hz, here we hide it so that it may be used as a regular user variable. 2280 */ 2281 #if defined(hz) 2282 #undef hz 2283 #endif 2284 2285 #include <limits.h> 2286 2287 /* The number of bits in a byte */ 2288 2289 #define PETSC_BITS_PER_BYTE CHAR_BIT 2290 2291 #if defined(PETSC_HAVE_SYS_TYPES_H) 2292 #include <sys/types.h> 2293 #endif 2294 2295 /*MC 2296 2297 PETSC_VERSION - This manual page provides information about how PETSc documents and uses its version information. This information is available to both C/C++ 2298 and Fortran compilers when petscsys.h is included. 2299 2300 The current PETSc version and the API for accessing it are defined in petscversion.h 2301 2302 The complete version number is given as the triple PETSC_VERSION_MAJOR.PETSC_VERSION_MINOR.PETSC_VERSION_SUBMINOR (in short hand x.y.z) 2303 2304 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 2305 where only a change in the major version number (x) indicates a change in the API. 2306 2307 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 2308 2309 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), 2310 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 2311 version number (x.y.z). 2312 2313 `PETSC_RELEASE_DATE` is the date the x.y version was released (i.e. the version before any patch releases) 2314 2315 `PETSC_VERSION_DATE` is the date the x.y.z version was released 2316 2317 `PETSC_VERSION_GIT` is the last git commit to the repository given in the form vx.y.z-wwwww 2318 2319 `PETSC_VERSION_DATE_GIT` is the date of the last git commit to the repository 2320 2321 `PETSC_VERSION_()` is deprecated and will eventually be removed. 2322 2323 Level: intermediate 2324 2325 M*/ 2326 2327 PETSC_EXTERN PetscErrorCode PetscGetArchType(char[], size_t); 2328 PETSC_EXTERN PetscErrorCode PetscGetHostName(char[], size_t); 2329 PETSC_EXTERN PetscErrorCode PetscGetUserName(char[], size_t); 2330 PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[], size_t); 2331 PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]); 2332 PETSC_EXTERN PetscErrorCode PetscGetDate(char[], size_t); 2333 PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t); 2334 PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt *, PetscInt *, PetscInt *, PetscInt *); 2335 2336 PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt, const PetscInt[], PetscBool *); 2337 PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt, const PetscMPIInt[], PetscBool *); 2338 PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt, const PetscReal[], PetscBool *); 2339 PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt, PetscInt[]); 2340 PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt, PetscInt[]); 2341 PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt *, PetscInt[]); 2342 PETSC_EXTERN PetscErrorCode PetscSortedCheckDupsInt(PetscInt, const PetscInt[], PetscBool *); 2343 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt *, PetscInt[]); 2344 PETSC_EXTERN PetscErrorCode PetscCheckDupsInt(PetscInt, const PetscInt[], PetscBool *); 2345 PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt *); 2346 PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt *); 2347 PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt, const PetscInt[], PetscInt[]); 2348 PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt, const char *[], PetscInt[]); 2349 PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt, PetscInt[], PetscInt[]); 2350 PETSC_EXTERN PetscErrorCode PetscSortIntWithCountArray(PetscCount, PetscInt[], PetscCount[]); 2351 PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt, PetscInt[], PetscInt[], PetscInt[]); 2352 PETSC_EXTERN PetscErrorCode PetscSortIntWithIntCountArrayPair(PetscCount, PetscInt[], PetscInt[], PetscCount[]); 2353 PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt, PetscMPIInt[]); 2354 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt *, PetscMPIInt[]); 2355 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt, PetscMPIInt[], PetscMPIInt[]); 2356 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt, PetscMPIInt[], PetscInt[]); 2357 PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt, PetscInt[], PetscScalar[]); 2358 PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt, PetscInt[], void *, size_t, void *); 2359 PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt, PetscReal[]); 2360 PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt, PetscReal[], PetscInt[]); 2361 PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt, const PetscReal[], PetscInt[]); 2362 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt *, PetscReal[]); 2363 PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal, PetscInt, const PetscReal[], PetscReal, PetscInt *); 2364 PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt, PetscInt, PetscScalar[], PetscInt[]); 2365 PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt, PetscInt, PetscReal[], PetscInt[]); 2366 PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt, const PetscBool[], const PetscInt[], PetscInt *, PetscInt **, PetscInt **, PetscInt **, PetscInt **); 2367 PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt, const PetscInt[], const PetscInt[], PetscInt, const PetscInt[], const PetscInt[], PetscInt *, PetscInt **, PetscInt **); 2368 PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscInt *, PetscInt **); 2369 PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt, const PetscMPIInt[], PetscInt, const PetscMPIInt[], PetscInt *, PetscMPIInt **); 2370 PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *); 2371 2372 PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt, void *, size_t, int (*)(const void *, const void *, void *), void *); 2373 PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt, PetscInt[]); 2374 PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt, PetscMPIInt[]); 2375 PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt, PetscReal[]); 2376 PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt, void *, size_t, void *, size_t, int (*)(const void *, const void *, void *), void *); 2377 PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt, PetscInt[], PetscInt[]); 2378 PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt, PetscMPIInt[], PetscMPIInt[]); 2379 PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt, PetscReal[], PetscInt[]); 2380 2381 PETSC_EXTERN PetscErrorCode PetscSetDisplay(void); 2382 PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[], size_t); 2383 2384 /*J 2385 PetscRandomType - String with the name of a PETSc randomizer 2386 2387 Level: beginner 2388 2389 Notes: 2390 To use `PETSCSPRNG` or `PETSCRANDOM123` you must have ./configure PETSc 2391 with the option --download-sprng or --download-random123 2392 2393 .seealso: `PetscRandomSetType()`, `PetscRandom`, `PetscRandomCreate()` 2394 J*/ 2395 typedef const char *PetscRandomType; 2396 #define PETSCRAND "rand" 2397 #define PETSCRAND48 "rand48" 2398 #define PETSCSPRNG "sprng" 2399 #define PETSCRANDER48 "rander48" 2400 #define PETSCRANDOM123 "random123" 2401 #define PETSCCURAND "curand" 2402 2403 /* Logging support */ 2404 PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID; 2405 2406 PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void); 2407 2408 /* Dynamic creation and loading functions */ 2409 PETSC_EXTERN PetscFunctionList PetscRandomList; 2410 2411 PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[], PetscErrorCode (*)(PetscRandom)); 2412 PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType); 2413 PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom); 2414 PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType *); 2415 PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom, PetscObject, const char[]); 2416 PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom, PetscViewer); 2417 2418 PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm, PetscRandom *); 2419 PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom, PetscScalar *); 2420 PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom, PetscReal *); 2421 PETSC_EXTERN PetscErrorCode PetscRandomGetValues(PetscRandom, PetscInt, PetscScalar *); 2422 PETSC_EXTERN PetscErrorCode PetscRandomGetValuesReal(PetscRandom, PetscInt, PetscReal *); 2423 PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom, PetscScalar *, PetscScalar *); 2424 PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom, PetscScalar, PetscScalar); 2425 PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom, unsigned long); 2426 PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom, unsigned long *); 2427 PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom); 2428 PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom *); 2429 2430 PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[], char[], size_t); 2431 PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[], char[], size_t); 2432 PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[], size_t); 2433 PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[], char[]); 2434 PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[], size_t); 2435 PETSC_EXTERN PetscErrorCode PetscTestFile(const char[], char, PetscBool *); 2436 PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[], char, PetscBool *); 2437 PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]); 2438 PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]); 2439 PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]); 2440 2441 static inline PetscBool PetscBinaryBigEndian(void) { 2442 long _petsc_v = 1; 2443 return ((char *)&_petsc_v)[0] ? PETSC_FALSE : PETSC_TRUE; 2444 } 2445 2446 PETSC_EXTERN PetscErrorCode PetscBinaryRead(int, void *, PetscInt, PetscInt *, PetscDataType); 2447 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm, int, void *, PetscInt, PetscInt *, PetscDataType); 2448 PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int, const void *, PetscInt, PetscDataType); 2449 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm, int, const void *, PetscInt, PetscDataType); 2450 PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[], PetscFileMode, int *); 2451 PETSC_EXTERN PetscErrorCode PetscBinaryClose(int); 2452 PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm, PetscBool *); 2453 PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm, PetscBool *); 2454 PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm, char[], size_t); 2455 PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm, const char[], char[], size_t, PetscBool *); 2456 PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm, const char[], char[], size_t, PetscBool *); 2457 #if defined(PETSC_USE_SOCKET_VIEWER) 2458 PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[], int, int *); 2459 #endif 2460 2461 PETSC_EXTERN PetscErrorCode PetscBinarySeek(int, off_t, PetscBinarySeekType, off_t *); 2462 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm, int, off_t, PetscBinarySeekType, off_t *); 2463 PETSC_EXTERN PetscErrorCode PetscByteSwap(void *, PetscDataType, PetscInt); 2464 2465 PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]); 2466 PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[], PetscBool); 2467 PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void); 2468 PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char *); 2469 PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void); 2470 PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void); 2471 PETSC_EXTERN PetscErrorCode PetscWaitOnError(void); 2472 2473 PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt *); 2474 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], PetscMPIInt **, PetscMPIInt **); 2475 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscMPIInt **, PetscMPIInt **, PetscMPIInt **); 2476 PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscInt ***, MPI_Request **); 2477 PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscMPIInt[], const PetscMPIInt[], PetscScalar ***, MPI_Request **); 2478 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); 2479 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); 2480 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); 2481 2482 PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[]; 2483 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm, PetscBuildTwoSidedType); 2484 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm, PetscBuildTwoSidedType *); 2485 2486 PETSC_EXTERN PetscErrorCode PetscSSEIsEnabled(MPI_Comm, PetscBool *, PetscBool *); 2487 2488 PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject); 2489 2490 PETSC_EXTERN const char *const PetscSubcommTypes[]; 2491 2492 struct _n_PetscSubcomm { 2493 MPI_Comm parent; /* parent communicator */ 2494 MPI_Comm dupparent; /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */ 2495 MPI_Comm child; /* the sub-communicator */ 2496 PetscMPIInt n; /* num of subcommunicators under the parent communicator */ 2497 PetscMPIInt color; /* color of processors belong to this communicator */ 2498 PetscMPIInt *subsize; /* size of subcommunicator[color] */ 2499 PetscSubcommType type; 2500 char *subcommprefix; 2501 }; 2502 2503 static inline MPI_Comm PetscSubcommParent(PetscSubcomm scomm) { 2504 return scomm->parent; 2505 } 2506 static inline MPI_Comm PetscSubcommChild(PetscSubcomm scomm) { 2507 return scomm->child; 2508 } 2509 static inline MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm) { 2510 return scomm->dupparent; 2511 } 2512 PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm, PetscSubcomm *); 2513 PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm *); 2514 PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm, PetscInt); 2515 PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm, PetscSubcommType); 2516 PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm, PetscMPIInt, PetscMPIInt); 2517 PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm, PetscViewer); 2518 PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm); 2519 PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm, const char[]); 2520 PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm, MPI_Comm *); 2521 PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm, MPI_Comm *); 2522 PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm, MPI_Comm *); 2523 2524 PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt, PetscHeap *); 2525 PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap, PetscInt, PetscInt); 2526 PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap, PetscInt *, PetscInt *); 2527 PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap, PetscInt *, PetscInt *); 2528 PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap, PetscInt, PetscInt); 2529 PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap); 2530 PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap *); 2531 PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap, PetscViewer); 2532 2533 PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer); 2534 PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm, PetscShmComm *); 2535 PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm, PetscMPIInt, PetscMPIInt *); 2536 PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm, PetscMPIInt, PetscMPIInt *); 2537 PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm, MPI_Comm *); 2538 2539 /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */ 2540 PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm, PetscInt, PetscOmpCtrl *); 2541 PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl, MPI_Comm *, MPI_Comm *, PetscBool *); 2542 PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl *); 2543 PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl); 2544 PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl); 2545 PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl); 2546 2547 PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t, size_t, PetscSegBuffer *); 2548 PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer *); 2549 PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer, size_t, void *); 2550 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer, void *); 2551 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer, void *); 2552 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer, void *); 2553 PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer, size_t *); 2554 PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer, size_t); 2555 2556 /* Type-safe wrapper to encourage use of PETSC_RESTRICT. Does not use PetscFunctionBegin because the error handling 2557 * prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as 2558 * possible. */ 2559 static inline PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg, size_t count, PetscInt *PETSC_RESTRICT *slot) { 2560 return PetscSegBufferGet(seg, count, (void **)slot); 2561 } 2562 2563 extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton; 2564 PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted *); 2565 PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted *); 2566 PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted, const char *, const char *, PetscBool *); 2567 2568 #include <stdarg.h> 2569 PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char *, size_t, const char[], size_t *, va_list); 2570 PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE *, const char[], va_list); 2571 2572 PETSC_EXTERN PetscSegBuffer PetscCitationsList; 2573 2574 /*@C 2575 PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code. 2576 2577 Not Collective - only what is registered on rank 0 of `PETSC_COMM_WORLD` will be printed 2578 2579 Input Parameters: 2580 + cite - the bibtex item, formated to displayed on multiple lines nicely 2581 - set - a boolean variable initially set to `PETSC_FALSE`; this is used to insure only a single registration of the citation 2582 2583 Level: intermediate 2584 2585 Not available from Fortran 2586 2587 Options Database: 2588 . -citations [filename] - print out the bibtex entries for the given computation 2589 @*/ 2590 static inline PetscErrorCode PetscCitationsRegister(const char cit[], PetscBool *set) { 2591 size_t len; 2592 char *vstring; 2593 2594 PetscFunctionBegin; 2595 if (set && *set) PetscFunctionReturn(0); 2596 PetscCall(PetscStrlen(cit, &len)); 2597 PetscCall(PetscSegBufferGet(PetscCitationsList, len, &vstring)); 2598 PetscCall(PetscArraycpy(vstring, cit, len)); 2599 if (set) *set = PETSC_TRUE; 2600 PetscFunctionReturn(0); 2601 } 2602 2603 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Google has discontinued its URL shortener service") PetscErrorCode PetscURLShorten(const char[], char[], size_t c); 2604 PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm, char[], char[], size_t); 2605 PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm, const char[], char[], size_t); 2606 PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm, const char[], const char[]); 2607 2608 PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm, char[], char[], size_t); 2609 PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm, const char[], char[], char[], size_t); 2610 2611 PETSC_EXTERN PetscErrorCode PetscGlobusGetTransfers(MPI_Comm, const char[], char[], size_t); 2612 2613 PETSC_EXTERN PetscErrorCode PetscTextBelt(MPI_Comm, const char[], const char[], PetscBool *); 2614 PETSC_EXTERN PetscErrorCode PetscTellMyCell(MPI_Comm, const char[], const char[], PetscBool *); 2615 2616 PETSC_EXTERN PetscErrorCode PetscPullJSONValue(const char[], const char[], char[], size_t, PetscBool *); 2617 PETSC_EXTERN PetscErrorCode PetscPushJSONValue(char[], const char[], const char[], size_t); 2618 2619 #if defined(PETSC_USE_DEBUG) 2620 static inline unsigned int PetscStrHash(const char *str) { 2621 unsigned int c, hash = 5381; 2622 2623 while ((c = (unsigned int)*str++)) hash = ((hash << 5) + hash) + c; /* hash * 33 + c */ 2624 return hash; 2625 } 2626 2627 /*MC 2628 MPIU_Allreduce - a PETSc replacement for `MPI_Allreduce()` that tries to determine if the call from all the MPI ranks occur from the 2629 same place in the PETSc code. This helps to detect bugs where different MPI ranks follow different code paths 2630 resulting in inconsistent and incorrect calls to `MPI_Allreduce()`. 2631 2632 Collective 2633 2634 Synopsis: 2635 #include <petscsys.h> 2636 PetscErrorCode MPIU_Allreduce(void *indata,void *outdata,PetscMPIInt count,MPI_Datatype datatype, MPI_Op op, MPI_Comm comm); 2637 2638 Input Parameters: 2639 + a - pointer to the input data to be reduced 2640 . c - the number of MPI data items in a and b 2641 . d - the MPI datatype, for example `MPI_INT` 2642 . e - the MPI operation, for example `MPI_SUM` 2643 - fcomm - the MPI communicator on which the operation occurs 2644 2645 Output Parameter: 2646 . b - the reduced values 2647 2648 Notes: 2649 In optimized mode this directly calls `MPI_Allreduce()` 2650 2651 This is defined as a macro that can return error codes internally so it cannot be used in a subroutine that returns void. 2652 2653 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 2654 2655 Level: developer 2656 2657 .seealso: `MPI_Allreduce()` 2658 M*/ 2659 #define MPIU_Allreduce(a, b, c, d, e, fcomm) \ 2660 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)));) 2661 #else 2662 #define MPIU_Allreduce(a, b, c, d, e, fcomm) PetscMacroReturnStandard(PetscCallMPI(MPI_Allreduce((a), (b), (c), d, e, (fcomm)))) 2663 #endif 2664 2665 #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) 2666 PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint, PetscMPIInt, MPI_Info, MPI_Comm, void *, MPI_Win *); 2667 PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win, PetscMPIInt, MPI_Aint *, PetscMPIInt *, void *); 2668 #endif 2669 2670 /* this is a vile hack */ 2671 #if defined(PETSC_HAVE_NECMPI) 2672 #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) 2673 #define MPI_Type_free(a) (*(a) = MPI_DATATYPE_NULL, 0); 2674 #endif 2675 #endif 2676 2677 /* 2678 List of external packages and queries on it 2679 */ 2680 PETSC_EXTERN PetscErrorCode PetscHasExternalPackage(const char[], PetscBool *); 2681 2682 /* 2683 OpenMP support 2684 */ 2685 #if defined(_OPENMP) 2686 #define PetscPragmaOMP(...) _Pragma(PetscStringize(omp __VA_ARGS__)) 2687 #else // no OpenMP so no threads 2688 #define PetscPragmaOMP(...) 2689 #endif 2690 2691 /* this cannot go here because it may be in a different shared library */ 2692 PETSC_EXTERN PetscErrorCode PCMPIServerBegin(void); 2693 PETSC_EXTERN PetscErrorCode PCMPIServerEnd(void); 2694 PETSC_EXTERN PetscErrorCode PCMPICommsDestroy(void); 2695 #endif 2696