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