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