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