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