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