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