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