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