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