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