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