1 /* 2 This is the main PETSc include file (for C and C++). It is included by all 3 other PETSc include files, so it almost never has to be specifically included. 4 */ 5 #if !defined(PETSCSYS_H) 6 #define PETSCSYS_H 7 /* ========================================================================== */ 8 /* 9 petscconf.h is contained in ${PETSC_ARCH}/include/petscconf.h it is 10 found automatically by the compiler due to the -I${PETSC_DIR}/${PETSC_ARCH}/include. 11 For --prefix installs the ${PETSC_ARCH}/ does not exist and petscconf.h is in the same 12 directory as the other PETSc include files. 13 */ 14 #include <petscconf.h> 15 #include <petscfix.h> 16 17 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP) || defined(PETSC_HAVE_KOKKOS) 18 #define PETSC_HAVE_DEVICE 19 #endif 20 21 #if defined(PETSC_DESIRE_FEATURE_TEST_MACROS) 22 /* 23 Feature test macros must be included before headers defined by IEEE Std 1003.1-2001 24 We only turn these in PETSc source files that require them by setting PETSC_DESIRE_FEATURE_TEST_MACROS 25 */ 26 # if defined(PETSC__POSIX_C_SOURCE_200112L) && !defined(_POSIX_C_SOURCE) 27 # define _POSIX_C_SOURCE 200112L 28 # endif 29 # if defined(PETSC__BSD_SOURCE) && !defined(_BSD_SOURCE) 30 # define _BSD_SOURCE 31 # endif 32 # if defined(PETSC__DEFAULT_SOURCE) && !defined(_DEFAULT_SOURCE) 33 # define _DEFAULT_SOURCE 34 # endif 35 # if defined(PETSC__GNU_SOURCE) && !defined(_GNU_SOURCE) 36 # define _GNU_SOURCE 37 # endif 38 #endif 39 40 #include <petscsystypes.h> 41 42 /* ========================================================================== */ 43 /* 44 This facilitates using the C version of PETSc from C++ and the C++ version from C. 45 */ 46 #if defined(__cplusplus) 47 # define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_CXX 48 #else 49 # define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_C 50 #endif 51 52 /* ========================================================================== */ 53 /* 54 Since PETSc manages its own extern "C" handling users should never include PETSc include 55 files within extern "C". This will generate a compiler error if a user does put the include 56 file within an extern "C". 57 */ 58 #if defined(__cplusplus) 59 void assert_never_put_petsc_headers_inside_an_extern_c(int); void assert_never_put_petsc_headers_inside_an_extern_c(double); 60 #endif 61 62 #if defined(__cplusplus) 63 # define PETSC_RESTRICT PETSC_CXX_RESTRICT 64 #else 65 # define PETSC_RESTRICT PETSC_C_RESTRICT 66 #endif 67 68 #if defined(__cplusplus) 69 # define PETSC_INLINE PETSC_CXX_INLINE 70 #else 71 # define PETSC_INLINE PETSC_C_INLINE 72 #endif 73 74 #define PETSC_STATIC_INLINE static PETSC_INLINE 75 76 #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES) /* For Win32 shared libraries */ 77 # define PETSC_DLLEXPORT __declspec(dllexport) 78 # define PETSC_DLLIMPORT __declspec(dllimport) 79 # define PETSC_VISIBILITY_INTERNAL 80 #elif defined(PETSC_USE_VISIBILITY_CXX) && defined(__cplusplus) 81 # define PETSC_DLLEXPORT __attribute__((visibility ("default"))) 82 # define PETSC_DLLIMPORT __attribute__((visibility ("default"))) 83 # define PETSC_VISIBILITY_INTERNAL __attribute__((visibility ("hidden"))) 84 #elif defined(PETSC_USE_VISIBILITY_C) && !defined(__cplusplus) 85 # define PETSC_DLLEXPORT __attribute__((visibility ("default"))) 86 # define PETSC_DLLIMPORT __attribute__((visibility ("default"))) 87 # define PETSC_VISIBILITY_INTERNAL __attribute__((visibility ("hidden"))) 88 #else 89 # define PETSC_DLLEXPORT 90 # define PETSC_DLLIMPORT 91 # define PETSC_VISIBILITY_INTERNAL 92 #endif 93 94 #if defined(petsc_EXPORTS) /* CMake defines this when building the shared library */ 95 # define PETSC_VISIBILITY_PUBLIC PETSC_DLLEXPORT 96 #else /* Win32 users need this to import symbols from petsc.dll */ 97 # define PETSC_VISIBILITY_PUBLIC PETSC_DLLIMPORT 98 #endif 99 100 /* 101 Functions tagged with PETSC_EXTERN in the header files are 102 always defined as extern "C" when compiled with C++ so they may be 103 used from C and are always visible in the shared libraries 104 */ 105 #if defined(__cplusplus) 106 # define PETSC_EXTERN extern "C" PETSC_VISIBILITY_PUBLIC 107 # define PETSC_EXTERN_TYPEDEF extern "C" 108 # define PETSC_INTERN extern "C" PETSC_VISIBILITY_INTERNAL 109 #else 110 # define PETSC_EXTERN extern PETSC_VISIBILITY_PUBLIC 111 # define PETSC_EXTERN_TYPEDEF 112 # define PETSC_INTERN extern PETSC_VISIBILITY_INTERNAL 113 #endif 114 115 #include <petscversion.h> 116 #define PETSC_AUTHOR_INFO " The PETSc Team\n petsc-maint@mcs.anl.gov\n https://www.mcs.anl.gov/petsc/\n" 117 118 /* ========================================================================== */ 119 120 /* 121 Defines the interface to MPI allowing the use of all MPI functions. 122 123 PETSc does not use the C++ binding of MPI at ALL. The following flag 124 makes sure the C++ bindings are not included. The C++ bindings REQUIRE 125 putting mpi.h before ANY C++ include files, we cannot control this 126 with all PETSc users. Users who want to use the MPI C++ bindings can include 127 mpicxx.h directly in their code 128 */ 129 #if !defined(MPICH_SKIP_MPICXX) 130 # define MPICH_SKIP_MPICXX 1 131 #endif 132 #if !defined(OMPI_SKIP_MPICXX) 133 # define OMPI_SKIP_MPICXX 1 134 #endif 135 #if defined(PETSC_HAVE_MPIUNI) 136 # include <petsc/mpiuni/mpi.h> 137 #else 138 # include <mpi.h> 139 #endif 140 141 /*MC 142 PetscDefined - determine whether a boolean macro is defined 143 144 Notes: 145 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 /* 380 Declare extern C stuff after including external header files 381 */ 382 383 PETSC_EXTERN const char *const PetscBools[]; 384 385 /* 386 Defines elementary mathematics functions and constants. 387 */ 388 #include <petscmath.h> 389 390 PETSC_EXTERN const char *const PetscCopyModes[]; 391 392 /*MC 393 PETSC_IGNORE - same as NULL, means PETSc will ignore this argument 394 395 Level: beginner 396 397 Note: 398 Accepted by many PETSc functions to not set a parameter and instead use 399 some default 400 401 Fortran Notes: 402 This macro does not exist in Fortran; you must use PETSC_NULL_INTEGER, 403 PETSC_NULL_DOUBLE_PRECISION etc 404 405 .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_DETERMINE 406 407 M*/ 408 #define PETSC_IGNORE NULL 409 410 /* This is deprecated */ 411 #define PETSC_NULL NULL 412 413 /*MC 414 PETSC_DECIDE - standard way of passing in integer or floating point parameter 415 where you wish PETSc to use the default. 416 417 Level: beginner 418 419 .seealso: PETSC_DEFAULT, PETSC_IGNORE, PETSC_DETERMINE 420 421 M*/ 422 #define PETSC_DECIDE -1 423 424 /*MC 425 PETSC_DETERMINE - standard way of passing in integer or floating point parameter 426 where you wish PETSc to compute the required value. 427 428 Level: beginner 429 430 Developer Note: 431 I would like to use const PetscInt PETSC_DETERMINE = PETSC_DECIDE; but for 432 some reason this is not allowed by the standard even though PETSC_DECIDE is a constant value. 433 434 .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_IGNORE, VecSetSizes() 435 436 M*/ 437 #define PETSC_DETERMINE PETSC_DECIDE 438 439 /*MC 440 PETSC_DEFAULT - standard way of passing in integer or floating point parameter 441 where you wish PETSc to use the default. 442 443 Level: beginner 444 445 Fortran Notes: 446 You need to use PETSC_DEFAULT_INTEGER or PETSC_DEFAULT_REAL. 447 448 .seealso: PETSC_DECIDE, PETSC_IGNORE, PETSC_DETERMINE 449 450 M*/ 451 #define PETSC_DEFAULT -2 452 453 /*MC 454 PETSC_COMM_WORLD - the equivalent of the MPI_COMM_WORLD communicator which represents 455 all the processes that PETSc knows about. 456 457 Level: beginner 458 459 Notes: 460 By default PETSC_COMM_WORLD and MPI_COMM_WORLD are identical unless you wish to 461 run PETSc on ONLY a subset of MPI_COMM_WORLD. In that case create your new (smaller) 462 communicator, call it, say comm, and set PETSC_COMM_WORLD = comm BEFORE calling 463 PetscInitialize(), but after MPI_Init() has been called. 464 465 The value of PETSC_COMM_WORLD should never be USED/accessed before PetscInitialize() 466 is called because it may not have a valid value yet. 467 468 .seealso: PETSC_COMM_SELF 469 470 M*/ 471 PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD; 472 473 /*MC 474 PETSC_COMM_SELF - This is always MPI_COMM_SELF 475 476 Level: beginner 477 478 Notes: 479 Do not USE/access or set this variable before PetscInitialize() has been called. 480 481 .seealso: PETSC_COMM_WORLD 482 483 M*/ 484 #define PETSC_COMM_SELF MPI_COMM_SELF 485 486 /*MC 487 PETSC_MPI_THREAD_REQUIRED - the required threading support used if PETSc initializes 488 MPI with MPI_Init_thread. 489 490 Level: beginner 491 492 Notes: 493 By default PETSC_MPI_THREAD_REQUIRED equals MPI_THREAD_FUNNELED. 494 495 .seealso: PetscInitialize() 496 497 M*/ 498 PETSC_EXTERN PetscMPIInt PETSC_MPI_THREAD_REQUIRED; 499 500 PETSC_EXTERN PetscBool PetscBeganMPI; 501 PETSC_EXTERN PetscBool PetscErrorHandlingInitialized; 502 PETSC_EXTERN PetscBool PetscInitializeCalled; 503 PETSC_EXTERN PetscBool PetscFinalizeCalled; 504 PETSC_EXTERN PetscBool PetscViennaCLSynchronize; 505 506 PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm)); 507 PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*); 508 PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm*); 509 510 #if defined(PETSC_HAVE_CUDA) 511 PETSC_EXTERN PetscBool PetscCUDASynchronize; 512 PETSC_EXTERN PetscErrorCode PetscCUDAInitialize(MPI_Comm,PetscInt); 513 PETSC_EXTERN PetscErrorCode PetscCUDAInitializeCheck(void); 514 #endif 515 516 #if defined(PETSC_HAVE_HIP) 517 PETSC_EXTERN PetscBool PetscHIPSynchronize; 518 PETSC_EXTERN PetscErrorCode PetscHIPInitialize(MPI_Comm,PetscInt); 519 PETSC_EXTERN PetscErrorCode PetscHIPInitializeCheck(void); 520 #endif 521 522 #if defined(PETSC_HAVE_ELEMENTAL) 523 PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void); 524 PETSC_EXTERN PetscErrorCode PetscElementalInitialized(PetscBool*); 525 PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void); 526 #endif 527 528 /*MC 529 PetscMalloc - Allocates memory, One should use PetscNew(), PetscMalloc1() or PetscCalloc1() usually instead of this 530 531 Synopsis: 532 #include <petscsys.h> 533 PetscErrorCode PetscMalloc(size_t m,void **result) 534 535 Not Collective 536 537 Input Parameter: 538 . m - number of bytes to allocate 539 540 Output Parameter: 541 . result - memory allocated 542 543 Level: beginner 544 545 Notes: 546 Memory is always allocated at least double aligned 547 548 It is safe to allocate size 0 and pass the resulting pointer (which may or may not be NULL) to PetscFree(). 549 550 .seealso: PetscFree(), PetscNew() 551 552 M*/ 553 #define PetscMalloc(a,b) ((*PetscTrMalloc)((a),PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(void**)(b))) 554 555 /*MC 556 PetscRealloc - Rellocates memory 557 558 Synopsis: 559 #include <petscsys.h> 560 PetscErrorCode PetscRealloc(size_t m,void **result) 561 562 Not Collective 563 564 Input Parameters: 565 + m - number of bytes to allocate 566 - result - previous memory 567 568 Output Parameter: 569 . result - new memory allocated 570 571 Level: developer 572 573 Notes: 574 Memory is always allocated at least double aligned 575 576 .seealso: PetscMalloc(), PetscFree(), PetscNew() 577 578 M*/ 579 #define PetscRealloc(a,b) ((*PetscTrRealloc)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__,(void**)(b))) 580 581 /*MC 582 PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment 583 584 Synopsis: 585 #include <petscsys.h> 586 void *PetscAddrAlign(void *addr) 587 588 Not Collective 589 590 Input Parameters: 591 . addr - address to align (any pointer type) 592 593 Level: developer 594 595 .seealso: PetscMallocAlign() 596 597 M*/ 598 #define PetscAddrAlign(a) (void*)((((PETSC_UINTPTR_T)(a))+(PETSC_MEMALIGN-1)) & ~(PETSC_MEMALIGN-1)) 599 600 /*MC 601 PetscMalloc1 - Allocates an array of memory aligned to PETSC_MEMALIGN 602 603 Synopsis: 604 #include <petscsys.h> 605 PetscErrorCode PetscMalloc1(size_t m1,type **r1) 606 607 Not Collective 608 609 Input Parameter: 610 . m1 - number of elements to allocate (may be zero) 611 612 Output Parameter: 613 . r1 - memory allocated 614 615 Note: 616 This uses the sizeof() of the memory type requested to determine the total memory to be allocated, therefore you should not 617 multiply the number of elements requested by the sizeof() the type. For example use 618 $ PetscInt *id; 619 $ PetscMalloc1(10,&id); 620 not 621 $ PetscInt *id; 622 $ PetscMalloc1(10*sizeof(PetscInt),&id); 623 624 Does not zero the memory allocated, use PetscCalloc1() to obtain memory that has been zeroed. 625 626 Level: beginner 627 628 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc1(), PetscMalloc2() 629 630 M*/ 631 #define PetscMalloc1(m1,r1) PetscMallocA(1,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1)) 632 633 /*MC 634 PetscCalloc1 - Allocates a cleared (zeroed) array of memory aligned to PETSC_MEMALIGN 635 636 Synopsis: 637 #include <petscsys.h> 638 PetscErrorCode PetscCalloc1(size_t m1,type **r1) 639 640 Not Collective 641 642 Input Parameter: 643 . m1 - number of elements to allocate in 1st chunk (may be zero) 644 645 Output Parameter: 646 . r1 - memory allocated 647 648 Notes: 649 See PetsMalloc1() for more details on usage. 650 651 Level: beginner 652 653 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc1(), PetscCalloc2() 654 655 M*/ 656 #define PetscCalloc1(m1,r1) PetscMallocA(1,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1)) 657 658 /*MC 659 PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN 660 661 Synopsis: 662 #include <petscsys.h> 663 PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2) 664 665 Not Collective 666 667 Input Parameter: 668 + m1 - number of elements to allocate in 1st chunk (may be zero) 669 - m2 - number of elements to allocate in 2nd chunk (may be zero) 670 671 Output Parameter: 672 + r1 - memory allocated in first chunk 673 - r2 - memory allocated in second chunk 674 675 Level: developer 676 677 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc1(), PetscCalloc2() 678 679 M*/ 680 #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)) 681 682 /*MC 683 PetscCalloc2 - Allocates 2 cleared (zeroed) arrays of memory both aligned to PETSC_MEMALIGN 684 685 Synopsis: 686 #include <petscsys.h> 687 PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2) 688 689 Not Collective 690 691 Input Parameter: 692 + m1 - number of elements to allocate in 1st chunk (may be zero) 693 - m2 - number of elements to allocate in 2nd chunk (may be zero) 694 695 Output Parameter: 696 + r1 - memory allocated in first chunk 697 - r2 - memory allocated in second chunk 698 699 Level: developer 700 701 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc1(), PetscMalloc2() 702 703 M*/ 704 #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)) 705 706 /*MC 707 PetscMalloc3 - Allocates 3 arrays of memory, all aligned to PETSC_MEMALIGN 708 709 Synopsis: 710 #include <petscsys.h> 711 PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3) 712 713 Not Collective 714 715 Input Parameter: 716 + m1 - number of elements to allocate in 1st chunk (may be zero) 717 . m2 - number of elements to allocate in 2nd chunk (may be zero) 718 - m3 - number of elements to allocate in 3rd chunk (may be zero) 719 720 Output Parameter: 721 + r1 - memory allocated in first chunk 722 . r2 - memory allocated in second chunk 723 - r3 - memory allocated in third chunk 724 725 Level: developer 726 727 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc3(), PetscFree3() 728 729 M*/ 730 #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)) 731 732 /*MC 733 PetscCalloc3 - Allocates 3 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN 734 735 Synopsis: 736 #include <petscsys.h> 737 PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3) 738 739 Not Collective 740 741 Input Parameter: 742 + m1 - number of elements to allocate in 1st chunk (may be zero) 743 . m2 - number of elements to allocate in 2nd chunk (may be zero) 744 - m3 - number of elements to allocate in 3rd chunk (may be zero) 745 746 Output Parameter: 747 + r1 - memory allocated in first chunk 748 . r2 - memory allocated in second chunk 749 - r3 - memory allocated in third chunk 750 751 Level: developer 752 753 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc2(), PetscMalloc3(), PetscFree3() 754 755 M*/ 756 #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)) 757 758 /*MC 759 PetscMalloc4 - Allocates 4 arrays of memory, all aligned to PETSC_MEMALIGN 760 761 Synopsis: 762 #include <petscsys.h> 763 PetscErrorCode PetscMalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4) 764 765 Not Collective 766 767 Input Parameter: 768 + m1 - number of elements to allocate in 1st chunk (may be zero) 769 . m2 - number of elements to allocate in 2nd chunk (may be zero) 770 . m3 - number of elements to allocate in 3rd chunk (may be zero) 771 - m4 - number of elements to allocate in 4th chunk (may be zero) 772 773 Output Parameter: 774 + r1 - memory allocated in first chunk 775 . r2 - memory allocated in second chunk 776 . r3 - memory allocated in third chunk 777 - r4 - memory allocated in fourth chunk 778 779 Level: developer 780 781 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc4(), PetscFree4() 782 783 M*/ 784 #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)) 785 786 /*MC 787 PetscCalloc4 - Allocates 4 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN 788 789 Synopsis: 790 #include <petscsys.h> 791 PetscErrorCode PetscCalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4) 792 793 Not Collective 794 795 Input Parameters: 796 + m1 - number of elements to allocate in 1st chunk (may be zero) 797 . m2 - number of elements to allocate in 2nd chunk (may be zero) 798 . m3 - number of elements to allocate in 3rd chunk (may be zero) 799 - m4 - number of elements to allocate in 4th chunk (may be zero) 800 801 Output Parameters: 802 + r1 - memory allocated in first chunk 803 . r2 - memory allocated in second chunk 804 . r3 - memory allocated in third chunk 805 - r4 - memory allocated in fourth chunk 806 807 Level: developer 808 809 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc4(), PetscFree4() 810 811 M*/ 812 #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)) 813 814 /*MC 815 PetscMalloc5 - Allocates 5 arrays of memory, all aligned to PETSC_MEMALIGN 816 817 Synopsis: 818 #include <petscsys.h> 819 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) 820 821 Not Collective 822 823 Input Parameters: 824 + m1 - number of elements to allocate in 1st chunk (may be zero) 825 . m2 - number of elements to allocate in 2nd chunk (may be zero) 826 . m3 - number of elements to allocate in 3rd chunk (may be zero) 827 . m4 - number of elements to allocate in 4th chunk (may be zero) 828 - m5 - number of elements to allocate in 5th chunk (may be zero) 829 830 Output Parameters: 831 + r1 - memory allocated in first chunk 832 . r2 - memory allocated in second chunk 833 . r3 - memory allocated in third chunk 834 . r4 - memory allocated in fourth chunk 835 - r5 - memory allocated in fifth chunk 836 837 Level: developer 838 839 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc5(), PetscFree5() 840 841 M*/ 842 #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)) 843 844 /*MC 845 PetscCalloc5 - Allocates 5 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN 846 847 Synopsis: 848 #include <petscsys.h> 849 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) 850 851 Not Collective 852 853 Input Parameters: 854 + m1 - number of elements to allocate in 1st chunk (may be zero) 855 . m2 - number of elements to allocate in 2nd chunk (may be zero) 856 . m3 - number of elements to allocate in 3rd chunk (may be zero) 857 . m4 - number of elements to allocate in 4th chunk (may be zero) 858 - m5 - number of elements to allocate in 5th chunk (may be zero) 859 860 Output Parameters: 861 + r1 - memory allocated in first chunk 862 . r2 - memory allocated in second chunk 863 . r3 - memory allocated in third chunk 864 . r4 - memory allocated in fourth chunk 865 - r5 - memory allocated in fifth chunk 866 867 Level: developer 868 869 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc5(), PetscFree5() 870 871 M*/ 872 #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)) 873 874 /*MC 875 PetscMalloc6 - Allocates 6 arrays of memory, all aligned to PETSC_MEMALIGN 876 877 Synopsis: 878 #include <petscsys.h> 879 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) 880 881 Not Collective 882 883 Input Parameters: 884 + m1 - number of elements to allocate in 1st chunk (may be zero) 885 . m2 - number of elements to allocate in 2nd chunk (may be zero) 886 . m3 - number of elements to allocate in 3rd chunk (may be zero) 887 . m4 - number of elements to allocate in 4th chunk (may be zero) 888 . m5 - number of elements to allocate in 5th chunk (may be zero) 889 - m6 - number of elements to allocate in 6th chunk (may be zero) 890 891 Output Parameteasr: 892 + r1 - memory allocated in first chunk 893 . r2 - memory allocated in second chunk 894 . r3 - memory allocated in third chunk 895 . r4 - memory allocated in fourth chunk 896 . r5 - memory allocated in fifth chunk 897 - r6 - memory allocated in sixth chunk 898 899 Level: developer 900 901 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc6(), PetscFree3(), PetscFree4(), PetscFree5(), PetscFree6() 902 903 M*/ 904 #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)) 905 906 /*MC 907 PetscCalloc6 - Allocates 6 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN 908 909 Synopsis: 910 #include <petscsys.h> 911 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) 912 913 Not Collective 914 915 Input Parameters: 916 + m1 - number of elements to allocate in 1st chunk (may be zero) 917 . m2 - number of elements to allocate in 2nd chunk (may be zero) 918 . m3 - number of elements to allocate in 3rd chunk (may be zero) 919 . m4 - number of elements to allocate in 4th chunk (may be zero) 920 . m5 - number of elements to allocate in 5th chunk (may be zero) 921 - m6 - number of elements to allocate in 6th chunk (may be zero) 922 923 Output Parameters: 924 + r1 - memory allocated in first chunk 925 . r2 - memory allocated in second chunk 926 . r3 - memory allocated in third chunk 927 . r4 - memory allocated in fourth chunk 928 . r5 - memory allocated in fifth chunk 929 - r6 - memory allocated in sixth chunk 930 931 Level: developer 932 933 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscMalloc6(), PetscFree6() 934 935 M*/ 936 #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)) 937 938 /*MC 939 PetscMalloc7 - Allocates 7 arrays of memory, all aligned to PETSC_MEMALIGN 940 941 Synopsis: 942 #include <petscsys.h> 943 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) 944 945 Not Collective 946 947 Input Parameters: 948 + m1 - number of elements to allocate in 1st chunk (may be zero) 949 . m2 - number of elements to allocate in 2nd chunk (may be zero) 950 . m3 - number of elements to allocate in 3rd chunk (may be zero) 951 . m4 - number of elements to allocate in 4th chunk (may be zero) 952 . m5 - number of elements to allocate in 5th chunk (may be zero) 953 . m6 - number of elements to allocate in 6th chunk (may be zero) 954 - m7 - number of elements to allocate in 7th chunk (may be zero) 955 956 Output Parameters: 957 + r1 - memory allocated in first chunk 958 . r2 - memory allocated in second chunk 959 . r3 - memory allocated in third chunk 960 . r4 - memory allocated in fourth chunk 961 . r5 - memory allocated in fifth chunk 962 . r6 - memory allocated in sixth chunk 963 - r7 - memory allocated in seventh chunk 964 965 Level: developer 966 967 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc7(), PetscFree7() 968 969 M*/ 970 #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)) 971 972 /*MC 973 PetscCalloc7 - Allocates 7 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN 974 975 Synopsis: 976 #include <petscsys.h> 977 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) 978 979 Not Collective 980 981 Input Parameters: 982 + m1 - number of elements to allocate in 1st chunk (may be zero) 983 . m2 - number of elements to allocate in 2nd chunk (may be zero) 984 . m3 - number of elements to allocate in 3rd chunk (may be zero) 985 . m4 - number of elements to allocate in 4th chunk (may be zero) 986 . m5 - number of elements to allocate in 5th chunk (may be zero) 987 . m6 - number of elements to allocate in 6th chunk (may be zero) 988 - m7 - number of elements to allocate in 7th chunk (may be zero) 989 990 Output Parameters: 991 + r1 - memory allocated in first chunk 992 . r2 - memory allocated in second chunk 993 . r3 - memory allocated in third chunk 994 . r4 - memory allocated in fourth chunk 995 . r5 - memory allocated in fifth chunk 996 . r6 - memory allocated in sixth chunk 997 - r7 - memory allocated in seventh chunk 998 999 Level: developer 1000 1001 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscMalloc7(), PetscFree7() 1002 1003 M*/ 1004 #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)) 1005 1006 /*MC 1007 PetscNew - Allocates memory of a particular type, zeros the memory! Aligned to PETSC_MEMALIGN 1008 1009 Synopsis: 1010 #include <petscsys.h> 1011 PetscErrorCode PetscNew(type **result) 1012 1013 Not Collective 1014 1015 Output Parameter: 1016 . result - memory allocated, sized to match pointer type 1017 1018 Level: beginner 1019 1020 .seealso: PetscFree(), PetscMalloc(), PetscNewLog(), PetscCalloc1(), PetscMalloc1() 1021 1022 M*/ 1023 #define PetscNew(b) PetscCalloc1(1,(b)) 1024 1025 /*MC 1026 PetscNewLog - Allocates memory of a type matching pointer, zeros the memory! Aligned to PETSC_MEMALIGN. Associates the memory allocated 1027 with the given object using PetscLogObjectMemory(). 1028 1029 Synopsis: 1030 #include <petscsys.h> 1031 PetscErrorCode PetscNewLog(PetscObject obj,type **result) 1032 1033 Not Collective 1034 1035 Input Parameter: 1036 . obj - object memory is logged to 1037 1038 Output Parameter: 1039 . result - memory allocated, sized to match pointer type 1040 1041 Level: developer 1042 1043 .seealso: PetscFree(), PetscMalloc(), PetscNew(), PetscLogObjectMemory(), PetscCalloc1(), PetscMalloc1() 1044 1045 M*/ 1046 #define PetscNewLog(o,b) (PetscNew((b)) || PetscLogObjectMemory((PetscObject)o,sizeof(**(b)))) 1047 1048 /*MC 1049 PetscFree - Frees memory 1050 1051 Synopsis: 1052 #include <petscsys.h> 1053 PetscErrorCode PetscFree(void *memory) 1054 1055 Not Collective 1056 1057 Input Parameter: 1058 . memory - memory to free (the pointer is ALWAYS set to NULL upon success) 1059 1060 Level: beginner 1061 1062 Notes: 1063 Do not free memory obtained with PetscMalloc2(), PetscCalloc2() etc, they must be freed with PetscFree2() etc. 1064 1065 It is safe to call PetscFree() on a NULL pointer. 1066 1067 .seealso: PetscNew(), PetscMalloc(), PetscNewLog(), PetscMalloc1(), PetscCalloc1() 1068 1069 M*/ 1070 #define PetscFree(a) ((*PetscTrFree)((void*)(a),__LINE__,PETSC_FUNCTION_NAME,__FILE__) || ((a) = NULL,0)) 1071 1072 /*MC 1073 PetscFree2 - Frees 2 chunks of memory obtained with PetscMalloc2() 1074 1075 Synopsis: 1076 #include <petscsys.h> 1077 PetscErrorCode PetscFree2(void *memory1,void *memory2) 1078 1079 Not Collective 1080 1081 Input Parameters: 1082 + memory1 - memory to free 1083 - memory2 - 2nd memory to free 1084 1085 Level: developer 1086 1087 Notes: 1088 Memory must have been obtained with PetscMalloc2() 1089 1090 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree() 1091 1092 M*/ 1093 #define PetscFree2(m1,m2) PetscFreeA(2,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2)) 1094 1095 /*MC 1096 PetscFree3 - Frees 3 chunks of memory obtained with PetscMalloc3() 1097 1098 Synopsis: 1099 #include <petscsys.h> 1100 PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3) 1101 1102 Not Collective 1103 1104 Input Parameters: 1105 + memory1 - memory to free 1106 . memory2 - 2nd memory to free 1107 - memory3 - 3rd memory to free 1108 1109 Level: developer 1110 1111 Notes: 1112 Memory must have been obtained with PetscMalloc3() 1113 1114 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3() 1115 1116 M*/ 1117 #define PetscFree3(m1,m2,m3) PetscFreeA(3,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3)) 1118 1119 /*MC 1120 PetscFree4 - Frees 4 chunks of memory obtained with PetscMalloc4() 1121 1122 Synopsis: 1123 #include <petscsys.h> 1124 PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4) 1125 1126 Not Collective 1127 1128 Input Parameters: 1129 + m1 - memory to free 1130 . m2 - 2nd memory to free 1131 . m3 - 3rd memory to free 1132 - m4 - 4th memory to free 1133 1134 Level: developer 1135 1136 Notes: 1137 Memory must have been obtained with PetscMalloc4() 1138 1139 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4() 1140 1141 M*/ 1142 #define PetscFree4(m1,m2,m3,m4) PetscFreeA(4,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4)) 1143 1144 /*MC 1145 PetscFree5 - Frees 5 chunks of memory obtained with PetscMalloc5() 1146 1147 Synopsis: 1148 #include <petscsys.h> 1149 PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5) 1150 1151 Not Collective 1152 1153 Input Parameters: 1154 + m1 - memory to free 1155 . m2 - 2nd memory to free 1156 . m3 - 3rd memory to free 1157 . m4 - 4th memory to free 1158 - m5 - 5th memory to free 1159 1160 Level: developer 1161 1162 Notes: 1163 Memory must have been obtained with PetscMalloc5() 1164 1165 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5() 1166 1167 M*/ 1168 #define PetscFree5(m1,m2,m3,m4,m5) PetscFreeA(5,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5)) 1169 1170 /*MC 1171 PetscFree6 - Frees 6 chunks of memory obtained with PetscMalloc6() 1172 1173 Synopsis: 1174 #include <petscsys.h> 1175 PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6) 1176 1177 Not Collective 1178 1179 Input Parameters: 1180 + m1 - memory to free 1181 . m2 - 2nd memory to free 1182 . m3 - 3rd memory to free 1183 . m4 - 4th memory to free 1184 . m5 - 5th memory to free 1185 - m6 - 6th memory to free 1186 1187 1188 Level: developer 1189 1190 Notes: 1191 Memory must have been obtained with PetscMalloc6() 1192 1193 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6() 1194 1195 M*/ 1196 #define PetscFree6(m1,m2,m3,m4,m5,m6) PetscFreeA(6,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5),&(m6)) 1197 1198 /*MC 1199 PetscFree7 - Frees 7 chunks of memory obtained with PetscMalloc7() 1200 1201 Synopsis: 1202 #include <petscsys.h> 1203 PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7) 1204 1205 Not Collective 1206 1207 Input Parameters: 1208 + m1 - memory to free 1209 . m2 - 2nd memory to free 1210 . m3 - 3rd memory to free 1211 . m4 - 4th memory to free 1212 . m5 - 5th memory to free 1213 . m6 - 6th memory to free 1214 - m7 - 7th memory to free 1215 1216 1217 Level: developer 1218 1219 Notes: 1220 Memory must have been obtained with PetscMalloc7() 1221 1222 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(), 1223 PetscMalloc7() 1224 1225 M*/ 1226 #define PetscFree7(m1,m2,m3,m4,m5,m6,m7) PetscFreeA(7,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5),&(m6),&(m7)) 1227 1228 PETSC_EXTERN PetscErrorCode PetscMallocA(int,PetscBool,int,const char *,const char *,size_t,void *,...); 1229 PETSC_EXTERN PetscErrorCode PetscFreeA(int,int,const char *,const char *,void *,...); 1230 PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,PetscBool,int,const char[],const char[],void**); 1231 PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]); 1232 PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t,int,const char[],const char[],void**); 1233 PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool); 1234 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 **)); 1235 PETSC_EXTERN PetscErrorCode PetscMallocClear(void); 1236 1237 /* 1238 Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair. 1239 */ 1240 PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void); 1241 PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void); 1242 #if defined(PETSC_HAVE_CUDA) 1243 PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void); 1244 PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void); 1245 #endif 1246 1247 #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE 1248 #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION 1249 1250 /* 1251 Routines for tracing memory corruption/bleeding with default PETSc memory allocation 1252 */ 1253 PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *); 1254 PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *); 1255 PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *); 1256 PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *); 1257 PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int); 1258 PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int,PetscLogDouble*); 1259 PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool,PetscBool); 1260 PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*,PetscBool*,PetscBool*); 1261 PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]); 1262 PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble); 1263 PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool*); 1264 PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool); 1265 PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool*); 1266 1267 PETSC_EXTERN const char *const PetscDataTypes[]; 1268 PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*); 1269 PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*); 1270 PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*); 1271 PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*); 1272 1273 /* 1274 Basic memory and string operations. These are usually simple wrappers 1275 around the basic Unix system calls, but a few of them have additional 1276 functionality and/or error checking. 1277 */ 1278 PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool *); 1279 PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*); 1280 PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***); 1281 PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **); 1282 PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool *); 1283 PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool *); 1284 PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *); 1285 PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *); 1286 PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]); 1287 PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]); 1288 PETSC_EXTERN PetscErrorCode PetscStrlcat(char[],const char[],size_t); 1289 PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t); 1290 PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]); 1291 PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]); 1292 PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]); 1293 PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]); 1294 PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]); 1295 PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]); 1296 PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*); 1297 PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*); 1298 PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*); 1299 PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]); 1300 PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***); 1301 PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***); 1302 PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt,const char *const*,char***); 1303 PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt,char***); 1304 PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t); 1305 1306 PETSC_EXTERN void PetscStrcmpNoError(const char[],const char[],PetscBool *); 1307 1308 PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*); 1309 PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]); 1310 PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*); 1311 1312 PETSC_EXTERN PetscErrorCode PetscStrInList(const char[],const char[],char,PetscBool*); 1313 PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*); 1314 PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*); 1315 1316 /* 1317 These are MPI operations for MPI_Allreduce() etc 1318 */ 1319 PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP; 1320 #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 1321 PETSC_EXTERN MPI_Op MPIU_SUM; 1322 #else 1323 #define MPIU_SUM MPI_SUM 1324 #endif 1325 #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16) 1326 PETSC_EXTERN MPI_Op MPIU_MAX; 1327 PETSC_EXTERN MPI_Op MPIU_MIN; 1328 #else 1329 #define MPIU_MAX MPI_MAX 1330 #define MPIU_MIN MPI_MIN 1331 #endif 1332 PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*); 1333 1334 PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm); 1335 PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm); 1336 1337 PETSC_EXTERN const char *const PetscFileModes[]; 1338 1339 /* 1340 Defines PETSc error handling. 1341 */ 1342 #include <petscerror.h> 1343 1344 #define PETSC_SMALLEST_CLASSID 1211211 1345 PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID; 1346 PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID; 1347 PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *); 1348 PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject,PetscObjectId*); 1349 PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject,PetscObjectId,PetscBool*); 1350 1351 1352 /* 1353 Routines that get memory usage information from the OS 1354 */ 1355 PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *); 1356 PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *); 1357 PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void); 1358 PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]); 1359 1360 PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal); 1361 1362 /* 1363 Initialization of PETSc 1364 */ 1365 PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]); 1366 PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]); 1367 PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void); 1368 PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *); 1369 PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *); 1370 PETSC_EXTERN PetscErrorCode PetscFinalize(void); 1371 PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void); 1372 PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***); 1373 PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***); 1374 PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **); 1375 1376 PETSC_EXTERN PetscErrorCode PetscEnd(void); 1377 PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void); 1378 1379 PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]); 1380 PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void); 1381 PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void); 1382 PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]); 1383 1384 PETSC_EXTERN PetscErrorCode PetscMonitorCompare(PetscErrorCode (*)(void),void *,PetscErrorCode (*)(void**),PetscErrorCode (*)(void),void *,PetscErrorCode (*)(void**),PetscBool *); 1385 1386 1387 /* 1388 These are so that in extern C code we can caste function pointers to non-extern C 1389 function pointers. Since the regular C++ code expects its function pointers to be C++ 1390 */ 1391 PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void); 1392 PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void); 1393 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void); 1394 1395 /* 1396 Functions that can act on any PETSc object. 1397 */ 1398 PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*); 1399 PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *); 1400 PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *); 1401 PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]); 1402 PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []); 1403 PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]); 1404 PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]); 1405 PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]); 1406 PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt); 1407 PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*); 1408 PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt); 1409 PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject); 1410 PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*); 1411 PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject); 1412 PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt*); 1413 PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject); 1414 PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]); 1415 PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject*); 1416 PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void)); 1417 #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d)) 1418 PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject); 1419 PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject); 1420 PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject); 1421 PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject,PetscObject); 1422 PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt*); 1423 1424 #include <petscviewertypes.h> 1425 #include <petscoptions.h> 1426 1427 PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer,PetscBool,PetscLogDouble); 1428 PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool*); 1429 1430 PETSC_EXTERN PetscErrorCode PetscObjectsListGetGlobalNumbering(MPI_Comm,PetscInt,PetscObject*,PetscInt*,PetscInt*); 1431 1432 PETSC_EXTERN PetscErrorCode PetscMemoryShowUsage(PetscViewer,const char[]); 1433 PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer,const char[]); 1434 PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer); 1435 PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer); 1436 #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr)) 1437 PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void)); 1438 PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]); 1439 PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]); 1440 PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]); 1441 PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]); 1442 PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]); 1443 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject); 1444 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void); 1445 PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,PetscObject,const char[]); 1446 PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject); 1447 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *); 1448 PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject,const char[],PetscBool *); 1449 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...); 1450 PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject,PetscBool*,const char[],...); 1451 PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void)); 1452 PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void); 1453 1454 #if defined(PETSC_HAVE_SAWS) 1455 PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void); 1456 PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject); 1457 PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool); 1458 PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject); 1459 PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject); 1460 PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject); 1461 PETSC_EXTERN void PetscStackSAWsGrantAccess(void); 1462 PETSC_EXTERN void PetscStackSAWsTakeAccess(void); 1463 PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void); 1464 PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void); 1465 1466 #else 1467 #define PetscSAWsBlock() 0 1468 #define PetscObjectSAWsViewOff(obj) 0 1469 #define PetscObjectSAWsSetBlock(obj,flg) 0 1470 #define PetscObjectSAWsBlock(obj) 0 1471 #define PetscObjectSAWsGrantAccess(obj) 0 1472 #define PetscObjectSAWsTakeAccess(obj) 0 1473 #define PetscStackViewSAWs() 0 1474 #define PetscStackSAWsViewOff() 0 1475 #define PetscStackSAWsTakeAccess() 1476 #define PetscStackSAWsGrantAccess() 1477 1478 #endif 1479 1480 PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *); 1481 PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *); 1482 PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **); 1483 1484 #if defined(PETSC_USE_DEBUG) 1485 PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**); 1486 #endif 1487 PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool); 1488 1489 PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*); 1490 PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*); 1491 PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*); 1492 PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject); 1493 PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]); 1494 PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *); 1495 1496 /* 1497 Dynamic library lists. Lists of names of routines in objects or in dynamic 1498 link libraries that will be loaded as needed. 1499 */ 1500 1501 #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr)) 1502 PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void)); 1503 PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*); 1504 #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr)) 1505 PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void)); 1506 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[],const char[]); 1507 PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *); 1508 PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer); 1509 PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*); 1510 1511 PETSC_EXTERN PetscDLLibrary PetscDLLibrariesLoaded; 1512 PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]); 1513 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]); 1514 PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **); 1515 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary); 1516 PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool *); 1517 PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *); 1518 PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary); 1519 1520 /* 1521 Useful utility routines 1522 */ 1523 PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*); 1524 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*); 1525 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm,PetscInt*,PetscInt*); 1526 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt); 1527 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt); 1528 PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject); 1529 PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*); 1530 PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm,PetscInt[],PetscInt[]); 1531 PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm,PetscReal[],PetscReal[]); 1532 1533 /*MC 1534 PetscNot - negates a logical type value and returns result as a PetscBool 1535 1536 Notes: 1537 This is useful in cases like 1538 $ int *a; 1539 $ PetscBool flag = PetscNot(a) 1540 where !a would not return a PetscBool because we cannot provide a cast from int to PetscBool in C. 1541 1542 Level: beginner 1543 1544 .seealso : PetscBool, PETSC_TRUE, PETSC_FALSE 1545 M*/ 1546 #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE) 1547 1548 /*MC 1549 PetscHelpPrintf - Prints help messages. 1550 1551 Synopsis: 1552 #include <petscsys.h> 1553 PetscErrorCode (*PetscHelpPrintf)(const char format[],...); 1554 1555 Not Collective 1556 1557 Input Parameters: 1558 . format - the usual printf() format string 1559 1560 Level: developer 1561 1562 Fortran Note: 1563 This routine is not supported in Fortran. 1564 1565 1566 .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf() 1567 M*/ 1568 PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...); 1569 1570 /* 1571 Defines PETSc profiling. 1572 */ 1573 #include <petsclog.h> 1574 1575 /* 1576 Simple PETSc parallel IO for ASCII printing 1577 */ 1578 PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]); 1579 PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**); 1580 PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*); 1581 PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...); 1582 PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...); 1583 PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...); 1584 PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...); 1585 PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[],size_t,const char*,PetscInt,const PetscReal[]); 1586 1587 PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...); 1588 PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...); 1589 PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...); 1590 1591 PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char*,size_t*); 1592 PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char*,char *); 1593 1594 #if defined(PETSC_HAVE_POPEN) 1595 PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **); 1596 PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*); 1597 PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]); 1598 #endif 1599 1600 PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...); 1601 PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...); 1602 PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*); 1603 PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]); 1604 PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**); 1605 PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**); 1606 PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]); 1607 1608 PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID; 1609 PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void**); 1610 PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void*); 1611 PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*); 1612 PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer*); 1613 PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer,PetscErrorCode (*)(void*)); 1614 PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void*); 1615 1616 /* 1617 For use in debuggers 1618 */ 1619 PETSC_EXTERN PetscMPIInt PetscGlobalRank; 1620 PETSC_EXTERN PetscMPIInt PetscGlobalSize; 1621 PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer); 1622 PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer); 1623 PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer); 1624 1625 #include <stddef.h> 1626 #include <string.h> /* for memcpy, memset */ 1627 #include <stdlib.h> 1628 1629 #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__) 1630 #include <xmmintrin.h> 1631 #endif 1632 1633 /*@C 1634 PetscMemmove - Copies n bytes, beginning at location b, to the space 1635 beginning at location a. Copying between regions that overlap will 1636 take place correctly. Use PetscMemcpy() if the locations do not overlap 1637 1638 Not Collective 1639 1640 Input Parameters: 1641 + b - pointer to initial memory space 1642 . a - pointer to copy space 1643 - n - length (in bytes) of space to copy 1644 1645 Level: intermediate 1646 1647 Note: 1648 PetscArraymove() is preferred 1649 This routine is analogous to memmove(). 1650 1651 Developers Note: This is inlined for performance 1652 1653 .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscStrallocpy(), 1654 PetscArraymove() 1655 @*/ 1656 PETSC_STATIC_INLINE PetscErrorCode PetscMemmove(void *a,const void *b,size_t n) 1657 { 1658 PetscFunctionBegin; 1659 if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to null pointer"); 1660 if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer"); 1661 #if !defined(PETSC_HAVE_MEMMOVE) 1662 if (a < b) { 1663 if (a <= b - n) memcpy(a,b,n); 1664 else { 1665 memcpy(a,b,(int)(b - a)); 1666 PetscMemmove(b,b + (int)(b - a),n - (int)(b - a)); 1667 } 1668 } else { 1669 if (b <= a - n) memcpy(a,b,n); 1670 else { 1671 memcpy(b + n,b + (n - (int)(a - b)),(int)(a - b)); 1672 PetscMemmove(a,b,n - (int)(a - b)); 1673 } 1674 } 1675 #else 1676 memmove((char*)(a),(char*)(b),n); 1677 #endif 1678 PetscFunctionReturn(0); 1679 } 1680 1681 /*@C 1682 PetscMemcpy - Copies n bytes, beginning at location b, to the space 1683 beginning at location a. The two memory regions CANNOT overlap, use 1684 PetscMemmove() in that case. 1685 1686 Not Collective 1687 1688 Input Parameters: 1689 + b - pointer to initial memory space 1690 - n - length (in bytes) of space to copy 1691 1692 Output Parameter: 1693 . a - pointer to copy space 1694 1695 Level: intermediate 1696 1697 Compile Option: 1698 PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used 1699 for memory copies on double precision values. 1700 PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used 1701 for memory copies on double precision values. 1702 PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used 1703 for memory copies on double precision values. 1704 1705 Note: 1706 Prefer PetscArraycpy() 1707 This routine is analogous to memcpy(). 1708 Not available from Fortran 1709 1710 Developer Note: this is inlined for fastest performance 1711 1712 .seealso: PetscMemzero(), PetscMemcmp(), PetscArrayzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy() 1713 1714 @*/ 1715 PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n) 1716 { 1717 #if defined(PETSC_USE_DEBUG) 1718 size_t al = (size_t) a,bl = (size_t) b; 1719 size_t nl = (size_t) n; 1720 PetscFunctionBegin; 1721 if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer"); 1722 if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer"); 1723 #else 1724 PetscFunctionBegin; 1725 #endif 1726 if (a != b && n > 0) { 1727 #if defined(PETSC_USE_DEBUG) 1728 if ((al > bl && (al - bl) < nl) || (bl - al) < nl) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\ 1729 or make sure your copy regions and lengths are correct. \n\ 1730 Length (bytes) %ld first address %ld second address %ld",nl,al,bl); 1731 #endif 1732 #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY)) 1733 if (!(a % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) { 1734 size_t len = n/sizeof(PetscScalar); 1735 #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) 1736 PetscBLASInt one = 1,blen; 1737 PetscErrorCode ierr; 1738 ierr = PetscBLASIntCast(len,&blen);CHKERRQ(ierr); 1739 PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one)); 1740 #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY) 1741 fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a); 1742 #else 1743 size_t i; 1744 PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a; 1745 for (i=0; i<len; i++) y[i] = x[i]; 1746 #endif 1747 } else { 1748 memcpy((char*)(a),(char*)(b),n); 1749 } 1750 #else 1751 memcpy((char*)(a),(char*)(b),n); 1752 #endif 1753 } 1754 PetscFunctionReturn(0); 1755 } 1756 1757 /*@C 1758 PetscMemzero - Zeros the specified memory. 1759 1760 Not Collective 1761 1762 Input Parameters: 1763 + a - pointer to beginning memory location 1764 - n - length (in bytes) of memory to initialize 1765 1766 Level: intermediate 1767 1768 Compile Option: 1769 PETSC_PREFER_BZERO - on certain machines (the IBM RS6000) the bzero() routine happens 1770 to be faster than the memset() routine. This flag causes the bzero() routine to be used. 1771 1772 Not available from Fortran 1773 Prefer PetscArrayzero() 1774 1775 Developer Note: this is inlined for fastest performance 1776 1777 .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy() 1778 @*/ 1779 PETSC_STATIC_INLINE PetscErrorCode PetscMemzero(void *a,size_t n) 1780 { 1781 if (n > 0) { 1782 #if defined(PETSC_USE_DEBUG) 1783 if (!a) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer with %zu bytes",n); 1784 #endif 1785 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) 1786 if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) { 1787 size_t i,len = n/sizeof(PetscScalar); 1788 PetscScalar *x = (PetscScalar*)a; 1789 for (i=0; i<len; i++) x[i] = 0.0; 1790 } else { 1791 #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO) 1792 if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) { 1793 PetscInt len = n/sizeof(PetscScalar); 1794 fortranzero_(&len,(PetscScalar*)a); 1795 } else { 1796 #endif 1797 #if defined(PETSC_PREFER_BZERO) 1798 bzero((char *)a,n); 1799 #else 1800 memset((char*)a,0,n); 1801 #endif 1802 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO) 1803 } 1804 #endif 1805 } 1806 return 0; 1807 } 1808 1809 /*MC 1810 PetscArraycmp - Compares two arrays in memory. 1811 1812 Synopsis: 1813 #include <petscsys.h> 1814 PetscErrorCode PetscArraycmp(const anytype *str1,const anytype *str2,size_t cnt,PetscBool *e) 1815 1816 Not Collective 1817 1818 Input Parameters: 1819 + str1 - First array 1820 . str2 - Second array 1821 - cnt - Count of the array, not in bytes, but number of entries in the arrays 1822 1823 Output Parameters: 1824 . e - PETSC_TRUE if equal else PETSC_FALSE. 1825 1826 Level: intermediate 1827 1828 Note: 1829 This routine is a preferred replacement to PetscMemcmp() 1830 The arrays must be of the same type 1831 1832 .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy(), 1833 PetscArraymove() 1834 M*/ 1835 #define PetscArraycmp(str1,str2,cnt,e) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcmp(str1,str2,(size_t)(cnt)*sizeof(*(str1)),e)) 1836 1837 /*MC 1838 PetscArraymove - Copies from one array in memory to another, the arrays may overlap. Use PetscArraycpy() when the arrays 1839 do not overlap 1840 1841 Synopsis: 1842 #include <petscsys.h> 1843 PetscErrorCode PetscArraymove(anytype *str1,const anytype *str2,size_t cnt) 1844 1845 Not Collective 1846 1847 Input Parameters: 1848 + str1 - First array 1849 . str2 - Second array 1850 - cnt - Count of the array, not in bytes, but number of entries in the arrays 1851 1852 Level: intermediate 1853 1854 Note: 1855 This routine is a preferred replacement to PetscMemmove() 1856 The arrays must be of the same type 1857 1858 .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycpy(), PetscMemmove(), PetscArraycmp(), PetscStrallocpy() 1859 M*/ 1860 #define PetscArraymove(str1,str2,cnt) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemmove(str1,str2,(size_t)(cnt)*sizeof(*(str1)))) 1861 1862 /*MC 1863 PetscArraycpy - Copies from one array in memory to another 1864 1865 Synopsis: 1866 #include <petscsys.h> 1867 PetscErrorCode PetscArraycpy(anytype *str1,const anytype *str2,size_t cnt) 1868 1869 Not Collective 1870 1871 Input Parameters: 1872 + str1 - First array (destination) 1873 . str2 - Second array (source) 1874 - cnt - Count of the array, not in bytes, but number of entries in the arrays 1875 1876 Level: intermediate 1877 1878 Note: 1879 This routine is a preferred replacement to PetscMemcpy() 1880 The arrays must be of the same type 1881 1882 .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraymove(), PetscMemmove(), PetscArraycmp(), PetscStrallocpy() 1883 M*/ 1884 #define PetscArraycpy(str1,str2,cnt) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcpy(str1,str2,(size_t)(cnt)*sizeof(*(str1)))) 1885 1886 /*MC 1887 PetscArrayzero - Zeros an array in memory. 1888 1889 Synopsis: 1890 #include <petscsys.h> 1891 PetscErrorCode PetscArrayzero(anytype *str1,size_t cnt) 1892 1893 Not Collective 1894 1895 Input Parameters: 1896 + str1 - array 1897 - cnt - Count of the array, not in bytes, but number of entries in the array 1898 1899 Level: intermediate 1900 1901 Note: 1902 This routine is a preferred replacement to PetscMemzero() 1903 1904 .seealso: PetscMemcpy(), PetscMemcmp(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy(), PetscArraymove() 1905 M*/ 1906 #define PetscArrayzero(str1,cnt) PetscMemzero(str1,(size_t)(cnt)*sizeof(*(str1))) 1907 1908 /*MC 1909 PetscPrefetchBlock - Prefetches a block of memory 1910 1911 Synopsis: 1912 #include <petscsys.h> 1913 void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t) 1914 1915 Not Collective 1916 1917 Input Parameters: 1918 + a - pointer to first element to fetch (any type but usually PetscInt or PetscScalar) 1919 . n - number of elements to fetch 1920 . rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors) 1921 - t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note 1922 1923 Level: developer 1924 1925 Notes: 1926 The last two arguments (rw and t) must be compile-time constants. 1927 1928 Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality. Not all architectures offer 1929 equivalent locality hints, but the following macros are always defined to their closest analogue. 1930 + PETSC_PREFETCH_HINT_NTA - Non-temporal. Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched). 1931 . 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. 1932 . PETSC_PREFETCH_HINT_T1 - Fetch to level 2 and higher (not L1). 1933 - PETSC_PREFETCH_HINT_T2 - Fetch to high-level cache only. (On many systems, T0 and T1 are equivalent.) 1934 1935 This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid 1936 address). 1937 1938 M*/ 1939 #define PetscPrefetchBlock(a,n,rw,t) do { \ 1940 const char *_p = (const char*)(a),*_end = (const char*)((a)+(n)); \ 1941 for (; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \ 1942 } while (0) 1943 1944 /* 1945 Determine if some of the kernel computation routines use 1946 Fortran (rather than C) for the numerical calculations. On some machines 1947 and compilers (like complex numbers) the Fortran version of the routines 1948 is faster than the C/C++ versions. The flag --with-fortran-kernels 1949 should be used with ./configure to turn these on. 1950 */ 1951 #if defined(PETSC_USE_FORTRAN_KERNELS) 1952 1953 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL) 1954 #define PETSC_USE_FORTRAN_KERNEL_MULTCRL 1955 #endif 1956 1957 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM) 1958 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM 1959 #endif 1960 1961 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ) 1962 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ 1963 #endif 1964 1965 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ) 1966 #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ 1967 #endif 1968 1969 #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM) 1970 #define PETSC_USE_FORTRAN_KERNEL_NORM 1971 #endif 1972 1973 #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY) 1974 #define PETSC_USE_FORTRAN_KERNEL_MAXPY 1975 #endif 1976 1977 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 1978 #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ 1979 #endif 1980 1981 #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ) 1982 #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ 1983 #endif 1984 1985 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ) 1986 #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ 1987 #endif 1988 1989 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ) 1990 #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ 1991 #endif 1992 1993 #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT) 1994 #define PETSC_USE_FORTRAN_KERNEL_MDOT 1995 #endif 1996 1997 #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY) 1998 #define PETSC_USE_FORTRAN_KERNEL_XTIMESY 1999 #endif 2000 2001 #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX) 2002 #define PETSC_USE_FORTRAN_KERNEL_AYPX 2003 #endif 2004 2005 #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY) 2006 #define PETSC_USE_FORTRAN_KERNEL_WAXPY 2007 #endif 2008 2009 #endif 2010 2011 /* 2012 Macros for indicating code that should be compiled with a C interface, 2013 rather than a C++ interface. Any routines that are dynamically loaded 2014 (such as the PCCreate_XXX() routines) must be wrapped so that the name 2015 mangler does not change the functions symbol name. This just hides the 2016 ugly extern "C" {} wrappers. 2017 */ 2018 #if defined(__cplusplus) 2019 # define EXTERN_C_BEGIN extern "C" { 2020 # define EXTERN_C_END } 2021 #else 2022 # define EXTERN_C_BEGIN 2023 # define EXTERN_C_END 2024 #endif 2025 2026 /* --------------------------------------------------------------------*/ 2027 2028 /*MC 2029 MPI_Comm - the basic object used by MPI to determine which processes are involved in a 2030 communication 2031 2032 Level: beginner 2033 2034 Note: This manual page is a place-holder because MPICH does not have a manual page for MPI_Comm 2035 2036 .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF 2037 M*/ 2038 2039 #if defined(PETSC_HAVE_MPIIO) 2040 PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*); 2041 PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*); 2042 PETSC_EXTERN PetscErrorCode MPIU_File_write_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*); 2043 PETSC_EXTERN PetscErrorCode MPIU_File_read_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*); 2044 PETSC_EXTERN PetscErrorCode MPIU_File_write_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*); 2045 PETSC_EXTERN PetscErrorCode MPIU_File_read_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*); 2046 #endif 2047 2048 /* the following petsc_static_inline require petscerror.h */ 2049 2050 /* Limit MPI to 32-bits */ 2051 #define PETSC_MPI_INT_MAX 2147483647 2052 #define PETSC_MPI_INT_MIN -2147483647 2053 /* Limit BLAS to 32-bits */ 2054 #define PETSC_BLAS_INT_MAX 2147483647 2055 #define PETSC_BLAS_INT_MIN -2147483647 2056 2057 /*@C 2058 PetscIntCast - casts a PetscInt64 (which is 64 bits in size) to a PetscInt (which may be 32 bits in size), generates an 2059 error if the PetscInt is not large enough to hold the number. 2060 2061 Not Collective 2062 2063 Input Parameter: 2064 . a - the PetscInt64 value 2065 2066 Output Parameter: 2067 . b - the resulting PetscInt value 2068 2069 Level: advanced 2070 2071 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 2072 2073 Not available from Fortran 2074 2075 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscBLASIntCast() 2076 @*/ 2077 PETSC_STATIC_INLINE PetscErrorCode PetscIntCast(PetscInt64 a,PetscInt *b) 2078 { 2079 PetscFunctionBegin; 2080 #if !defined(PETSC_USE_64BIT_INDICES) 2081 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"); 2082 #endif 2083 *b = (PetscInt)(a); 2084 PetscFunctionReturn(0); 2085 } 2086 2087 /*@C 2088 PetscBLASIntCast - casts a PetscInt (which may be 64 bits in size) to a PetscBLASInt (which may be 32 bits in size), generates an 2089 error if the PetscBLASInt is not large enough to hold the number. 2090 2091 Not Collective 2092 2093 Input Parameter: 2094 . a - the PetscInt value 2095 2096 Output Parameter: 2097 . b - the resulting PetscBLASInt value 2098 2099 Level: advanced 2100 2101 Not available from Fortran 2102 2103 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscIntCast() 2104 @*/ 2105 PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b) 2106 { 2107 PetscFunctionBegin; 2108 #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES) 2109 *b = 0; 2110 if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK"); 2111 #endif 2112 *b = (PetscBLASInt)(a); 2113 PetscFunctionReturn(0); 2114 } 2115 2116 /*@C 2117 PetscMPIIntCast - casts a PetscInt (which may be 64 bits in size) to a PetscMPIInt (which may be 32 bits in size), generates an 2118 error if the PetscMPIInt is not large enough to hold the number. 2119 2120 Not Collective 2121 2122 Input Parameter: 2123 . a - the PetscInt value 2124 2125 Output Parameter: 2126 . b - the resulting PetscMPIInt value 2127 2128 Level: advanced 2129 2130 Not available from Fortran 2131 2132 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntCast() 2133 @*/ 2134 PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b) 2135 { 2136 PetscFunctionBegin; 2137 #if defined(PETSC_USE_64BIT_INDICES) 2138 *b = 0; 2139 if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for MPI"); 2140 #endif 2141 *b = (PetscMPIInt)(a); 2142 PetscFunctionReturn(0); 2143 } 2144 2145 #define PetscInt64Mult(a,b) ((PetscInt64)(a))*((PetscInt64)(b)) 2146 2147 /*@C 2148 2149 PetscRealIntMultTruncate - Computes the product of a positive PetscReal and a positive PetscInt and truncates the value to slightly less than the maximal possible value 2150 2151 Not Collective 2152 2153 Input Parameter: 2154 + a - the PetscReal value 2155 - b - the second value 2156 2157 Returns: 2158 the result as a PetscInt value 2159 2160 Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64 2161 Use PetscIntMultTruncate() to compute the product of two positive PetscInt and truncate to fit a PetscInt 2162 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 2163 2164 Developers Note: 2165 We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks. 2166 2167 This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able. 2168 2169 Not available from Fortran 2170 2171 Level: advanced 2172 2173 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult() 2174 @*/ 2175 PETSC_STATIC_INLINE PetscInt PetscRealIntMultTruncate(PetscReal a,PetscInt b) 2176 { 2177 PetscInt64 r; 2178 2179 r = (PetscInt64) (a*(PetscReal)b); 2180 if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100; 2181 return (PetscInt) r; 2182 } 2183 2184 /*@C 2185 2186 PetscIntMultTruncate - Computes the product of two positive PetscInt and truncates the value to slightly less than the maximal possible value 2187 2188 Not Collective 2189 2190 Input Parameter: 2191 + a - the PetscInt value 2192 - b - the second value 2193 2194 Returns: 2195 the result as a PetscInt value 2196 2197 Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64 2198 Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt 2199 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 2200 2201 Not available from Fortran 2202 2203 Developers Note: We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks. 2204 2205 This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able. 2206 2207 Level: advanced 2208 2209 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult() 2210 @*/ 2211 PETSC_STATIC_INLINE PetscInt PetscIntMultTruncate(PetscInt a,PetscInt b) 2212 { 2213 PetscInt64 r; 2214 2215 r = PetscInt64Mult(a,b); 2216 if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100; 2217 return (PetscInt) r; 2218 } 2219 2220 /*@C 2221 2222 PetscIntSumTruncate - Computes the sum of two positive PetscInt and truncates the value to slightly less than the maximal possible value 2223 2224 Not Collective 2225 2226 Input Parameter: 2227 + a - the PetscInt value 2228 - b - the second value 2229 2230 Returns: 2231 the result as a PetscInt value 2232 2233 Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64 2234 Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt 2235 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 2236 2237 This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able. 2238 2239 Not available from Fortran 2240 2241 Level: advanced 2242 2243 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult() 2244 @*/ 2245 PETSC_STATIC_INLINE PetscInt PetscIntSumTruncate(PetscInt a,PetscInt b) 2246 { 2247 PetscInt64 r; 2248 2249 r = ((PetscInt64)a) + ((PetscInt64)b); 2250 if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100; 2251 return (PetscInt) r; 2252 } 2253 2254 /*@C 2255 2256 PetscIntMultError - Computes the product of two positive PetscInt and generates an error with overflow. 2257 2258 Not Collective 2259 2260 Input Parameter: 2261 + a - the PetscInt value 2262 - b - the second value 2263 2264 Output Parameter:ma 2265 . result - the result as a PetscInt value, or NULL if you do not want the result, you just want to check if it overflows 2266 2267 Use PetscInt64Mult() to compute the product of two 32 bit PetscInt and store in a PetscInt64 2268 Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt 2269 2270 Not available from Fortran 2271 2272 Developers Note: We currently assume that PetscInt addition does not overflow, this is obviously wrong but requires many more checks. 2273 2274 Level: advanced 2275 2276 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64(), PetscIntSumError() 2277 @*/ 2278 PETSC_STATIC_INLINE PetscErrorCode PetscIntMultError(PetscInt a,PetscInt b,PetscInt *result) 2279 { 2280 PetscInt64 r; 2281 2282 PetscFunctionBegin; 2283 r = PetscInt64Mult(a,b); 2284 #if !defined(PETSC_USE_64BIT_INDICES) 2285 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); 2286 #endif 2287 if (result) *result = (PetscInt) r; 2288 PetscFunctionReturn(0); 2289 } 2290 2291 /*@C 2292 2293 PetscIntSumError - Computes the sum of two positive PetscInt and generates an error with overflow. 2294 2295 Not Collective 2296 2297 Input Parameter: 2298 + a - the PetscInt value 2299 - b - the second value 2300 2301 Output Parameter:ma 2302 . c - the result as a PetscInt value, or NULL if you do not want the result, you just want to check if it overflows 2303 2304 Use PetscInt64Mult() to compute the product of two 32 bit PetscInt and store in a PetscInt64 2305 Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt 2306 2307 Not available from Fortran 2308 2309 Level: advanced 2310 2311 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult() 2312 @*/ 2313 PETSC_STATIC_INLINE PetscErrorCode PetscIntSumError(PetscInt a,PetscInt b,PetscInt *result) 2314 { 2315 PetscInt64 r; 2316 2317 PetscFunctionBegin; 2318 r = ((PetscInt64)a) + ((PetscInt64)b); 2319 #if !defined(PETSC_USE_64BIT_INDICES) 2320 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); 2321 #endif 2322 if (result) *result = (PetscInt) r; 2323 PetscFunctionReturn(0); 2324 } 2325 2326 /* 2327 The IBM include files define hz, here we hide it so that it may be used as a regular user variable. 2328 */ 2329 #if defined(hz) 2330 # undef hz 2331 #endif 2332 2333 #include <limits.h> 2334 2335 /* The number of bits in a byte */ 2336 2337 #define PETSC_BITS_PER_BYTE CHAR_BIT 2338 2339 #if defined(PETSC_HAVE_SYS_TYPES_H) 2340 # include <sys/types.h> 2341 #endif 2342 2343 /*MC 2344 2345 PETSC_VERSION - This manual page provides information about how PETSc documents and uses its version information. This information is available to both C/C++ 2346 and Fortran compilers when petscsys.h is included. 2347 2348 2349 The current PETSc version and the API for accessing it are defined in petscversion.h 2350 2351 The complete version number is given as the triple PETSC_VERSION_MAJOR.PETSC_VERSION_MINOR.PETSC_VERSION_SUBMINOR (in short hand x.y.z) 2352 2353 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 2354 where only a change in the major version number (x) indicates a change in the API. 2355 2356 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 2357 2358 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), 2359 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 2360 version number (x.y.z). 2361 2362 PETSC_RELEASE_DATE is the date the x.y version was released (i.e. the version before any patch releases) 2363 2364 PETSC_VERSION_DATE is the date the x.y.z version was released 2365 2366 PETSC_VERSION_GIT is the last git commit to the repository given in the form vx.y.z-wwwww 2367 2368 PETSC_VERSION_DATE_GIT is the date of the last git commit to the repository 2369 2370 Level: intermediate 2371 2372 PETSC_VERSION_() and PETSC_VERSION_PATCH are deprecated and will eventually be removed. For several releases PETSC_VERSION_PATCH is always 0 2373 2374 M*/ 2375 2376 /*MC 2377 2378 UsingFortran - To use PETSc with Fortran you must use both PETSc include files and modules. At the beginning 2379 of every function and module definition you need something like 2380 2381 $ 2382 $#include "petsc/finclude/petscXXX.h" 2383 $ use petscXXX 2384 2385 You can declare PETSc variables using either of the following. 2386 2387 $ XXX variablename 2388 $ type(tXXX) variablename 2389 2390 For example, 2391 2392 $#include "petsc/finclude/petscvec.h" 2393 $ use petscvec 2394 $ 2395 $ Vec b 2396 $ type(tVec) x 2397 2398 Level: beginner 2399 2400 M*/ 2401 2402 PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t); 2403 PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t); 2404 PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t); 2405 PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t); 2406 PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]); 2407 PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t); 2408 PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t); 2409 PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt*,PetscInt*,PetscInt*,PetscInt*); 2410 2411 PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt,const PetscInt[],PetscBool*); 2412 PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt,const PetscMPIInt[],PetscBool*); 2413 PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt,const PetscReal[],PetscBool*); 2414 PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]); 2415 PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt,PetscInt[]); 2416 PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt*,PetscInt[]); 2417 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]); 2418 PETSC_EXTERN PetscErrorCode PetscCheckDupsInt(PetscInt,const PetscInt[],PetscBool*); 2419 PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*); 2420 PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt*); 2421 PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]); 2422 PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]); 2423 PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]); 2424 PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt[],PetscInt[],PetscInt[]); 2425 PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]); 2426 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]); 2427 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]); 2428 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt,PetscMPIInt[],PetscInt[]); 2429 PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]); 2430 PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt,PetscInt[],void*,size_t,void*); 2431 PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]); 2432 PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt,PetscReal[],PetscInt[]); 2433 PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]); 2434 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]); 2435 PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal,PetscInt,const PetscReal[],PetscReal,PetscInt*); 2436 PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]); 2437 PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]); 2438 PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**); 2439 PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt[],const PetscInt[],PetscInt,const PetscInt[],const PetscInt[],PetscInt*,PetscInt**,PetscInt**); 2440 PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt*,PetscInt**); 2441 PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt,const PetscMPIInt[],PetscInt,const PetscMPIInt[],PetscInt*,PetscMPIInt**); 2442 PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *); 2443 2444 PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt,void*,size_t,int(*)(const void*,const void*)); 2445 PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt,PetscInt[]); 2446 PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt,PetscMPIInt[]); 2447 PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt,PetscReal[]); 2448 PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt,void*,size_t,void*,size_t,int(*)(const void*,const void*)); 2449 PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt,PetscInt[],PetscInt[]); 2450 PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt,PetscMPIInt[],PetscMPIInt[]); 2451 PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt,PetscReal[],PetscInt[]); 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