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