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