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