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