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