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