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