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