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