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