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 /* 1176 Variable type where we stash PETSc object pointers in Fortran. 1177 On most machines size(pointer) == sizeof(long) - except Microsoft Windows where its sizeof(long long). 1178 Note we could just use size_t here, but then PETSC_VIEWER_*_*_FORTRAN would need to be given as positive numbers. 1179 */ 1180 1181 #if (PETSC_SIZEOF_VOID_P) == (PETSC_SIZEOF_LONG) 1182 #define PetscFortranAddr long 1183 #elif (PETSC_SIZEOF_VOID_P) == (PETSC_SIZEOF_LONG_LONG) 1184 #define PetscFortranAddr long long 1185 #else 1186 #error "Unknown size for PetscFortranAddr! Send us a bugreport at petsc-maint@mcs.anl.gov" 1187 #endif 1188 1189 /*E 1190 PetscDataType - Used for handling different basic data types. 1191 1192 Level: beginner 1193 1194 Developer comment: It would be nice if we could always just use MPI Datatypes, why can we not? 1195 1196 .seealso: PetscBinaryRead(), PetscBinaryWrite(), PetscDataTypeToMPIDataType(), 1197 PetscDataTypeGetSize() 1198 1199 E*/ 1200 typedef enum {PETSC_INT = 0,PETSC_DOUBLE = 1,PETSC_COMPLEX = 2, PETSC_LONG = 3 ,PETSC_SHORT = 4,PETSC_FLOAT = 5, 1201 PETSC_CHAR = 6,PETSC_BIT_LOGICAL = 7,PETSC_ENUM = 8,PETSC_BOOL=9, PETSC_LONG_DOUBLE = 10, PETSC_QD_DD = 11} PetscDataType; 1202 extern const char *PetscDataTypes[]; 1203 1204 #if defined(PETSC_USE_COMPLEX) 1205 #define PETSC_SCALAR PETSC_COMPLEX 1206 #else 1207 #if defined(PETSC_USE_SCALAR_SINGLE) 1208 #define PETSC_SCALAR PETSC_FLOAT 1209 #elif defined(PETSC_USE_SCALAR_LONG_DOUBLE) 1210 #define PETSC_SCALAR PETSC_LONG_DOUBLE 1211 #elif defined(PETSC_USE_SCALAR_INT) 1212 #define PETSC_SCALAR PETSC_INT 1213 #elif defined(PETSC_USE_SCALAR_QD_DD) 1214 #define PETSC_SCALAR PETSC_QD_DD 1215 #else 1216 #define PETSC_SCALAR PETSC_DOUBLE 1217 #endif 1218 #endif 1219 #if defined(PETSC_USE_SCALAR_SINGLE) 1220 #define PETSC_REAL PETSC_FLOAT 1221 #elif defined(PETSC_USE_SCALAR_LONG_DOUBLE) 1222 #define PETSC_REAL PETSC_LONG_DOUBLE 1223 #elif defined(PETSC_USE_SCALAR_INT) 1224 #define PETSC_REAL PETSC_INT 1225 #elif defined(PETSC_USE_SCALAR_QD_DD) 1226 #define PETSC_REAL PETSC_QD_DD 1227 #else 1228 #define PETSC_REAL PETSC_DOUBLE 1229 #endif 1230 #define PETSC_FORTRANADDR PETSC_LONG 1231 1232 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*); 1233 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*); 1234 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscDataTypeGetSize(PetscDataType,size_t*); 1235 1236 /* 1237 Basic memory and string operations. These are usually simple wrappers 1238 around the basic Unix system calls, but a few of them have additional 1239 functionality and/or error checking. 1240 */ 1241 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscBitMemcpy(void*,PetscInt,const void*,PetscInt,PetscInt,PetscDataType); 1242 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscMemmove(void*,void *,size_t); 1243 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscMemcmp(const void*,const void*,size_t,PetscBool *); 1244 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrlen(const char[],size_t*); 1245 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrToArray(const char[],int*,char ***); 1246 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrToArrayDestroy(int,char **); 1247 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrcmp(const char[],const char[],PetscBool *); 1248 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrgrt(const char[],const char[],PetscBool *); 1249 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrcasecmp(const char[],const char[],PetscBool *); 1250 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrncmp(const char[],const char[],size_t,PetscBool *); 1251 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrcpy(char[],const char[]); 1252 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrcat(char[],const char[]); 1253 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrncat(char[],const char[],size_t); 1254 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrncpy(char[],const char[],size_t); 1255 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrchr(const char[],char,char *[]); 1256 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrtolower(char[]); 1257 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrrchr(const char[],char,char *[]); 1258 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrstr(const char[],const char[],char *[]); 1259 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrrstr(const char[],const char[],char *[]); 1260 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrallocpy(const char[],char *[]); 1261 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStrreplace(MPI_Comm,const char[],char[],size_t); 1262 1263 /*S 1264 PetscToken - 'Token' used for managing tokenizing strings 1265 1266 Level: intermediate 1267 1268 .seealso: PetscTokenCreate(), PetscTokenFind(), PetscTokenDestroy() 1269 S*/ 1270 typedef struct _p_PetscToken* PetscToken; 1271 1272 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscTokenCreate(const char[],const char,PetscToken*); 1273 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscTokenFind(PetscToken,char *[]); 1274 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscTokenDestroy(PetscToken); 1275 1276 /* 1277 These are MPI operations for MPI_Allreduce() etc 1278 */ 1279 EXTERN PETSCSYS_DLLEXPORT MPI_Op PetscMaxSum_Op; 1280 #if defined(PETSC_USE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX) 1281 EXTERN PETSCSYS_DLLEXPORT MPI_Op MPIU_SUM; 1282 #else 1283 #define MPIU_SUM MPI_SUM 1284 #endif 1285 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*); 1286 1287 /*S 1288 PetscObject - any PETSc object, PetscViewer, Mat, Vec, KSP etc 1289 1290 Level: beginner 1291 1292 Note: This is the base class from which all objects appear. 1293 1294 .seealso: PetscObjectDestroy(), PetscObjectView(), PetscObjectGetName(), PetscObjectSetName(), PetscObjectReference(), PetscObjectDereferenc() 1295 S*/ 1296 typedef struct _p_PetscObject* PetscObject; 1297 1298 /*S 1299 PetscFList - Linked list of functions, possibly stored in dynamic libraries, accessed 1300 by string name 1301 1302 Level: advanced 1303 1304 .seealso: PetscFListAdd(), PetscFListDestroy() 1305 S*/ 1306 typedef struct _n_PetscFList *PetscFList; 1307 1308 /*E 1309 PetscFileMode - Access mode for a file. 1310 1311 Level: beginner 1312 1313 FILE_MODE_READ - open a file at its beginning for reading 1314 1315 FILE_MODE_WRITE - open a file at its beginning for writing (will create if the file does not exist) 1316 1317 FILE_MODE_APPEND - open a file at end for writing 1318 1319 FILE_MODE_UPDATE - open a file for updating, meaning for reading and writing 1320 1321 FILE_MODE_APPEND_UPDATE - open a file for updating, meaning for reading and writing, at the end 1322 1323 .seealso: PetscViewerFileSetMode() 1324 E*/ 1325 typedef enum {FILE_MODE_READ, FILE_MODE_WRITE, FILE_MODE_APPEND, FILE_MODE_UPDATE, FILE_MODE_APPEND_UPDATE} PetscFileMode; 1326 1327 #include "petscviewer.h" 1328 #include "petscoptions.h" 1329 1330 #define PETSC_SMALLEST_CLASSID 1211211 1331 extern PETSCSYS_DLLEXPORT PetscClassId PETSC_LARGEST_CLASSID; 1332 extern PETSCSYS_DLLEXPORT PetscClassId PETSC_OBJECT_CLASSID; 1333 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscClassIdRegister(const char[],PetscClassId *); 1334 1335 /* 1336 Routines that get memory usage information from the OS 1337 */ 1338 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscMemoryGetCurrentUsage(PetscLogDouble *); 1339 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscMemoryGetMaximumUsage(PetscLogDouble *); 1340 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscMemorySetGetMaximumUsage(void); 1341 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscMemoryShowUsage(PetscViewer,const char[]); 1342 1343 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscInfoAllow(PetscBool ,const char []); 1344 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetTime(PetscLogDouble*); 1345 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetCPUTime(PetscLogDouble*); 1346 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSleep(PetscReal); 1347 1348 /* 1349 Initialization of PETSc 1350 */ 1351 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscInitialize(int*,char***,const char[],const char[]); 1352 PetscPolymorphicSubroutine(PetscInitialize,(int *argc,char ***args),(argc,args,PETSC_NULL,PETSC_NULL)) 1353 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscInitializeNoArguments(void); 1354 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscInitialized(PetscBool *); 1355 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFinalized(PetscBool *); 1356 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFinalize(void); 1357 EXTERN PetscErrorCode PetscInitializeFortran(void); 1358 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetArgs(int*,char ***); 1359 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetArguments(char ***); 1360 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFreeArguments(char **); 1361 1362 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscEnd(void); 1363 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSysInitializePackage(const char[]); 1364 1365 extern MPI_Comm PETSC_COMM_LOCAL_WORLD; 1366 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscOpenMPMerge(PetscMPIInt,PetscErrorCode (*)(void*),void*); 1367 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscOpenMPSpawn(PetscMPIInt); 1368 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscOpenMPFinalize(void); 1369 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscOpenMPRun(MPI_Comm,PetscErrorCode (*)(MPI_Comm,void *),void*); 1370 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscOpenMPRunCtx(MPI_Comm,PetscErrorCode (*)(MPI_Comm,void*,void *),void*); 1371 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscOpenMPFree(MPI_Comm,void*); 1372 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscOpenMPMalloc(MPI_Comm,size_t,void**); 1373 1374 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscPythonInitialize(const char[],const char[]); 1375 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscPythonFinalize(void); 1376 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscPythonPrintError(void); 1377 1378 /* 1379 These are so that in extern C code we can caste function pointers to non-extern C 1380 function pointers. Since the regular C++ code expects its function pointers to be 1381 C++. 1382 */ 1383 typedef void (**PetscVoidStarFunction)(void); 1384 typedef void (*PetscVoidFunction)(void); 1385 typedef PetscErrorCode (*PetscErrorCodeFunction)(void); 1386 1387 /* 1388 PetscTryMethod - Queries an object for a method, if it exists then calls it. 1389 These are intended to be used only inside PETSc functions. 1390 1391 Level: developer 1392 1393 .seealso: PetscUseMethod() 1394 */ 1395 #define PetscTryMethod(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 } 1400 1401 /* 1402 PetscUseMethod - Queries an object for a method, if it exists then calls it, otherwise generates an error. 1403 These are intended to be used only inside PETSc functions. 1404 1405 Level: developer 1406 1407 .seealso: PetscTryMethod() 1408 */ 1409 #define PetscUseMethod(obj,A,B,C) \ 1410 0;{ PetscErrorCode (*f)B, __ierr; \ 1411 __ierr = PetscObjectQueryFunction((PetscObject)obj,A,(PetscVoidStarFunction)&f);CHKERRQ(__ierr); \ 1412 if (f) {__ierr = (*f)C;CHKERRQ(__ierr);}\ 1413 else SETERRQ1(((PetscObject)obj)->comm,PETSC_ERR_SUP,"Cannot locate function %s in object",A); \ 1414 } 1415 1416 /* 1417 Functions that can act on any PETSc object. 1418 */ 1419 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectCreate(MPI_Comm,PetscObject*); 1420 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectCreateGeneric(MPI_Comm, PetscClassId, const char [], PetscObject *); 1421 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectDestroy(PetscObject); 1422 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectGetComm(PetscObject,MPI_Comm *); 1423 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectGetClassId(PetscObject,PetscClassId *); 1424 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectSetType(PetscObject,const char []); 1425 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectSetPrecision(PetscObject,PetscPrecision); 1426 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectGetType(PetscObject,const char *[]); 1427 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectSetName(PetscObject,const char[]); 1428 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectGetName(PetscObject,const char*[]); 1429 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer,const char[]); 1430 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectSetTabLevel(PetscObject,PetscInt); 1431 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectGetTabLevel(PetscObject,PetscInt*); 1432 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt); 1433 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectReference(PetscObject); 1434 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectGetReference(PetscObject,PetscInt*); 1435 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectDereference(PetscObject); 1436 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectGetNewTag(PetscObject,PetscMPIInt *); 1437 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectView(PetscObject,PetscViewer); 1438 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectCompose(PetscObject,const char[],PetscObject); 1439 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectQuery(PetscObject,const char[],PetscObject *); 1440 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectComposeFunction(PetscObject,const char[],const char[],void (*)(void)); 1441 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectSetFromOptions(PetscObject); 1442 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectSetUp(PetscObject); 1443 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscCommGetNewTag(MPI_Comm,PetscMPIInt *); 1444 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectAddOptionsHandler(PetscObject,PetscErrorCode (*)(PetscObject,void*),PetscErrorCode (*)(PetscObject,void*),void*); 1445 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectProcessOptionsHandlers(PetscObject); 1446 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectDestroyOptionsHandlers(PetscObject); 1447 1448 /*MC 1449 PetscObjectComposeFunctionDynamic - Associates a function with a given PETSc object. 1450 1451 Synopsis: 1452 PetscErrorCode PetscObjectComposeFunctionDynamic(PetscObject obj,const char name[],const char fname[],void *ptr) 1453 1454 Logically Collective on PetscObject 1455 1456 Input Parameters: 1457 + obj - the PETSc object; this must be cast with a (PetscObject), for example, 1458 PetscObjectCompose((PetscObject)mat,...); 1459 . name - name associated with the child function 1460 . fname - name of the function 1461 - ptr - function pointer (or PETSC_NULL if using dynamic libraries) 1462 1463 Level: advanced 1464 1465 1466 Notes: 1467 To remove a registered routine, pass in a PETSC_NULL rname and fnc(). 1468 1469 PetscObjectComposeFunctionDynamic() can be used with any PETSc object (such as 1470 Mat, Vec, KSP, SNES, etc.) or any user-provided object. 1471 1472 The composed function must be wrapped in a EXTERN_C_BEGIN/END for this to 1473 work in C++/complex with dynamic link libraries (./configure options --with-shared-libraries --with-dynamic-loading) 1474 enabled. 1475 1476 Concepts: objects^composing functions 1477 Concepts: composing functions 1478 Concepts: functions^querying 1479 Concepts: objects^querying 1480 Concepts: querying objects 1481 1482 .seealso: PetscObjectQueryFunction() 1483 M*/ 1484 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1485 #define PetscObjectComposeFunctionDynamic(a,b,c,d) PetscObjectComposeFunction(a,b,c,0) 1486 #else 1487 #define PetscObjectComposeFunctionDynamic(a,b,c,d) PetscObjectComposeFunction(a,b,c,(PetscVoidFunction)(d)) 1488 #endif 1489 1490 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectQueryFunction(PetscObject,const char[],void (**)(void)); 1491 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectSetOptionsPrefix(PetscObject,const char[]); 1492 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectAppendOptionsPrefix(PetscObject,const char[]); 1493 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectPrependOptionsPrefix(PetscObject,const char[]); 1494 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectGetOptionsPrefix(PetscObject,const char*[]); 1495 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectAMSPublish(PetscObject); 1496 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectUnPublish(PetscObject); 1497 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectChangeTypeName(PetscObject,const char[]); 1498 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectRegisterDestroy(PetscObject); 1499 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectRegisterDestroyAll(void); 1500 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscObjectName(PetscObject); 1501 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscTypeCompare(PetscObject,const char[],PetscBool *); 1502 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRegisterFinalize(PetscErrorCode (*)(void)); 1503 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRegisterFinalizeAll(void); 1504 1505 /* 1506 Defines PETSc error handling. 1507 */ 1508 #include "petscerror.h" 1509 1510 /*S 1511 PetscOList - Linked list of PETSc objects, each accessable by string name 1512 1513 Level: developer 1514 1515 Notes: Used by PetscObjectCompose() and PetscObjectQuery() 1516 1517 .seealso: PetscOListAdd(), PetscOListDestroy(), PetscOListFind(), PetscObjectCompose(), PetscObjectQuery() 1518 S*/ 1519 typedef struct _n_PetscOList *PetscOList; 1520 1521 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscOListDestroy(PetscOList); 1522 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscOListFind(PetscOList,const char[],PetscObject*); 1523 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscOListReverseFind(PetscOList,PetscObject,char**); 1524 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscOListAdd(PetscOList *,const char[],PetscObject); 1525 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscOListDuplicate(PetscOList,PetscOList *); 1526 1527 /* 1528 Dynamic library lists. Lists of names of routines in objects or in dynamic 1529 link libraries that will be loaded as needed. 1530 */ 1531 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFListAdd(PetscFList*,const char[],const char[],void (*)(void)); 1532 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFListDestroy(PetscFList*); 1533 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFListFind(PetscFList,MPI_Comm,const char[],void (**)(void)); 1534 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFList,const char[]); 1535 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1536 #define PetscFListAddDynamic(a,b,p,c) PetscFListAdd(a,b,p,0) 1537 #else 1538 #define PetscFListAddDynamic(a,b,p,c) PetscFListAdd(a,b,p,(void (*)(void))c) 1539 #endif 1540 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFListDuplicate(PetscFList,PetscFList *); 1541 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFListView(PetscFList,PetscViewer); 1542 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFListConcat(const char [],const char [],char []); 1543 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFListGet(PetscFList,char ***,int*); 1544 1545 /*S 1546 PetscDLLibrary - Linked list of dynamics libraries to search for functions 1547 1548 Level: advanced 1549 1550 --with-shared-libraries --with-dynamic-loading must be used with ./configure to use dynamic libraries 1551 1552 .seealso: PetscDLLibraryOpen() 1553 S*/ 1554 typedef struct _n_PetscDLLibrary *PetscDLLibrary; 1555 extern PETSCSYS_DLLEXPORT PetscDLLibrary DLLibrariesLoaded; 1556 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]); 1557 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]); 1558 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **); 1559 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscDLLibraryPrintPath(PetscDLLibrary); 1560 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool *); 1561 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *); 1562 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscDLLibraryClose(PetscDLLibrary); 1563 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscDLLibraryCCAAppend(MPI_Comm,PetscDLLibrary *,const char[]); 1564 1565 /* 1566 PetscFwk support. Needs to be documented. 1567 Logically it is an extension of PetscDLLXXX, PetscObjectCompose, etc. 1568 */ 1569 #include "petscfwk.h" 1570 1571 /* 1572 Useful utility routines 1573 */ 1574 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*); 1575 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*); 1576 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt); 1577 PetscPolymorphicSubroutine(PetscSequentialPhaseBegin,(MPI_Comm comm),(comm,1)) 1578 PetscPolymorphicSubroutine(PetscSequentialPhaseBegin,(void),(PETSC_COMM_WORLD,1)) 1579 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt); 1580 PetscPolymorphicSubroutine(PetscSequentialPhaseEnd,(MPI_Comm comm),(comm,1)) 1581 PetscPolymorphicSubroutine(PetscSequentialPhaseEnd,(void),(PETSC_COMM_WORLD,1)) 1582 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscBarrier(PetscObject); 1583 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscMPIDump(FILE*); 1584 1585 /* 1586 PetscNot - negates a logical type value and returns result as a PetscBool 1587 1588 Notes: This is useful in cases like 1589 $ int *a; 1590 $ PetscBool flag = PetscNot(a) 1591 where !a does not return a PetscBool because we cannot provide a cast from int to PetscBool in C. 1592 */ 1593 #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE) 1594 1595 /* 1596 Defines basic graphics available from PETSc. 1597 */ 1598 #include "petscdraw.h" 1599 1600 /* 1601 Defines the base data structures for all PETSc objects 1602 */ 1603 #include "private/petscimpl.h" 1604 1605 /* 1606 Defines PETSc profiling. 1607 */ 1608 #include "petsclog.h" 1609 1610 /* 1611 For locking, unlocking and destroying AMS memories associated with PETSc objects. ams.h is included in petscviewer.h 1612 */ 1613 #if defined(PETSC_HAVE_AMS) 1614 extern PetscBool PetscAMSPublishAll; 1615 #define PetscObjectTakeAccess(obj) ((((PetscObject)(obj))->amem == -1) ? 0 : AMS_Memory_take_access(((PetscObject)(obj))->amem)) 1616 #define PetscObjectGrantAccess(obj) ((((PetscObject)(obj))->amem == -1) ? 0 : AMS_Memory_grant_access(((PetscObject)(obj))->amem)) 1617 #define PetscObjectDepublish(obj) ((((PetscObject)(obj))->amem == -1) ? 0 : AMS_Memory_destroy(((PetscObject)(obj))->amem));((PetscObject)(obj))->amem = -1; 1618 #else 1619 #define PetscObjectTakeAccess(obj) 0 1620 #define PetscObjectGrantAccess(obj) 0 1621 #define PetscObjectDepublish(obj) 0 1622 #endif 1623 1624 /* 1625 Simple PETSc parallel IO for ASCII printing 1626 */ 1627 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFixFilename(const char[],char[]); 1628 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFOpen(MPI_Comm,const char[],const char[],FILE**); 1629 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFClose(MPI_Comm,FILE*); 1630 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFPrintf(MPI_Comm,FILE*,const char[],...); 1631 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscPrintf(MPI_Comm,const char[],...); 1632 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSNPrintf(char*,size_t,const char [],...); 1633 1634 1635 1636 /* These are used internally by PETSc ASCII IO routines*/ 1637 #include <stdarg.h> 1638 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list); 1639 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT (*PetscVFPrintf)(FILE*,const char[],va_list); 1640 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscVFPrintfDefault(FILE*,const char[],va_list); 1641 1642 #if defined(PETSC_HAVE_MATLAB_ENGINE) 1643 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscVFPrintf_Matlab(FILE*,const char[],va_list); 1644 #endif 1645 1646 /*MC 1647 PetscErrorPrintf - Prints error messages. 1648 1649 Synopsis: 1650 PetscErrorCode (*PetscErrorPrintf)(const char format[],...); 1651 1652 Not Collective 1653 1654 Input Parameters: 1655 . format - the usual printf() format string 1656 1657 Options Database Keys: 1658 + -error_output_stdout - cause error messages to be printed to stdout instead of the 1659 (default) stderr 1660 - -error_output_none to turn off all printing of error messages (does not change the way the 1661 error is handled.) 1662 1663 Notes: Use 1664 $ PetscErrorPrintf = PetscErrorPrintfNone; to turn off all printing of error messages (does not change the way the 1665 $ error is handled.) and 1666 $ PetscErrorPrintf = PetscErrorPrintfDefault; to turn it back on 1667 $ of you can use your own function 1668 1669 Use 1670 PETSC_STDERR = FILE* obtained from a file open etc. to have stderr printed to the file. 1671 PETSC_STDOUT = FILE* obtained from a file open etc. to have stdout printed to the file. 1672 1673 Use 1674 PetscPushErrorHandler() to provide your own error handler that determines what kind of messages to print 1675 1676 Level: developer 1677 1678 Fortran Note: 1679 This routine is not supported in Fortran. 1680 1681 Concepts: error messages^printing 1682 Concepts: printing^error messages 1683 1684 .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscHelpPrintf(), PetscPrintf(), PetscErrorHandlerPush(), PetscVFPrintf(), PetscHelpPrintf() 1685 M*/ 1686 EXTERN PETSCSYS_DLLEXPORT PetscErrorCode (*PetscErrorPrintf)(const char[],...); 1687 1688 /*MC 1689 PetscHelpPrintf - Prints help messages. 1690 1691 Synopsis: 1692 PetscErrorCode (*PetscHelpPrintf)(const char format[],...); 1693 1694 Not Collective 1695 1696 Input Parameters: 1697 . format - the usual printf() format string 1698 1699 Level: developer 1700 1701 Fortran Note: 1702 This routine is not supported in Fortran. 1703 1704 Concepts: help messages^printing 1705 Concepts: printing^help messages 1706 1707 .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf() 1708 M*/ 1709 EXTERN PETSCSYS_DLLEXPORT PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...); 1710 1711 EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...); 1712 EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...); 1713 EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...); 1714 1715 #if defined(PETSC_HAVE_POPEN) 1716 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **); 1717 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscPClose(MPI_Comm,FILE*); 1718 #endif 1719 1720 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSynchronizedPrintf(MPI_Comm,const char[],...); 1721 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...); 1722 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSynchronizedFlush(MPI_Comm); 1723 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]); 1724 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**); 1725 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStartJava(MPI_Comm,const char[],const char[],FILE**); 1726 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetPetscDir(const char*[]); 1727 1728 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscPopUpSelect(MPI_Comm,const char*,const char*,int,const char**,int*); 1729 1730 /*S 1731 PetscContainer - Simple PETSc object that contains a pointer to any required data 1732 1733 Level: advanced 1734 1735 .seealso: PetscObject, PetscContainerCreate() 1736 S*/ 1737 extern PetscClassId PETSCSYS_DLLEXPORT PETSC_CONTAINER_CLASSID; 1738 typedef struct _p_PetscContainer* PetscContainer; 1739 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscContainerGetPointer(PetscContainer,void **); 1740 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscContainerSetPointer(PetscContainer,void *); 1741 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscContainerDestroy(PetscContainer); 1742 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscContainerCreate(MPI_Comm,PetscContainer *); 1743 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void*)); 1744 1745 /* 1746 For use in debuggers 1747 */ 1748 extern PETSCSYS_DLLEXPORT PetscMPIInt PetscGlobalRank; 1749 extern PETSCSYS_DLLEXPORT PetscMPIInt PetscGlobalSize; 1750 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscIntView(PetscInt,const PetscInt[],PetscViewer); 1751 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRealView(PetscInt,const PetscReal[],PetscViewer); 1752 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscScalarView(PetscInt,const PetscScalar[],PetscViewer); 1753 1754 #if defined(PETSC_HAVE_MEMORY_H) 1755 #include <memory.h> 1756 #endif 1757 #if defined(PETSC_HAVE_STDLIB_H) 1758 #include <stdlib.h> 1759 #endif 1760 #if defined(PETSC_HAVE_STRINGS_H) 1761 #include <strings.h> 1762 #endif 1763 #if defined(PETSC_HAVE_STRING_H) 1764 #include <string.h> 1765 #endif 1766 1767 1768 #if defined(PETSC_HAVE_XMMINTRIN_H) 1769 #include <xmmintrin.h> 1770 #endif 1771 #if defined(PETSC_HAVE_STDINT_H) 1772 #include <stdint.h> 1773 #endif 1774 1775 /*@C 1776 PetscMemcpy - Copies n bytes, beginning at location b, to the space 1777 beginning at location a. The two memory regions CANNOT overlap, use 1778 PetscMemmove() in that case. 1779 1780 Not Collective 1781 1782 Input Parameters: 1783 + b - pointer to initial memory space 1784 - n - length (in bytes) of space to copy 1785 1786 Output Parameter: 1787 . a - pointer to copy space 1788 1789 Level: intermediate 1790 1791 Compile Option: 1792 PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used 1793 for memory copies on double precision values. 1794 PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used 1795 for memory copies on double precision values. 1796 PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used 1797 for memory copies on double precision values. 1798 1799 Note: 1800 This routine is analogous to memcpy(). 1801 1802 Developer Note: this is inlined for fastest performance 1803 1804 Concepts: memory^copying 1805 Concepts: copying^memory 1806 1807 .seealso: PetscMemmove() 1808 1809 @*/ 1810 PETSC_STATIC_INLINE PetscErrorCode PETSCSYS_DLLEXPORT PetscMemcpy(void *a,const void *b,size_t n) 1811 { 1812 #if defined(PETSC_USE_DEBUG) 1813 unsigned long al = (unsigned long) a,bl = (unsigned long) b; 1814 unsigned long nl = (unsigned long) n; 1815 if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer"); 1816 if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer"); 1817 #endif 1818 PetscFunctionBegin; 1819 if (a != b) { 1820 #if defined(PETSC_USE_DEBUG) 1821 if ((al > bl && (al - bl) < nl) || (bl - al) < nl) { 1822 SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\ 1823 or make sure your copy regions and lengths are correct. \n\ 1824 Length (bytes) %ld first address %ld second address %ld",nl,al,bl); 1825 } 1826 #endif 1827 #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY)) 1828 if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) { 1829 size_t len = n/sizeof(PetscScalar); 1830 #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) 1831 PetscBLASInt one = 1,blen = PetscBLASIntCast(len); 1832 BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one); 1833 #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY) 1834 fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a); 1835 #else 1836 size_t i; 1837 PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a; 1838 for (i=0; i<len; i++) y[i] = x[i]; 1839 #endif 1840 } else { 1841 memcpy((char*)(a),(char*)(b),n); 1842 } 1843 #elif defined(PETSC_HAVE__INTEL_FAST_MEMCPY) 1844 _intel_fast_memcpy((char*)(a),(char*)(b),n); 1845 #else 1846 memcpy((char*)(a),(char*)(b),n); 1847 #endif 1848 } 1849 PetscFunctionReturn(0); 1850 } 1851 1852 /*@C 1853 PetscMemzero - Zeros the specified memory. 1854 1855 Not Collective 1856 1857 Input Parameters: 1858 + a - pointer to beginning memory location 1859 - n - length (in bytes) of memory to initialize 1860 1861 Level: intermediate 1862 1863 Compile Option: 1864 PETSC_PREFER_BZERO - on certain machines (the IBM RS6000) the bzero() routine happens 1865 to be faster than the memset() routine. This flag causes the bzero() routine to be used. 1866 1867 Developer Note: this is inlined for fastest performance 1868 1869 Concepts: memory^zeroing 1870 Concepts: zeroing^memory 1871 1872 .seealso: PetscMemcpy() 1873 @*/ 1874 PETSC_STATIC_INLINE PetscErrorCode PETSCSYS_DLLEXPORT PetscMemzero(void *a,size_t n) 1875 { 1876 if (n > 0) { 1877 #if defined(PETSC_USE_DEBUG) 1878 if (!a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer"); 1879 #endif 1880 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) 1881 if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) { 1882 size_t i,len = n/sizeof(PetscScalar); 1883 PetscScalar *x = (PetscScalar*)a; 1884 for (i=0; i<len; i++) x[i] = 0.0; 1885 } else { 1886 #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO) 1887 if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) { 1888 PetscInt len = n/sizeof(PetscScalar); 1889 fortranzero_(&len,(PetscScalar*)a); 1890 } else { 1891 #endif 1892 #if defined(PETSC_PREFER_BZERO) 1893 bzero((char *)a,n); 1894 #elif defined (PETSC_HAVE__INTEL_FAST_MEMSET) 1895 _intel_fast_memset((char*)a,0,n); 1896 #else 1897 memset((char*)a,0,n); 1898 #endif 1899 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO) 1900 } 1901 #endif 1902 } 1903 return 0; 1904 } 1905 1906 /*MC 1907 PetscPrefetchBlock - Prefetches a block of memory 1908 1909 Synopsis: 1910 void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t) 1911 1912 Not Collective 1913 1914 Input Parameters: 1915 + a - pointer to first element to fetch (any type but usually PetscInt or PetscScalar) 1916 . n - number of elements to fetch 1917 . rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors) 1918 - t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note 1919 1920 Level: developer 1921 1922 Notes: 1923 The last two arguments (rw and t) must be compile-time constants. 1924 1925 Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality. Not all architectures offer 1926 equivalent locality hints, but the following macros are always defined to their closest analogue. 1927 + PETSC_PREFETCH_HINT_NTA - Non-temporal. Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched). 1928 . 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. 1929 . PETSC_PREFETCH_HINT_T1 - Fetch to level 2 and higher (not L1). 1930 - PETSC_PREFETCH_HINT_T2 - Fetch to high-level cache only. (On many systems, T0 and T1 are equivalent.) 1931 1932 This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid 1933 address). 1934 1935 Concepts: memory 1936 M*/ 1937 #define PetscPrefetchBlock(a,n,rw,t) do { \ 1938 const char *_p = (const char*)(a),*_end = (const char*)((a)+(n)); \ 1939 for ( ; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \ 1940 } while (0) 1941 1942 /* 1943 Allows accessing Matlab Engine 1944 */ 1945 #include "petscmatlab.h" 1946 1947 /* 1948 Determine if some of the kernel computation routines use 1949 Fortran (rather than C) for the numerical calculations. On some machines 1950 and compilers (like complex numbers) the Fortran version of the routines 1951 is faster than the C/C++ versions. The flag --with-fortran-kernels 1952 should be used with ./configure to turn these on. 1953 */ 1954 #if defined(PETSC_USE_FORTRAN_KERNELS) 1955 1956 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL) 1957 #define PETSC_USE_FORTRAN_KERNEL_MULTCRL 1958 #endif 1959 1960 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM) 1961 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM 1962 #endif 1963 1964 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ) 1965 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ 1966 #endif 1967 1968 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ) 1969 #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ 1970 #endif 1971 1972 #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM) 1973 #define PETSC_USE_FORTRAN_KERNEL_NORM 1974 #endif 1975 1976 #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY) 1977 #define PETSC_USE_FORTRAN_KERNEL_MAXPY 1978 #endif 1979 1980 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 1981 #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ 1982 #endif 1983 1984 #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ) 1985 #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ 1986 #endif 1987 1988 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ) 1989 #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ 1990 #endif 1991 1992 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ) 1993 #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ 1994 #endif 1995 1996 #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT) 1997 #define PETSC_USE_FORTRAN_KERNEL_MDOT 1998 #endif 1999 2000 #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY) 2001 #define PETSC_USE_FORTRAN_KERNEL_XTIMESY 2002 #endif 2003 2004 #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX) 2005 #define PETSC_USE_FORTRAN_KERNEL_AYPX 2006 #endif 2007 2008 #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY) 2009 #define PETSC_USE_FORTRAN_KERNEL_WAXPY 2010 #endif 2011 2012 #endif 2013 2014 /* 2015 Macros for indicating code that should be compiled with a C interface, 2016 rather than a C++ interface. Any routines that are dynamically loaded 2017 (such as the PCCreate_XXX() routines) must be wrapped so that the name 2018 mangler does not change the functions symbol name. This just hides the 2019 ugly extern "C" {} wrappers. 2020 */ 2021 #if defined(__cplusplus) 2022 #define EXTERN_C_BEGIN extern "C" { 2023 #define EXTERN_C_END } 2024 #else 2025 #define EXTERN_C_BEGIN 2026 #define EXTERN_C_END 2027 #endif 2028 2029 /* --------------------------------------------------------------------*/ 2030 2031 /*MC 2032 MPI_Comm - the basic object used by MPI to determine which processes are involved in a 2033 communication 2034 2035 Level: beginner 2036 2037 Note: This manual page is a place-holder because MPICH does not have a manual page for MPI_Comm 2038 2039 .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF 2040 M*/ 2041 2042 /*MC 2043 PetscScalar - PETSc type that represents either a double precision real number, a double precision 2044 complex number, a single precision real number, a long double or an int - if the code is configured 2045 with --with-scalar-type=real,complex --with-precision=single,double,longdouble,int,matsingle 2046 2047 2048 Level: beginner 2049 2050 .seealso: PetscReal, PassiveReal, PassiveScalar, MPIU_SCALAR, PetscInt 2051 M*/ 2052 2053 /*MC 2054 PetscReal - PETSc type that represents a real number version of PetscScalar 2055 2056 Level: beginner 2057 2058 .seealso: PetscScalar, PassiveReal, PassiveScalar 2059 M*/ 2060 2061 /*MC 2062 PassiveScalar - PETSc type that represents a PetscScalar 2063 Level: beginner 2064 2065 This is the same as a PetscScalar except in code that is automatically differentiated it is 2066 treated as a constant (not an indendent or dependent variable) 2067 2068 .seealso: PetscReal, PassiveReal, PetscScalar 2069 M*/ 2070 2071 /*MC 2072 PassiveReal - PETSc type that represents a PetscReal 2073 2074 Level: beginner 2075 2076 This is the same as a PetscReal except in code that is automatically differentiated it is 2077 treated as a constant (not an indendent or dependent variable) 2078 2079 .seealso: PetscScalar, PetscReal, PassiveScalar 2080 M*/ 2081 2082 /*MC 2083 MPIU_SCALAR - MPI datatype corresponding to PetscScalar 2084 2085 Level: beginner 2086 2087 Note: In MPI calls that require an MPI datatype that matches a PetscScalar or array of PetscScalars 2088 pass this value 2089 2090 .seealso: PetscReal, PassiveReal, PassiveScalar, PetscScalar, MPIU_INT 2091 M*/ 2092 2093 #if defined(PETSC_HAVE_MPIIO) 2094 #if !defined(PETSC_WORDS_BIGENDIAN) 2095 extern PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*); 2096 extern PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*); 2097 #else 2098 #define MPIU_File_write_all(a,b,c,d,e) MPI_File_write_all(a,b,c,d,e) 2099 #define MPIU_File_read_all(a,b,c,d,e) MPI_File_read_all(a,b,c,d,e) 2100 #endif 2101 #endif 2102 2103 /* the following petsc_static_inline require petscerror.h */ 2104 2105 /* Limit MPI to 32-bits */ 2106 #define PETSC_MPI_INT_MAX 2147483647 2107 #define PETSC_MPI_INT_MIN -2147483647 2108 /* Limit BLAS to 32-bits */ 2109 #define PETSC_BLAS_INT_MAX 2147483647 2110 #define PETSC_BLAS_INT_MIN -2147483647 2111 /* On 32 bit systems HDF5 is limited by size of integer, because hsize_t is defined as size_t */ 2112 #define PETSC_HDF5_INT_MAX 2147483647 2113 #define PETSC_HDF5_INT_MIN -2147483647 2114 2115 #if defined(PETSC_USE_64BIT_INDICES) 2116 #define PetscMPIIntCheck(a) if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Message too long for MPI") 2117 #define PetscBLASIntCheck(a) if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK") 2118 #define PetscMPIIntCast(a) (a);PetscMPIIntCheck(a) 2119 #define PetscBLASIntCast(a) (a);PetscBLASIntCheck(a) 2120 2121 #if (PETSC_SIZEOF_SIZE_T == 4) 2122 #define PetscHDF5IntCheck(a) if ((a) > PETSC_HDF5_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for HDF5") 2123 #define PetscHDF5IntCast(a) (a);PetscHDF5IntCheck(a) 2124 #else 2125 #define PetscHDF5IntCheck(a) 2126 #define PetscHDF5IntCast(a) a 2127 #endif 2128 2129 #else 2130 #define PetscMPIIntCheck(a) 2131 #define PetscBLASIntCheck(a) 2132 #define PetscHDF5IntCheck(a) 2133 #define PetscMPIIntCast(a) a 2134 #define PetscBLASIntCast(a) a 2135 #define PetscHDF5IntCast(a) a 2136 #endif 2137 2138 2139 /* 2140 The IBM include files define hz, here we hide it so that it may be used 2141 as a regular user variable. 2142 */ 2143 #if defined(hz) 2144 #undef hz 2145 #endif 2146 2147 /* For arrays that contain filenames or paths */ 2148 2149 2150 #if defined(PETSC_HAVE_LIMITS_H) 2151 #include <limits.h> 2152 #endif 2153 #if defined(PETSC_HAVE_SYS_PARAM_H) 2154 #include <sys/param.h> 2155 #endif 2156 #if defined(PETSC_HAVE_SYS_TYPES_H) 2157 #include <sys/types.h> 2158 #endif 2159 #if defined(MAXPATHLEN) 2160 # define PETSC_MAX_PATH_LEN MAXPATHLEN 2161 #elif defined(MAX_PATH) 2162 # define PETSC_MAX_PATH_LEN MAX_PATH 2163 #elif defined(_MAX_PATH) 2164 # define PETSC_MAX_PATH_LEN _MAX_PATH 2165 #else 2166 # define PETSC_MAX_PATH_LEN 4096 2167 #endif 2168 2169 /* Special support for C++ */ 2170 #include "petscsys.hh" 2171 2172 2173 /*MC 2174 2175 UsingFortran - Fortran can be used with PETSc in four distinct approaches 2176 2177 $ 1) classic Fortran 77 style 2178 $#include "finclude/petscXXX.h" to work with material from the XXX component of PETSc 2179 $ XXX variablename 2180 $ You cannot use this approach if you wish to use the Fortran 90 specific PETSc routines 2181 $ which end in F90; such as VecGetArrayF90() 2182 $ 2183 $ 2) classic Fortran 90 style 2184 $#include "finclude/petscXXX.h" 2185 $#include "finclude/petscXXX.h90" to work with material from the XXX component of PETSc 2186 $ XXX variablename 2187 $ 2188 $ 3) Using Fortran modules 2189 $#include "finclude/petscXXXdef.h" 2190 $ use petscXXXX 2191 $ XXX variablename 2192 $ 2193 $ 4) Use Fortran modules and Fortran data types for PETSc types 2194 $#include "finclude/petscXXXdef.h" 2195 $ use petscXXXX 2196 $ type(XXX) variablename 2197 $ To use this approach you must ./configure PETSc with the additional 2198 $ option --with-fortran-datatypes You cannot use the type(XXX) declaration approach without using Fortran modules 2199 2200 Finally if you absolutely do not want to use any #include you can use either 2201 2202 $ 3a) skip the #include BUT you cannot use any PETSc data type names like Vec, Mat, PetscInt, PetscErrorCode etc 2203 $ and you must declare the variables as integer, for example 2204 $ integer variablename 2205 $ 2206 $ 4a) skip the #include, you use the object types like type(Vec) type(Mat) but cannot use the data type 2207 $ names like PetscErrorCode, PetscInt etc. again for those you must use integer 2208 2209 We recommend either 2 or 3. Approaches 2 and 3 provide type checking for most PETSc function calls; 4 has type checking 2210 for only a few PETSc functions. 2211 2212 Fortran type checking with interfaces is strick, this means you cannot pass a scalar value when an array value 2213 is expected (even though it is legal Fortran). For example when setting a single value in a matrix with MatSetValues() 2214 you cannot have something like 2215 $ PetscInt row,col 2216 $ PetscScalar val 2217 $ ... 2218 $ call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr) 2219 You must instead have 2220 $ PetscInt row(1),col(1) 2221 $ PetscScalar val(1) 2222 $ ... 2223 $ call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr) 2224 2225 2226 See the example src/vec/vec/examples/tutorials/ex20f90.F90 for an example that can use all four approaches 2227 2228 Developer Notes: The finclude/petscXXXdef.h contain all the #defines (would be typedefs in C code) these 2229 automatically include their predecessors; for example finclude/petscvecdef.h includes finclude/petscisdef.h 2230 2231 The finclude/petscXXXX.h contain all the parameter statements for that package. These automatically include 2232 their finclude/petscXXXdef.h file but DO NOT automatically include their predecessors; for example 2233 finclude/petscvec.h does NOT automatically include finclude/petscis.h 2234 2235 The finclude/ftn-custom/petscXXXdef.h90 are not intended to be used directly in code, they define the 2236 Fortran data type type(XXX) (for example type(Vec)) when PETSc is ./configure with the --with-fortran-datatypes option. 2237 2238 The finclude/ftn-custom/petscXXX.h90 (not included directly by code) contain interface definitions for 2239 the PETSc Fortran stubs that have different bindings then their C version (for example VecGetArrayF90). 2240 2241 The finclude/ftn-auto/petscXXX.h90 (not included directly by code) contain interface definitions generated 2242 automatically by "make allfortranstubs". 2243 2244 The finclude/petscXXX.h90 includes the custom finclude/ftn-custom/petscXXX.h90 and if ./configure 2245 was run with --with-fortran-interfaces it also includes the finclude/ftn-auto/petscXXX.h90 These DO NOT automatically 2246 include their predecessors 2247 2248 Level: beginner 2249 2250 M*/ 2251 2252 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetArchType(char[],size_t); 2253 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetHostName(char[],size_t); 2254 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetUserName(char[],size_t); 2255 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetProgramName(char[],size_t); 2256 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSetProgramName(const char[]); 2257 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetDate(char[],size_t); 2258 2259 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSortInt(PetscInt,PetscInt[]); 2260 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSortRemoveDupsInt(PetscInt*,PetscInt[]); 2261 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]); 2262 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]); 2263 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]); 2264 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]); 2265 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]); 2266 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSortReal(PetscInt,PetscReal[]); 2267 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]); 2268 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]); 2269 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]); 2270 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**); 2271 2272 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSetDisplay(void); 2273 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetDisplay(char[],size_t); 2274 2275 /*E 2276 PetscRandomType - String with the name of a PETSc randomizer 2277 with an optional dynamic library name, for example 2278 http://www.mcs.anl.gov/petsc/lib.a:myrandcreate() 2279 2280 Level: beginner 2281 2282 Notes: to use the SPRNG you must have ./configure PETSc 2283 with the option --download-sprng 2284 2285 .seealso: PetscRandomSetType(), PetscRandom 2286 E*/ 2287 #define PetscRandomType char* 2288 #define PETSCRAND "rand" 2289 #define PETSCRAND48 "rand48" 2290 #define PETSCSPRNG "sprng" 2291 2292 /* Logging support */ 2293 extern PETSCSYS_DLLEXPORT PetscClassId PETSC_RANDOM_CLASSID; 2294 2295 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomInitializePackage(const char[]); 2296 2297 /*S 2298 PetscRandom - Abstract PETSc object that manages generating random numbers 2299 2300 Level: intermediate 2301 2302 Concepts: random numbers 2303 2304 .seealso: PetscRandomCreate(), PetscRandomGetValue(), PetscRandomType 2305 S*/ 2306 typedef struct _p_PetscRandom* PetscRandom; 2307 2308 /* Dynamic creation and loading functions */ 2309 extern PetscFList PetscRandomList; 2310 extern PetscBool PetscRandomRegisterAllCalled; 2311 2312 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomRegisterAll(const char []); 2313 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomRegister(const char[],const char[],const char[],PetscErrorCode (*)(PetscRandom)); 2314 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomRegisterDestroy(void); 2315 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomSetType(PetscRandom, const PetscRandomType); 2316 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomSetFromOptions(PetscRandom); 2317 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomGetType(PetscRandom, const PetscRandomType*); 2318 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomViewFromOptions(PetscRandom,char*); 2319 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomView(PetscRandom,PetscViewer); 2320 2321 /*MC 2322 PetscRandomRegisterDynamic - Adds a new PetscRandom component implementation 2323 2324 Synopsis: 2325 PetscErrorCode PetscRandomRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(PetscRandom)) 2326 2327 Not Collective 2328 2329 Input Parameters: 2330 + name - The name of a new user-defined creation routine 2331 . path - The path (either absolute or relative) of the library containing this routine 2332 . func_name - The name of routine to create method context 2333 - create_func - The creation routine itself 2334 2335 Notes: 2336 PetscRandomRegisterDynamic() may be called multiple times to add several user-defined randome number generators 2337 2338 If dynamic libraries are used, then the fourth input argument (routine_create) is ignored. 2339 2340 Sample usage: 2341 .vb 2342 PetscRandomRegisterDynamic("my_rand","/home/username/my_lib/lib/libO/solaris/libmy.a", "MyPetscRandomtorCreate", MyPetscRandomtorCreate); 2343 .ve 2344 2345 Then, your random type can be chosen with the procedural interface via 2346 .vb 2347 PetscRandomCreate(MPI_Comm, PetscRandom *); 2348 PetscRandomSetType(PetscRandom,"my_random_name"); 2349 .ve 2350 or at runtime via the option 2351 .vb 2352 -random_type my_random_name 2353 .ve 2354 2355 Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 2356 2357 For an example of the code needed to interface your own random number generator see 2358 src/sys/random/impls/rand/rand.c 2359 2360 Level: advanced 2361 2362 .keywords: PetscRandom, register 2363 .seealso: PetscRandomRegisterAll(), PetscRandomRegisterDestroy(), PetscRandomRegister() 2364 M*/ 2365 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 2366 #define PetscRandomRegisterDynamic(a,b,c,d) PetscRandomRegister(a,b,c,0) 2367 #else 2368 #define PetscRandomRegisterDynamic(a,b,c,d) PetscRandomRegister(a,b,c,d) 2369 #endif 2370 2371 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomCreate(MPI_Comm,PetscRandom*); 2372 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomGetValue(PetscRandom,PetscScalar*); 2373 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomGetValueReal(PetscRandom,PetscReal*); 2374 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*); 2375 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar); 2376 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomSetSeed(PetscRandom,unsigned long); 2377 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomGetSeed(PetscRandom,unsigned long *); 2378 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomSeed(PetscRandom); 2379 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscRandomDestroy(PetscRandom); 2380 2381 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetFullPath(const char[],char[],size_t); 2382 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetRelativePath(const char[],char[],size_t); 2383 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetWorkingDirectory(char[],size_t); 2384 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetRealPath(const char[],char[]); 2385 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetHomeDirectory(char[],size_t); 2386 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscTestFile(const char[],char,PetscBool *); 2387 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscTestDirectory(const char[],char,PetscBool *); 2388 2389 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscBinaryRead(int,void*,PetscInt,PetscDataType); 2390 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscDataType); 2391 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscBinarySynchronizedWrite(MPI_Comm,int,void*,PetscInt,PetscDataType,PetscBool ); 2392 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscBinaryWrite(int,void*,PetscInt,PetscDataType,PetscBool ); 2393 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscBinaryOpen(const char[],PetscFileMode,int *); 2394 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscBinaryClose(int); 2395 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSharedTmp(MPI_Comm,PetscBool *); 2396 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSharedWorkingDirectory(MPI_Comm,PetscBool *); 2397 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGetTmp(MPI_Comm,char[],size_t); 2398 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *); 2399 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *); 2400 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscOpenSocket(char*,int,int*); 2401 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscWebServe(MPI_Comm,int); 2402 2403 /* 2404 In binary files variables are stored using the following lengths, 2405 regardless of how they are stored in memory on any one particular 2406 machine. Use these rather then sizeof() in computing sizes for 2407 PetscBinarySeek(). 2408 */ 2409 #define PETSC_BINARY_INT_SIZE (32/8) 2410 #define PETSC_BINARY_FLOAT_SIZE (32/8) 2411 #define PETSC_BINARY_CHAR_SIZE (8/8) 2412 #define PETSC_BINARY_SHORT_SIZE (16/8) 2413 #define PETSC_BINARY_DOUBLE_SIZE (64/8) 2414 #define PETSC_BINARY_SCALAR_SIZE sizeof(PetscScalar) 2415 2416 /*E 2417 PetscBinarySeekType - argument to PetscBinarySeek() 2418 2419 Level: advanced 2420 2421 .seealso: PetscBinarySeek(), PetscBinarySynchronizedSeek() 2422 E*/ 2423 typedef enum {PETSC_BINARY_SEEK_SET = 0,PETSC_BINARY_SEEK_CUR = 1,PETSC_BINARY_SEEK_END = 2} PetscBinarySeekType; 2424 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*); 2425 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*); 2426 2427 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSetDebugTerminal(const char[]); 2428 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSetDebugger(const char[],PetscBool ); 2429 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSetDefaultDebugger(void); 2430 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSetDebuggerFromString(char*); 2431 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscAttachDebugger(void); 2432 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscStopForDebugger(void); 2433 2434 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*); 2435 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**); 2436 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**); 2437 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**); 2438 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**); 2439 2440 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSSEIsEnabled(MPI_Comm,PetscBool *,PetscBool *); 2441 2442 /*E 2443 InsertMode - Whether entries are inserted or added into vectors or matrices 2444 2445 Level: beginner 2446 2447 .seealso: VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(), 2448 VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), 2449 MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd() 2450 E*/ 2451 typedef enum {NOT_SET_VALUES, INSERT_VALUES, ADD_VALUES, MAX_VALUES} InsertMode; 2452 2453 /*MC 2454 INSERT_VALUES - Put a value into a vector or matrix, overwrites any previous value 2455 2456 Level: beginner 2457 2458 .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(), 2459 VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), ADD_VALUES, 2460 MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES 2461 2462 M*/ 2463 2464 /*MC 2465 ADD_VALUES - Adds a value into a vector or matrix, if there previously was no value, just puts the 2466 value into that location 2467 2468 Level: beginner 2469 2470 .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(), 2471 VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), INSERT_VALUES, 2472 MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES 2473 2474 M*/ 2475 2476 /*MC 2477 MAX_VALUES - Puts the maximum of the scattered/gathered value and the current value into each location 2478 2479 Level: beginner 2480 2481 .seealso: InsertMode, VecScatterBegin(), VecScatterEnd(), ADD_VALUES, INSERT_VALUES 2482 2483 M*/ 2484 2485 /*S 2486 PetscSubcomm - Context of MPI subcommunicators, used by PCREDUNDANT 2487 2488 Level: advanced 2489 2490 Concepts: communicator, create 2491 S*/ 2492 typedef struct _n_PetscSubcomm* PetscSubcomm; 2493 2494 struct _n_PetscSubcomm { 2495 MPI_Comm parent; /* parent communicator */ 2496 MPI_Comm dupparent; /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */ 2497 MPI_Comm comm; /* this communicator */ 2498 PetscInt n; /* num of subcommunicators under the parent communicator */ 2499 PetscInt color; /* color of processors belong to this communicator */ 2500 }; 2501 2502 typedef enum {PETSC_SUBCOMM_GENERAL=0,PETSC_SUBCOMM_CONTIGUOUS=1,PETSC_SUBCOMM_INTERLACED=2} PetscSubcommType; 2503 extern const char *PetscSubcommTypes[]; 2504 2505 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSubcommCreate(MPI_Comm,PetscSubcomm*); 2506 EXTERN PetscErrorCode PETSCSYS_DLLEXPORT PetscSubcommDestroy(PetscSubcomm); 2507 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommSetNumber(PetscSubcomm,PetscInt); 2508 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommSetType(PetscSubcomm,const PetscSubcommType); 2509 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt,PetscMPIInt); 2510 2511 PETSC_EXTERN_CXX_END 2512 2513 #endif 2514