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