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) ((*(r2) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(PETSC_MEMALIGN-1),r1)) || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),0) \ 627 || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0)) 628 #endif 629 630 /*MC 631 PetscCalloc2 - Allocates 2 cleared (zeroed) arrays of memory both aligned to PETSC_MEMALIGN 632 633 Synopsis: 634 #include "petscsys.h" 635 PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2) 636 637 Not Collective 638 639 Input Parameter: 640 + m1 - number of elements to allocate in 1st chunk (may be zero) 641 - m2 - number of elements to allocate in 2nd chunk (may be zero) 642 643 Output Parameter: 644 + r1 - memory allocated in first chunk 645 - r2 - memory allocated in second chunk 646 647 Level: developer 648 649 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc1(), PetscMalloc2() 650 651 Concepts: memory allocation 652 M*/ 653 #define PetscCalloc2(m1,r1,m2,r2) (PetscMalloc2((m1),(r1),(m2),(r2)) || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2)))) 654 655 /*MC 656 PetscMalloc3 - Allocates 3 arrays of memory, all aligned to PETSC_MEMALIGN 657 658 Synopsis: 659 #include "petscsys.h" 660 PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3) 661 662 Not Collective 663 664 Input Parameter: 665 + m1 - number of elements to allocate in 1st chunk (may be zero) 666 . m2 - number of elements to allocate in 2nd chunk (may be zero) 667 - m3 - number of elements to allocate in 3rd chunk (may be zero) 668 669 Output Parameter: 670 + r1 - memory allocated in first chunk 671 . r2 - memory allocated in second chunk 672 - r3 - memory allocated in third chunk 673 674 Level: developer 675 676 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc3(), PetscFree3() 677 678 Concepts: memory allocation 679 680 M*/ 681 #if defined(PETSC_USE_DEBUG) 682 #define PetscMalloc3(m1,r1,m2,r2,m3,r3) (PetscMalloc((m1)*sizeof(**(r1)),r1) || PetscMalloc((m2)*sizeof(**(r2)),r2) || PetscMalloc((m3)*sizeof(**(r3)),r3)) 683 #else 684 #define PetscMalloc3(m1,r1,m2,r2,m3,r3) ((*(r2) = 0,*(r3) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(m3)*sizeof(**(r3))+2*(PETSC_MEMALIGN-1),r1)) \ 685 || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),0) \ 686 || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0)) 687 #endif 688 689 /*MC 690 PetscCalloc3 - Allocates 3 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN 691 692 Synopsis: 693 #include "petscsys.h" 694 PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3) 695 696 Not Collective 697 698 Input Parameter: 699 + m1 - number of elements to allocate in 1st chunk (may be zero) 700 . m2 - number of elements to allocate in 2nd chunk (may be zero) 701 - m3 - number of elements to allocate in 3rd chunk (may be zero) 702 703 Output Parameter: 704 + r1 - memory allocated in first chunk 705 . r2 - memory allocated in second chunk 706 - r3 - memory allocated in third chunk 707 708 Level: developer 709 710 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc2(), PetscMalloc3(), PetscFree3() 711 712 Concepts: memory allocation 713 M*/ 714 #define PetscCalloc3(m1,r1,m2,r2,m3,r3) \ 715 (PetscMalloc3((m1),(r1),(m2),(r2),(m3),(r3)) \ 716 || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3)))) 717 718 /*MC 719 PetscMalloc4 - Allocates 4 arrays of memory, all aligned to PETSC_MEMALIGN 720 721 Synopsis: 722 #include "petscsys.h" 723 PetscErrorCode PetscMalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4) 724 725 Not Collective 726 727 Input Parameter: 728 + m1 - number of elements to allocate in 1st chunk (may be zero) 729 . m2 - number of elements to allocate in 2nd chunk (may be zero) 730 . m3 - number of elements to allocate in 3rd chunk (may be zero) 731 - m4 - number of elements to allocate in 4th chunk (may be zero) 732 733 Output Parameter: 734 + r1 - memory allocated in first chunk 735 . r2 - memory allocated in second chunk 736 . r3 - memory allocated in third chunk 737 - r4 - memory allocated in fourth chunk 738 739 Level: developer 740 741 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc4(), PetscFree4() 742 743 Concepts: memory allocation 744 745 M*/ 746 #if defined(PETSC_USE_DEBUG) 747 #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)) 748 #else 749 #define PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4) \ 750 ((*(r2) = 0, *(r3) = 0, *(r4) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(m3)*sizeof(**(r3))+(m4)*sizeof(**(r4))+3*(PETSC_MEMALIGN-1),r1)) \ 751 || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),0) \ 752 || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0) || (!(m4) ? (*(r4) = 0,0) : 0)) 753 #endif 754 755 /*MC 756 PetscCalloc4 - Allocates 4 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN 757 758 Synopsis: 759 #include "petscsys.h" 760 PetscErrorCode PetscCalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4) 761 762 Not Collective 763 764 Input Parameter: 765 + m1 - number of elements to allocate in 1st chunk (may be zero) 766 . m2 - number of elements to allocate in 2nd chunk (may be zero) 767 . m3 - number of elements to allocate in 3rd chunk (may be zero) 768 - m4 - number of elements to allocate in 4th chunk (may be zero) 769 770 Output Parameter: 771 + r1 - memory allocated in first chunk 772 . r2 - memory allocated in second chunk 773 . r3 - memory allocated in third chunk 774 - r4 - memory allocated in fourth chunk 775 776 Level: developer 777 778 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc4(), PetscFree4() 779 780 Concepts: memory allocation 781 782 M*/ 783 #define PetscCalloc4(m1,r1,m2,r2,m3,r3,m4,r4) \ 784 (PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4) \ 785 || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \ 786 || PetscMemzero(*(r4),(m4)*sizeof(**(r4)))) 787 788 /*MC 789 PetscMalloc5 - Allocates 5 arrays of memory, all aligned to PETSC_MEMALIGN 790 791 Synopsis: 792 #include "petscsys.h" 793 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) 794 795 Not Collective 796 797 Input Parameter: 798 + m1 - number of elements to allocate in 1st chunk (may be zero) 799 . m2 - number of elements to allocate in 2nd chunk (may be zero) 800 . m3 - number of elements to allocate in 3rd chunk (may be zero) 801 . m4 - number of elements to allocate in 4th chunk (may be zero) 802 - m5 - number of elements to allocate in 5th chunk (may be zero) 803 804 Output Parameter: 805 + r1 - memory allocated in first chunk 806 . r2 - memory allocated in second chunk 807 . r3 - memory allocated in third chunk 808 . r4 - memory allocated in fourth chunk 809 - r5 - memory allocated in fifth chunk 810 811 Level: developer 812 813 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc5(), PetscFree5() 814 815 Concepts: memory allocation 816 817 M*/ 818 #if defined(PETSC_USE_DEBUG) 819 #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)) 820 #else 821 #define PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5) \ 822 ((*(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)) \ 823 || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),*(void**)(r5) = PetscAddrAlign(*(r4)+(m4)),0) \ 824 || (!(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)) 825 #endif 826 827 /*MC 828 PetscCalloc5 - Allocates 5 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN 829 830 Synopsis: 831 #include "petscsys.h" 832 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) 833 834 Not Collective 835 836 Input Parameter: 837 + m1 - number of elements to allocate in 1st chunk (may be zero) 838 . m2 - number of elements to allocate in 2nd chunk (may be zero) 839 . m3 - number of elements to allocate in 3rd chunk (may be zero) 840 . m4 - number of elements to allocate in 4th chunk (may be zero) 841 - m5 - number of elements to allocate in 5th chunk (may be zero) 842 843 Output Parameter: 844 + r1 - memory allocated in first chunk 845 . r2 - memory allocated in second chunk 846 . r3 - memory allocated in third chunk 847 . r4 - memory allocated in fourth chunk 848 - r5 - memory allocated in fifth chunk 849 850 Level: developer 851 852 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc5(), PetscFree5() 853 854 Concepts: memory allocation 855 856 M*/ 857 #define PetscCalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5) \ 858 (PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5) \ 859 || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \ 860 || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5)))) 861 862 /*MC 863 PetscMalloc6 - Allocates 6 arrays of memory, all aligned to PETSC_MEMALIGN 864 865 Synopsis: 866 #include "petscsys.h" 867 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) 868 869 Not Collective 870 871 Input Parameter: 872 + m1 - number of elements to allocate in 1st chunk (may be zero) 873 . m2 - number of elements to allocate in 2nd chunk (may be zero) 874 . m3 - number of elements to allocate in 3rd chunk (may be zero) 875 . m4 - number of elements to allocate in 4th chunk (may be zero) 876 . m5 - number of elements to allocate in 5th chunk (may be zero) 877 - m6 - number of elements to allocate in 6th chunk (may be zero) 878 879 Output Parameter: 880 + r1 - memory allocated in first chunk 881 . r2 - memory allocated in second chunk 882 . r3 - memory allocated in third chunk 883 . r4 - memory allocated in fourth chunk 884 . r5 - memory allocated in fifth chunk 885 - r6 - memory allocated in sixth chunk 886 887 Level: developer 888 889 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc6(), PetscFree3(), PetscFree4(), PetscFree5(), PetscFree6() 890 891 Concepts: memory allocation 892 893 M*/ 894 #if defined(PETSC_USE_DEBUG) 895 #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)) 896 #else 897 #define PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) \ 898 ((*(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)) \ 899 || (*(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) \ 900 || (!(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)) 901 #endif 902 903 /*MC 904 PetscCalloc6 - Allocates 6 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN 905 906 Synopsis: 907 #include "petscsys.h" 908 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) 909 910 Not Collective 911 912 Input Parameter: 913 + m1 - number of elements to allocate in 1st chunk (may be zero) 914 . m2 - number of elements to allocate in 2nd chunk (may be zero) 915 . m3 - number of elements to allocate in 3rd chunk (may be zero) 916 . m4 - number of elements to allocate in 4th chunk (may be zero) 917 . m5 - number of elements to allocate in 5th chunk (may be zero) 918 - m6 - number of elements to allocate in 6th chunk (may be zero) 919 920 Output Parameter: 921 + r1 - memory allocated in first chunk 922 . r2 - memory allocated in second chunk 923 . r3 - memory allocated in third chunk 924 . r4 - memory allocated in fourth chunk 925 . r5 - memory allocated in fifth chunk 926 - r6 - memory allocated in sixth chunk 927 928 Level: developer 929 930 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscMalloc6(), PetscFree6() 931 932 Concepts: memory allocation 933 M*/ 934 #define PetscCalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) \ 935 (PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) \ 936 || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \ 937 || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))) || PetscMemzero(*(r6),(m6)*sizeof(**(r6)))) 938 939 /*MC 940 PetscMalloc7 - Allocates 7 arrays of memory, all aligned to PETSC_MEMALIGN 941 942 Synopsis: 943 #include "petscsys.h" 944 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) 945 946 Not Collective 947 948 Input Parameter: 949 + m1 - number of elements to allocate in 1st chunk (may be zero) 950 . m2 - number of elements to allocate in 2nd chunk (may be zero) 951 . m3 - number of elements to allocate in 3rd chunk (may be zero) 952 . m4 - number of elements to allocate in 4th chunk (may be zero) 953 . m5 - number of elements to allocate in 5th chunk (may be zero) 954 . m6 - number of elements to allocate in 6th chunk (may be zero) 955 - m7 - number of elements to allocate in 7th chunk (may be zero) 956 957 Output Parameter: 958 + r1 - memory allocated in first chunk 959 . r2 - memory allocated in second chunk 960 . r3 - memory allocated in third chunk 961 . r4 - memory allocated in fourth chunk 962 . r5 - memory allocated in fifth chunk 963 . r6 - memory allocated in sixth chunk 964 - r7 - memory allocated in seventh chunk 965 966 Level: developer 967 968 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc7(), PetscFree7() 969 970 Concepts: memory allocation 971 972 M*/ 973 #if defined(PETSC_USE_DEBUG) 974 #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)) 975 #else 976 #define PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) \ 977 ((*(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)) \ 978 || (*(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) \ 979 || (!(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)) 980 #endif 981 982 /*MC 983 PetscCalloc7 - Allocates 7 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN 984 985 Synopsis: 986 #include "petscsys.h" 987 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) 988 989 Not Collective 990 991 Input Parameter: 992 + m1 - number of elements to allocate in 1st chunk (may be zero) 993 . m2 - number of elements to allocate in 2nd chunk (may be zero) 994 . m3 - number of elements to allocate in 3rd chunk (may be zero) 995 . m4 - number of elements to allocate in 4th chunk (may be zero) 996 . m5 - number of elements to allocate in 5th chunk (may be zero) 997 . m6 - number of elements to allocate in 6th chunk (may be zero) 998 - m7 - number of elements to allocate in 7th chunk (may be zero) 999 1000 Output Parameter: 1001 + r1 - memory allocated in first chunk 1002 . r2 - memory allocated in second chunk 1003 . r3 - memory allocated in third chunk 1004 . r4 - memory allocated in fourth chunk 1005 . r5 - memory allocated in fifth chunk 1006 . r6 - memory allocated in sixth chunk 1007 - r7 - memory allocated in seventh chunk 1008 1009 Level: developer 1010 1011 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscMalloc7(), PetscFree7() 1012 1013 Concepts: memory allocation 1014 M*/ 1015 #define PetscCalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) \ 1016 (PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) \ 1017 || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \ 1018 || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))) || PetscMemzero(*(r6),(m6)*sizeof(**(r6))) \ 1019 || PetscMemzero(*(r7),(m7)*sizeof(**(r7)))) 1020 1021 /*MC 1022 PetscNew - Allocates memory of a particular type, zeros the memory! Aligned to PETSC_MEMALIGN 1023 1024 Synopsis: 1025 #include "petscsys.h" 1026 PetscErrorCode PetscNew(type **result) 1027 1028 Not Collective 1029 1030 Output Parameter: 1031 . result - memory allocated, sized to match pointer type 1032 1033 Level: beginner 1034 1035 .seealso: PetscFree(), PetscMalloc(), PetscNewLog() 1036 1037 Concepts: memory allocation 1038 1039 M*/ 1040 #define PetscNew(b) PetscCalloc1(1,(b)) 1041 1042 /*MC 1043 PetscNewLog - Allocates memory of a type matching pointer, zeros the memory! Aligned to PETSC_MEMALIGN. Associates the memory allocated 1044 with the given object using PetscLogObjectMemory(). 1045 1046 Synopsis: 1047 #include "petscsys.h" 1048 PetscErrorCode PetscNewLog(PetscObject obj,type **result) 1049 1050 Not Collective 1051 1052 Input Parameter: 1053 . obj - object memory is logged to 1054 1055 Output Parameter: 1056 . result - memory allocated, sized to match pointer type 1057 1058 Level: developer 1059 1060 .seealso: PetscFree(), PetscMalloc(), PetscNew(), PetscLogObjectMemory() 1061 1062 Concepts: memory allocation 1063 1064 M*/ 1065 #define PetscNewLog(o,b) (PetscNew((b)) || ((o) ? PetscLogObjectMemory((PetscObject)o,sizeof(**(b))) : 0)) 1066 1067 /*MC 1068 PetscFree - Frees memory 1069 1070 Synopsis: 1071 #include "petscsys.h" 1072 PetscErrorCode PetscFree(void *memory) 1073 1074 Not Collective 1075 1076 Input Parameter: 1077 . memory - memory to free (the pointer is ALWAYS set to 0 upon sucess) 1078 1079 Level: beginner 1080 1081 Notes: Memory must have been obtained with PetscNew() or PetscMalloc() 1082 1083 .seealso: PetscNew(), PetscMalloc(), PetscFreeVoid() 1084 1085 Concepts: memory allocation 1086 1087 M*/ 1088 #define PetscFree(a) ((a) && ((*PetscTrFree)((void*)(a),__LINE__,PETSC_FUNCTION_NAME,__FILE__) || ((a) = 0,0))) 1089 1090 /*MC 1091 PetscFreeVoid - Frees memory 1092 1093 Synopsis: 1094 #include "petscsys.h" 1095 void PetscFreeVoid(void *memory) 1096 1097 Not Collective 1098 1099 Input Parameter: 1100 . memory - memory to free 1101 1102 Level: beginner 1103 1104 Notes: This is different from PetscFree() in that no error code is returned 1105 1106 .seealso: PetscFree(), PetscNew(), PetscMalloc() 1107 1108 Concepts: memory allocation 1109 1110 M*/ 1111 #define PetscFreeVoid(a) ((*PetscTrFree)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__),(a) = 0) 1112 1113 1114 /*MC 1115 PetscFree2 - Frees 2 chunks of memory obtained with PetscMalloc2() 1116 1117 Synopsis: 1118 #include "petscsys.h" 1119 PetscErrorCode PetscFree2(void *memory1,void *memory2) 1120 1121 Not Collective 1122 1123 Input Parameter: 1124 + memory1 - memory to free 1125 - memory2 - 2nd memory to free 1126 1127 Level: developer 1128 1129 Notes: Memory must have been obtained with PetscMalloc2() 1130 1131 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree() 1132 1133 Concepts: memory allocation 1134 1135 M*/ 1136 #if defined(PETSC_USE_DEBUG) 1137 #define PetscFree2(m1,m2) (PetscFree(m2) || PetscFree(m1)) 1138 #else 1139 #define PetscFree2(m1,m2) ((m1) ? ((m2)=0,PetscFree(m1)) : ((m1)=0,PetscFree(m2))) 1140 #endif 1141 1142 /*MC 1143 PetscFree3 - Frees 3 chunks of memory obtained with PetscMalloc3() 1144 1145 Synopsis: 1146 #include "petscsys.h" 1147 PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3) 1148 1149 Not Collective 1150 1151 Input Parameter: 1152 + memory1 - memory to free 1153 . memory2 - 2nd memory to free 1154 - memory3 - 3rd memory to free 1155 1156 Level: developer 1157 1158 Notes: Memory must have been obtained with PetscMalloc3() 1159 1160 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3() 1161 1162 Concepts: memory allocation 1163 1164 M*/ 1165 #if defined(PETSC_USE_DEBUG) 1166 #define PetscFree3(m1,m2,m3) (PetscFree(m3) || PetscFree(m2) || PetscFree(m1)) 1167 #else 1168 #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)))) 1169 #endif 1170 1171 /*MC 1172 PetscFree4 - Frees 4 chunks of memory obtained with PetscMalloc4() 1173 1174 Synopsis: 1175 #include "petscsys.h" 1176 PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4) 1177 1178 Not Collective 1179 1180 Input Parameter: 1181 + m1 - memory to free 1182 . m2 - 2nd memory to free 1183 . m3 - 3rd memory to free 1184 - m4 - 4th memory to free 1185 1186 Level: developer 1187 1188 Notes: Memory must have been obtained with PetscMalloc4() 1189 1190 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4() 1191 1192 Concepts: memory allocation 1193 1194 M*/ 1195 #if defined(PETSC_USE_DEBUG) 1196 #define PetscFree4(m1,m2,m3,m4) (PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1)) 1197 #else 1198 #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))))) 1199 #endif 1200 1201 /*MC 1202 PetscFree5 - Frees 5 chunks of memory obtained with PetscMalloc5() 1203 1204 Synopsis: 1205 #include "petscsys.h" 1206 PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5) 1207 1208 Not Collective 1209 1210 Input Parameter: 1211 + m1 - memory to free 1212 . m2 - 2nd memory to free 1213 . m3 - 3rd memory to free 1214 . m4 - 4th memory to free 1215 - m5 - 5th memory to free 1216 1217 Level: developer 1218 1219 Notes: Memory must have been obtained with PetscMalloc5() 1220 1221 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5() 1222 1223 Concepts: memory allocation 1224 1225 M*/ 1226 #if defined(PETSC_USE_DEBUG) 1227 #define PetscFree5(m1,m2,m3,m4,m5) (PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1)) 1228 #else 1229 #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)) : \ 1230 ((m4) ? ((m5)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m4)) : ((m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m5)))))) 1231 #endif 1232 1233 1234 /*MC 1235 PetscFree6 - Frees 6 chunks of memory obtained with PetscMalloc6() 1236 1237 Synopsis: 1238 #include "petscsys.h" 1239 PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6) 1240 1241 Not Collective 1242 1243 Input Parameter: 1244 + m1 - memory to free 1245 . m2 - 2nd memory to free 1246 . m3 - 3rd memory to free 1247 . m4 - 4th memory to free 1248 . m5 - 5th memory to free 1249 - m6 - 6th memory to free 1250 1251 1252 Level: developer 1253 1254 Notes: Memory must have been obtained with PetscMalloc6() 1255 1256 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6() 1257 1258 Concepts: memory allocation 1259 1260 M*/ 1261 #if defined(PETSC_USE_DEBUG) 1262 #define PetscFree6(m1,m2,m3,m4,m5,m6) (PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1)) 1263 #else 1264 #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)) : \ 1265 ((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)) : \ 1266 ((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))))))) 1267 #endif 1268 1269 /*MC 1270 PetscFree7 - Frees 7 chunks of memory obtained with PetscMalloc7() 1271 1272 Synopsis: 1273 #include "petscsys.h" 1274 PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7) 1275 1276 Not Collective 1277 1278 Input Parameter: 1279 + m1 - memory to free 1280 . m2 - 2nd memory to free 1281 . m3 - 3rd memory to free 1282 . m4 - 4th memory to free 1283 . m5 - 5th memory to free 1284 . m6 - 6th memory to free 1285 - m7 - 7th memory to free 1286 1287 1288 Level: developer 1289 1290 Notes: Memory must have been obtained with PetscMalloc7() 1291 1292 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(), 1293 PetscMalloc7() 1294 1295 Concepts: memory allocation 1296 1297 M*/ 1298 #if defined(PETSC_USE_DEBUG) 1299 #define PetscFree7(m1,m2,m3,m4,m5,m6,m7) (PetscFree(m7) || PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1)) 1300 #else 1301 #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)) : \ 1302 ((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)) : \ 1303 ((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)) : \ 1304 ((m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m7)))))))) 1305 #endif 1306 1307 PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,int,const char[],const char[],void**); 1308 PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]); 1309 PETSC_EXTERN PetscErrorCode PetscMallocSet(PetscErrorCode (*)(size_t,int,const char[],const char[],void**),PetscErrorCode (*)(void*,int,const char[],const char[])); 1310 PETSC_EXTERN PetscErrorCode PetscMallocClear(void); 1311 1312 /* 1313 PetscLogDouble variables are used to contain double precision numbers 1314 that are not used in the numerical computations, but rather in logging, 1315 timing etc. 1316 */ 1317 typedef double PetscLogDouble; 1318 #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE 1319 1320 /* 1321 Routines for tracing memory corruption/bleeding with default PETSc memory allocation 1322 */ 1323 PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *); 1324 PETSC_EXTERN PetscErrorCode PetscMallocDumpLog(FILE *); 1325 PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *); 1326 PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *); 1327 PETSC_EXTERN PetscErrorCode PetscMallocDebug(PetscBool); 1328 PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*); 1329 PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]); 1330 PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLog(void); 1331 PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLogThreshold(PetscLogDouble); 1332 PETSC_EXTERN PetscErrorCode PetscMallocGetDumpLog(PetscBool*); 1333 1334 /*E 1335 PetscDataType - Used for handling different basic data types. 1336 1337 Level: beginner 1338 1339 Developer comment: It would be nice if we could always just use MPI Datatypes, why can we not? 1340 1341 .seealso: PetscBinaryRead(), PetscBinaryWrite(), PetscDataTypeToMPIDataType(), 1342 PetscDataTypeGetSize() 1343 1344 E*/ 1345 typedef enum {PETSC_INT = 0,PETSC_DOUBLE = 1,PETSC_COMPLEX = 2, PETSC_LONG = 3 ,PETSC_SHORT = 4,PETSC_FLOAT = 5, 1346 PETSC_CHAR = 6,PETSC_BIT_LOGICAL = 7,PETSC_ENUM = 8,PETSC_BOOL=9, PETSC___FLOAT128 = 10,PETSC_OBJECT = 11, PETSC_FUNCTION = 12} PetscDataType; 1347 PETSC_EXTERN const char *const PetscDataTypes[]; 1348 1349 #if defined(PETSC_USE_COMPLEX) 1350 #define PETSC_SCALAR PETSC_COMPLEX 1351 #else 1352 #if defined(PETSC_USE_REAL_SINGLE) 1353 #define PETSC_SCALAR PETSC_FLOAT 1354 #elif defined(PETSC_USE_REAL___FLOAT128) 1355 #define PETSC_SCALAR PETSC___FLOAT128 1356 #else 1357 #define PETSC_SCALAR PETSC_DOUBLE 1358 #endif 1359 #endif 1360 #if defined(PETSC_USE_REAL_SINGLE) 1361 #define PETSC_REAL PETSC_FLOAT 1362 #elif defined(PETSC_USE_REAL___FLOAT128) 1363 #define PETSC_REAL PETSC___FLOAT128 1364 #else 1365 #define PETSC_REAL PETSC_DOUBLE 1366 #endif 1367 #define PETSC_FORTRANADDR PETSC_LONG 1368 1369 PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*); 1370 PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*); 1371 PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*); 1372 PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*); 1373 1374 /* 1375 Basic memory and string operations. These are usually simple wrappers 1376 around the basic Unix system calls, but a few of them have additional 1377 functionality and/or error checking. 1378 */ 1379 PETSC_EXTERN PetscErrorCode PetscBitMemcpy(void*,PetscInt,const void*,PetscInt,PetscInt,PetscDataType); 1380 PETSC_EXTERN PetscErrorCode PetscMemmove(void*,void *,size_t); 1381 PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool *); 1382 PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*); 1383 PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***); 1384 PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **); 1385 PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool *); 1386 PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool *); 1387 PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *); 1388 PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *); 1389 PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]); 1390 PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]); 1391 PETSC_EXTERN PetscErrorCode PetscStrncat(char[],const char[],size_t); 1392 PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t); 1393 PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]); 1394 PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]); 1395 PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]); 1396 PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]); 1397 PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]); 1398 PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]); 1399 PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*); 1400 PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*); 1401 PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*); 1402 PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]); 1403 PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***); 1404 PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***); 1405 PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t); 1406 1407 PETSC_EXTERN void PetscStrcmpNoError(const char[],const char[],PetscBool *); 1408 1409 /*S 1410 PetscToken - 'Token' used for managing tokenizing strings 1411 1412 Level: intermediate 1413 1414 .seealso: PetscTokenCreate(), PetscTokenFind(), PetscTokenDestroy() 1415 S*/ 1416 typedef struct _p_PetscToken* PetscToken; 1417 1418 PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*); 1419 PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]); 1420 PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*); 1421 1422 PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*); 1423 PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*); 1424 1425 /* 1426 These are MPI operations for MPI_Allreduce() etc 1427 */ 1428 PETSC_EXTERN MPI_Op PetscMaxSum_Op; 1429 #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) 1430 PETSC_EXTERN MPI_Op MPIU_SUM; 1431 #else 1432 #define MPIU_SUM MPI_SUM 1433 #endif 1434 #if defined(PETSC_USE_REAL___FLOAT128) 1435 PETSC_EXTERN MPI_Op MPIU_MAX; 1436 PETSC_EXTERN MPI_Op MPIU_MIN; 1437 #else 1438 #define MPIU_MAX MPI_MAX 1439 #define MPIU_MIN MPI_MIN 1440 #endif 1441 PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*); 1442 1443 PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm); 1444 PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm); 1445 1446 /*S 1447 PetscObject - any PETSc object, PetscViewer, Mat, Vec, KSP etc 1448 1449 Level: beginner 1450 1451 Note: This is the base class from which all PETSc objects are derived from. 1452 1453 .seealso: PetscObjectDestroy(), PetscObjectView(), PetscObjectGetName(), PetscObjectSetName(), PetscObjectReference(), PetscObjectDereferenc() 1454 S*/ 1455 typedef struct _p_PetscObject* PetscObject; 1456 1457 /*MC 1458 PetscObjectId - unique integer Id for a PetscObject 1459 1460 Level: developer 1461 1462 Notes: Unlike pointer values, object ids are never reused. 1463 1464 .seealso: PetscObjectState, PetscObjectGetId() 1465 M*/ 1466 typedef Petsc64bitInt PetscObjectId; 1467 1468 /*MC 1469 PetscObjectState - integer state for a PetscObject 1470 1471 Level: developer 1472 1473 Notes: 1474 Object state is always-increasing and (for objects that track state) can be used to determine if an object has 1475 changed since the last time you interacted with it. It is 64-bit so that it will not overflow for a very long time. 1476 1477 .seealso: PetscObjectId, PetscObjectStateGet(), PetscObjectStateIncrease(), PetscObjectStateSet() 1478 M*/ 1479 typedef Petsc64bitInt PetscObjectState; 1480 1481 /*S 1482 PetscFunctionList - Linked list of functions, possibly stored in dynamic libraries, accessed 1483 by string name 1484 1485 Level: advanced 1486 1487 .seealso: PetscFunctionListAdd(), PetscFunctionListDestroy(), PetscOpFlist 1488 S*/ 1489 typedef struct _n_PetscFunctionList *PetscFunctionList; 1490 1491 /*E 1492 PetscFileMode - Access mode for a file. 1493 1494 Level: beginner 1495 1496 FILE_MODE_READ - open a file at its beginning for reading 1497 1498 FILE_MODE_WRITE - open a file at its beginning for writing (will create if the file does not exist) 1499 1500 FILE_MODE_APPEND - open a file at end for writing 1501 1502 FILE_MODE_UPDATE - open a file for updating, meaning for reading and writing 1503 1504 FILE_MODE_APPEND_UPDATE - open a file for updating, meaning for reading and writing, at the end 1505 1506 .seealso: PetscViewerFileSetMode() 1507 E*/ 1508 typedef enum {FILE_MODE_READ, FILE_MODE_WRITE, FILE_MODE_APPEND, FILE_MODE_UPDATE, FILE_MODE_APPEND_UPDATE} PetscFileMode; 1509 extern const char *const PetscFileModes[]; 1510 1511 /* 1512 Defines PETSc error handling. 1513 */ 1514 #include <petscerror.h> 1515 1516 #define PETSC_SMALLEST_CLASSID 1211211 1517 PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID; 1518 PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID; 1519 PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *); 1520 1521 /* 1522 Routines that get memory usage information from the OS 1523 */ 1524 PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *); 1525 PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *); 1526 PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void); 1527 1528 PETSC_EXTERN PetscErrorCode PetscInfoAllow(PetscBool ,const char []); 1529 PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal); 1530 1531 /* 1532 Initialization of PETSc 1533 */ 1534 PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]); 1535 PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]); 1536 PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void); 1537 PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *); 1538 PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *); 1539 PETSC_EXTERN PetscErrorCode PetscFinalize(void); 1540 PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void); 1541 PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***); 1542 PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***); 1543 PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **); 1544 1545 PETSC_EXTERN PetscErrorCode PetscEnd(void); 1546 PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void); 1547 1548 PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]); 1549 PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void); 1550 PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void); 1551 PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]); 1552 1553 /* 1554 These are so that in extern C code we can caste function pointers to non-extern C 1555 function pointers. Since the regular C++ code expects its function pointers to be C++ 1556 */ 1557 PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void); 1558 PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void); 1559 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void); 1560 1561 /* 1562 Functions that can act on any PETSc object. 1563 */ 1564 PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*); 1565 PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *); 1566 PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *); 1567 PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]); 1568 PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []); 1569 PETSC_EXTERN PetscErrorCode PetscObjectSetPrecision(PetscObject,PetscPrecision); 1570 PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]); 1571 PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]); 1572 PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]); 1573 PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt); 1574 PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*); 1575 PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt); 1576 PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject); 1577 PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*); 1578 PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject); 1579 PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt *); 1580 PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject); 1581 PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]); 1582 PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject *); 1583 PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void)); 1584 #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d)) 1585 PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject); 1586 PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject); 1587 PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt *); 1588 PETSC_EXTERN PetscErrorCode PetscObjectAddOptionsHandler(PetscObject,PetscErrorCode (*)(PetscObject,void*),PetscErrorCode (*)(PetscObject,void*),void*); 1589 PETSC_EXTERN PetscErrorCode PetscObjectProcessOptionsHandlers(PetscObject); 1590 PETSC_EXTERN PetscErrorCode PetscObjectDestroyOptionsHandlers(PetscObject); 1591 PETSC_EXTERN PetscErrorCode PetscObjectsGetGlobalNumbering(MPI_Comm,PetscInt,PetscObject*,PetscInt*,PetscInt*); 1592 1593 #include <petscviewertypes.h> 1594 #include <petscoptions.h> 1595 1596 PETSC_EXTERN PetscErrorCode PetscMemoryShowUsage(PetscViewer,const char[]); 1597 PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer); 1598 PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer); 1599 #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr)) 1600 PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void)); 1601 PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]); 1602 PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]); 1603 PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]); 1604 PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]); 1605 PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]); 1606 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject); 1607 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void); 1608 PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,const char[],const char[]); 1609 PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject); 1610 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *); 1611 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...); 1612 PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void)); 1613 PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void); 1614 1615 #if defined(PETSC_HAVE_SAWS) 1616 PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void); 1617 PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject); 1618 PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool); 1619 PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject); 1620 PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject); 1621 PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject); 1622 PETSC_EXTERN void PetscStackSAWsGrantAccess(void); 1623 PETSC_EXTERN void PetscStackSAWsTakeAccess(void); 1624 PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void); 1625 PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void); 1626 1627 #else 1628 #define PetscSAWsBlock() 0 1629 #define PetscObjectSAWsViewOff(obj) 0 1630 #define PetscObjectSAWsSetBlock(obj,flg) 0 1631 #define PetscObjectSAWsBlock(obj) 0 1632 #define PetscObjectSAWsGrantAccess(obj) 0 1633 #define PetscObjectSAWsTakeAccess(obj) 0 1634 #define PetscStackViewSAWs() 0 1635 #define PetscStackSAWsViewOff() 0 1636 #define PetscStackSAWsTakeAccess() 1637 #define PetscStackSAWsGrantAccess() 1638 1639 #endif 1640 1641 typedef void* PetscDLHandle; 1642 typedef enum {PETSC_DL_DECIDE=0,PETSC_DL_NOW=1,PETSC_DL_LOCAL=2} PetscDLMode; 1643 PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *); 1644 PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *); 1645 PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **); 1646 1647 1648 #if defined(PETSC_USE_DEBUG) 1649 PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**); 1650 #endif 1651 PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool); 1652 1653 /*S 1654 PetscObjectList - Linked list of PETSc objects, each accessable by string name 1655 1656 Level: developer 1657 1658 Notes: Used by PetscObjectCompose() and PetscObjectQuery() 1659 1660 .seealso: PetscObjectListAdd(), PetscObjectListDestroy(), PetscObjectListFind(), PetscObjectCompose(), PetscObjectQuery(), PetscFunctionList 1661 S*/ 1662 typedef struct _n_PetscObjectList *PetscObjectList; 1663 1664 PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*); 1665 PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*); 1666 PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*); 1667 PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject); 1668 PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]); 1669 PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *); 1670 1671 /* 1672 Dynamic library lists. Lists of names of routines in objects or in dynamic 1673 link libraries that will be loaded as needed. 1674 */ 1675 1676 #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr)) 1677 PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void)); 1678 PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*); 1679 #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr)) 1680 PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void)); 1681 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[]); 1682 PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *); 1683 PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer); 1684 PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*); 1685 1686 /*S 1687 PetscDLLibrary - Linked list of dynamics libraries to search for functions 1688 1689 Level: advanced 1690 1691 .seealso: PetscDLLibraryOpen() 1692 S*/ 1693 typedef struct _n_PetscDLLibrary *PetscDLLibrary; 1694 PETSC_EXTERN PetscDLLibrary PetscDLLibrariesLoaded; 1695 PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]); 1696 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]); 1697 PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **); 1698 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary); 1699 PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool *); 1700 PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *); 1701 PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary); 1702 1703 /* 1704 Useful utility routines 1705 */ 1706 PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*); 1707 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*); 1708 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt); 1709 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt); 1710 PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject); 1711 PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*); 1712 1713 /* 1714 PetscNot - negates a logical type value and returns result as a PetscBool 1715 1716 Notes: This is useful in cases like 1717 $ int *a; 1718 $ PetscBool flag = PetscNot(a) 1719 where !a does not return a PetscBool because we cannot provide a cast from int to PetscBool in C. 1720 */ 1721 #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE) 1722 1723 #if defined(PETSC_HAVE_VALGRIND) 1724 # include <valgrind/valgrind.h> 1725 # define PETSC_RUNNING_ON_VALGRIND RUNNING_ON_VALGRIND 1726 #else 1727 # define PETSC_RUNNING_ON_VALGRIND PETSC_FALSE 1728 #endif 1729 1730 1731 /*MC 1732 PetscHelpPrintf - Prints help messages. 1733 1734 Synopsis: 1735 #include "petscsys.h" 1736 PetscErrorCode (*PetscHelpPrintf)(const char format[],...); 1737 1738 Not Collective 1739 1740 Input Parameters: 1741 . format - the usual printf() format string 1742 1743 Level: developer 1744 1745 Fortran Note: 1746 This routine is not supported in Fortran. 1747 1748 Concepts: help messages^printing 1749 Concepts: printing^help messages 1750 1751 .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf() 1752 M*/ 1753 PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...); 1754 1755 /* 1756 Defines PETSc profiling. 1757 */ 1758 #include <petsclog.h> 1759 1760 1761 1762 /* 1763 Simple PETSc parallel IO for ASCII printing 1764 */ 1765 PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]); 1766 PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**); 1767 PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*); 1768 PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...); 1769 PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...); 1770 PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...); 1771 PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...); 1772 1773 /* These are used internally by PETSc ASCII IO routines*/ 1774 #include <stdarg.h> 1775 PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list); 1776 PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE*,const char[],va_list); 1777 PETSC_EXTERN PetscErrorCode PetscVFPrintfDefault(FILE*,const char[],va_list); 1778 1779 #if defined(PETSC_HAVE_MATLAB_ENGINE) 1780 PETSC_EXTERN PetscErrorCode PetscVFPrintf_Matlab(FILE*,const char[],va_list); 1781 #endif 1782 1783 PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...); 1784 PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...); 1785 PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...); 1786 1787 #if defined(PETSC_HAVE_POPEN) 1788 PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **); 1789 PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*,int*); 1790 PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]); 1791 #endif 1792 1793 PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...); 1794 PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...); 1795 PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*); 1796 PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]); 1797 PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**); 1798 PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**); 1799 PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]); 1800 1801 PETSC_EXTERN PetscErrorCode PetscPopUpSelect(MPI_Comm,const char*,const char*,int,const char**,int*); 1802 1803 /*S 1804 PetscContainer - Simple PETSc object that contains a pointer to any required data 1805 1806 Level: advanced 1807 1808 .seealso: PetscObject, PetscContainerCreate() 1809 S*/ 1810 PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID; 1811 typedef struct _p_PetscContainer* PetscContainer; 1812 PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void **); 1813 PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void *); 1814 PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*); 1815 PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer *); 1816 PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void*)); 1817 1818 /* 1819 For use in debuggers 1820 */ 1821 PETSC_EXTERN PetscMPIInt PetscGlobalRank; 1822 PETSC_EXTERN PetscMPIInt PetscGlobalSize; 1823 PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer); 1824 PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer); 1825 PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer); 1826 1827 #include <stddef.h> 1828 #include <string.h> /* for memcpy, memset */ 1829 #if defined(PETSC_HAVE_STDLIB_H) 1830 #include <stdlib.h> 1831 #endif 1832 1833 #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__) 1834 #include <xmmintrin.h> 1835 #endif 1836 1837 #undef __FUNCT__ 1838 #define __FUNCT__ "PetscMemcpy" 1839 /*@C 1840 PetscMemcpy - Copies n bytes, beginning at location b, to the space 1841 beginning at location a. The two memory regions CANNOT overlap, use 1842 PetscMemmove() in that case. 1843 1844 Not Collective 1845 1846 Input Parameters: 1847 + b - pointer to initial memory space 1848 - n - length (in bytes) of space to copy 1849 1850 Output Parameter: 1851 . a - pointer to copy space 1852 1853 Level: intermediate 1854 1855 Compile Option: 1856 PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used 1857 for memory copies on double precision values. 1858 PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used 1859 for memory copies on double precision values. 1860 PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used 1861 for memory copies on double precision values. 1862 1863 Note: 1864 This routine is analogous to memcpy(). 1865 1866 Developer Note: this is inlined for fastest performance 1867 1868 Concepts: memory^copying 1869 Concepts: copying^memory 1870 1871 .seealso: PetscMemmove() 1872 1873 @*/ 1874 PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n) 1875 { 1876 #if defined(PETSC_USE_DEBUG) 1877 unsigned long al = (unsigned long) a,bl = (unsigned long) b; 1878 unsigned long nl = (unsigned long) n; 1879 PetscFunctionBegin; 1880 if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer"); 1881 if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer"); 1882 #else 1883 PetscFunctionBegin; 1884 #endif 1885 if (a != b) { 1886 #if defined(PETSC_USE_DEBUG) 1887 if ((al > bl && (al - bl) < nl) || (bl - al) < nl) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\ 1888 or make sure your copy regions and lengths are correct. \n\ 1889 Length (bytes) %ld first address %ld second address %ld",nl,al,bl); 1890 #endif 1891 #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY)) 1892 if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) { 1893 size_t len = n/sizeof(PetscScalar); 1894 #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) 1895 PetscBLASInt one = 1,blen; 1896 PetscErrorCode ierr; 1897 ierr = PetscBLASIntCast(len,&blen);CHKERRQ(ierr); 1898 PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one)); 1899 #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY) 1900 fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a); 1901 #else 1902 size_t i; 1903 PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a; 1904 for (i=0; i<len; i++) y[i] = x[i]; 1905 #endif 1906 } else { 1907 memcpy((char*)(a),(char*)(b),n); 1908 } 1909 #else 1910 memcpy((char*)(a),(char*)(b),n); 1911 #endif 1912 } 1913 PetscFunctionReturn(0); 1914 } 1915 1916 /*@C 1917 PetscMemzero - Zeros the specified memory. 1918 1919 Not Collective 1920 1921 Input Parameters: 1922 + a - pointer to beginning memory location 1923 - n - length (in bytes) of memory to initialize 1924 1925 Level: intermediate 1926 1927 Compile Option: 1928 PETSC_PREFER_BZERO - on certain machines (the IBM RS6000) the bzero() routine happens 1929 to be faster than the memset() routine. This flag causes the bzero() routine to be used. 1930 1931 Developer Note: this is inlined for fastest performance 1932 1933 Concepts: memory^zeroing 1934 Concepts: zeroing^memory 1935 1936 .seealso: PetscMemcpy() 1937 @*/ 1938 PETSC_STATIC_INLINE PetscErrorCode PetscMemzero(void *a,size_t n) 1939 { 1940 if (n > 0) { 1941 #if defined(PETSC_USE_DEBUG) 1942 if (!a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer"); 1943 #endif 1944 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) 1945 if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) { 1946 size_t i,len = n/sizeof(PetscScalar); 1947 PetscScalar *x = (PetscScalar*)a; 1948 for (i=0; i<len; i++) x[i] = 0.0; 1949 } else { 1950 #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO) 1951 if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) { 1952 PetscInt len = n/sizeof(PetscScalar); 1953 fortranzero_(&len,(PetscScalar*)a); 1954 } else { 1955 #endif 1956 #if defined(PETSC_PREFER_BZERO) 1957 bzero((char *)a,n); 1958 #else 1959 memset((char*)a,0,n); 1960 #endif 1961 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO) 1962 } 1963 #endif 1964 } 1965 return 0; 1966 } 1967 1968 /*MC 1969 PetscPrefetchBlock - Prefetches a block of memory 1970 1971 Synopsis: 1972 #include "petscsys.h" 1973 void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t) 1974 1975 Not Collective 1976 1977 Input Parameters: 1978 + a - pointer to first element to fetch (any type but usually PetscInt or PetscScalar) 1979 . n - number of elements to fetch 1980 . rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors) 1981 - t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note 1982 1983 Level: developer 1984 1985 Notes: 1986 The last two arguments (rw and t) must be compile-time constants. 1987 1988 Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality. Not all architectures offer 1989 equivalent locality hints, but the following macros are always defined to their closest analogue. 1990 + PETSC_PREFETCH_HINT_NTA - Non-temporal. Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched). 1991 . 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. 1992 . PETSC_PREFETCH_HINT_T1 - Fetch to level 2 and higher (not L1). 1993 - PETSC_PREFETCH_HINT_T2 - Fetch to high-level cache only. (On many systems, T0 and T1 are equivalent.) 1994 1995 This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid 1996 address). 1997 1998 Concepts: memory 1999 M*/ 2000 #define PetscPrefetchBlock(a,n,rw,t) do { \ 2001 const char *_p = (const char*)(a),*_end = (const char*)((a)+(n)); \ 2002 for ( ; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \ 2003 } while (0) 2004 2005 /* 2006 Determine if some of the kernel computation routines use 2007 Fortran (rather than C) for the numerical calculations. On some machines 2008 and compilers (like complex numbers) the Fortran version of the routines 2009 is faster than the C/C++ versions. The flag --with-fortran-kernels 2010 should be used with ./configure to turn these on. 2011 */ 2012 #if defined(PETSC_USE_FORTRAN_KERNELS) 2013 2014 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL) 2015 #define PETSC_USE_FORTRAN_KERNEL_MULTCRL 2016 #endif 2017 2018 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM) 2019 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM 2020 #endif 2021 2022 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ) 2023 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ 2024 #endif 2025 2026 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ) 2027 #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ 2028 #endif 2029 2030 #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM) 2031 #define PETSC_USE_FORTRAN_KERNEL_NORM 2032 #endif 2033 2034 #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY) 2035 #define PETSC_USE_FORTRAN_KERNEL_MAXPY 2036 #endif 2037 2038 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 2039 #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ 2040 #endif 2041 2042 #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ) 2043 #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ 2044 #endif 2045 2046 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ) 2047 #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ 2048 #endif 2049 2050 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ) 2051 #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ 2052 #endif 2053 2054 #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT) 2055 #define PETSC_USE_FORTRAN_KERNEL_MDOT 2056 #endif 2057 2058 #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY) 2059 #define PETSC_USE_FORTRAN_KERNEL_XTIMESY 2060 #endif 2061 2062 #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX) 2063 #define PETSC_USE_FORTRAN_KERNEL_AYPX 2064 #endif 2065 2066 #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY) 2067 #define PETSC_USE_FORTRAN_KERNEL_WAXPY 2068 #endif 2069 2070 #endif 2071 2072 /* 2073 Macros for indicating code that should be compiled with a C interface, 2074 rather than a C++ interface. Any routines that are dynamically loaded 2075 (such as the PCCreate_XXX() routines) must be wrapped so that the name 2076 mangler does not change the functions symbol name. This just hides the 2077 ugly extern "C" {} wrappers. 2078 */ 2079 #if defined(__cplusplus) 2080 #define EXTERN_C_BEGIN extern "C" { 2081 #define EXTERN_C_END } 2082 #else 2083 #define EXTERN_C_BEGIN 2084 #define EXTERN_C_END 2085 #endif 2086 2087 /* --------------------------------------------------------------------*/ 2088 2089 /*MC 2090 MPI_Comm - the basic object used by MPI to determine which processes are involved in a 2091 communication 2092 2093 Level: beginner 2094 2095 Note: This manual page is a place-holder because MPICH does not have a manual page for MPI_Comm 2096 2097 .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF 2098 M*/ 2099 2100 /*MC 2101 PetscScalar - PETSc type that represents either a double precision real number, a double precision 2102 complex number, a single precision real number, a long double or an int - if the code is configured 2103 with --with-scalar-type=real,complex --with-precision=single,double,longdouble,int,matsingle 2104 2105 2106 Level: beginner 2107 2108 .seealso: PetscReal, PassiveReal, PassiveScalar, MPIU_SCALAR, PetscInt 2109 M*/ 2110 2111 /*MC 2112 PetscComplex - PETSc type that represents a complex number with precision matching that of PetscReal. 2113 2114 Synopsis: 2115 #define PETSC_DESIRE_COMPLEX 2116 #include <petscsys.h> 2117 PetscComplex number = 1. + 2.*PETSC_i; 2118 2119 Level: beginner 2120 2121 Note: 2122 Complex numbers are automatically available if PETSc was configured --with-scalar-type=complex (in which case 2123 PetscComplex will match PetscScalar), otherwise the macro PETSC_DESIRE_COMPLEX must be defined before including any 2124 PETSc headers. If the compiler supports complex numbers, PetscComplex and associated variables and functions will be 2125 defined and PETSC_HAVE_COMPLEX will be set. 2126 2127 .seealso: PetscReal, PetscComplex, MPIU_COMPLEX, PetscInt, PETSC_i 2128 M*/ 2129 2130 /*MC 2131 PetscReal - PETSc type that represents a real number version of PetscScalar 2132 2133 Level: beginner 2134 2135 .seealso: PetscScalar, PassiveReal, PassiveScalar 2136 M*/ 2137 2138 /*MC 2139 PassiveScalar - PETSc type that represents a PetscScalar 2140 Level: beginner 2141 2142 This is the same as a PetscScalar except in code that is automatically differentiated it is 2143 treated as a constant (not an indendent or dependent variable) 2144 2145 .seealso: PetscReal, PassiveReal, PetscScalar 2146 M*/ 2147 2148 /*MC 2149 PassiveReal - PETSc type that represents a PetscReal 2150 2151 Level: beginner 2152 2153 This is the same as a PetscReal except in code that is automatically differentiated it is 2154 treated as a constant (not an indendent or dependent variable) 2155 2156 .seealso: PetscScalar, PetscReal, PassiveScalar 2157 M*/ 2158 2159 /*MC 2160 MPIU_SCALAR - MPI datatype corresponding to PetscScalar 2161 2162 Level: beginner 2163 2164 Note: In MPI calls that require an MPI datatype that matches a PetscScalar or array of PetscScalars 2165 pass this value 2166 2167 .seealso: PetscReal, PassiveReal, PassiveScalar, PetscScalar, MPIU_INT 2168 M*/ 2169 2170 #if defined(PETSC_HAVE_MPIIO) 2171 #if !defined(PETSC_WORDS_BIGENDIAN) 2172 PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*); 2173 PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*); 2174 #else 2175 #define MPIU_File_write_all(a,b,c,d,e) MPI_File_write_all(a,b,c,d,e) 2176 #define MPIU_File_read_all(a,b,c,d,e) MPI_File_read_all(a,b,c,d,e) 2177 #endif 2178 #endif 2179 2180 /* the following petsc_static_inline require petscerror.h */ 2181 2182 /* Limit MPI to 32-bits */ 2183 #define PETSC_MPI_INT_MAX 2147483647 2184 #define PETSC_MPI_INT_MIN -2147483647 2185 /* Limit BLAS to 32-bits */ 2186 #define PETSC_BLAS_INT_MAX 2147483647 2187 #define PETSC_BLAS_INT_MIN -2147483647 2188 2189 #undef __FUNCT__ 2190 #define __FUNCT__ "PetscBLASIntCast" 2191 PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b) 2192 { 2193 PetscFunctionBegin; 2194 #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES) 2195 if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK"); 2196 #endif 2197 *b = (PetscBLASInt)(a); 2198 PetscFunctionReturn(0); 2199 } 2200 2201 #undef __FUNCT__ 2202 #define __FUNCT__ "PetscMPIIntCast" 2203 PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b) 2204 { 2205 PetscFunctionBegin; 2206 #if defined(PETSC_USE_64BIT_INDICES) 2207 if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for MPI"); 2208 #endif 2209 *b = (PetscMPIInt)(a); 2210 PetscFunctionReturn(0); 2211 } 2212 2213 2214 /* 2215 The IBM include files define hz, here we hide it so that it may be used as a regular user variable. 2216 */ 2217 #if defined(hz) 2218 #undef hz 2219 #endif 2220 2221 /* For arrays that contain filenames or paths */ 2222 2223 2224 #if defined(PETSC_HAVE_LIMITS_H) 2225 #include <limits.h> 2226 #endif 2227 #if defined(PETSC_HAVE_SYS_PARAM_H) 2228 #include <sys/param.h> 2229 #endif 2230 #if defined(PETSC_HAVE_SYS_TYPES_H) 2231 #include <sys/types.h> 2232 #endif 2233 #if defined(MAXPATHLEN) 2234 # define PETSC_MAX_PATH_LEN MAXPATHLEN 2235 #elif defined(MAX_PATH) 2236 # define PETSC_MAX_PATH_LEN MAX_PATH 2237 #elif defined(_MAX_PATH) 2238 # define PETSC_MAX_PATH_LEN _MAX_PATH 2239 #else 2240 # define PETSC_MAX_PATH_LEN 4096 2241 #endif 2242 2243 /*MC 2244 2245 UsingFortran - Fortran can be used with PETSc in four distinct approaches 2246 2247 $ 1) classic Fortran 77 style 2248 $#include "finclude/petscXXX.h" to work with material from the XXX component of PETSc 2249 $ XXX variablename 2250 $ You cannot use this approach if you wish to use the Fortran 90 specific PETSc routines 2251 $ which end in F90; such as VecGetArrayF90() 2252 $ 2253 $ 2) classic Fortran 90 style 2254 $#include "finclude/petscXXX.h" 2255 $#include "finclude/petscXXX.h90" to work with material from the XXX component of PETSc 2256 $ XXX variablename 2257 $ 2258 $ 3) Using Fortran modules 2259 $#include "finclude/petscXXXdef.h" 2260 $ use petscXXXX 2261 $ XXX variablename 2262 $ 2263 $ 4) Use Fortran modules and Fortran data types for PETSc types 2264 $#include "finclude/petscXXXdef.h" 2265 $ use petscXXXX 2266 $ type(XXX) variablename 2267 $ To use this approach you must ./configure PETSc with the additional 2268 $ option --with-fortran-datatypes You cannot use the type(XXX) declaration approach without using Fortran modules 2269 2270 Finally if you absolutely do not want to use any #include you can use either 2271 2272 $ 3a) skip the #include BUT you cannot use any PETSc data type names like Vec, Mat, PetscInt, PetscErrorCode etc 2273 $ and you must declare the variables as integer, for example 2274 $ integer variablename 2275 $ 2276 $ 4a) skip the #include, you use the object types like type(Vec) type(Mat) but cannot use the data type 2277 $ names like PetscErrorCode, PetscInt etc. again for those you must use integer 2278 2279 We recommend either 2 or 3. Approaches 2 and 3 provide type checking for most PETSc function calls; 4 has type checking 2280 for only a few PETSc functions. 2281 2282 Fortran type checking with interfaces is strick, this means you cannot pass a scalar value when an array value 2283 is expected (even though it is legal Fortran). For example when setting a single value in a matrix with MatSetValues() 2284 you cannot have something like 2285 $ PetscInt row,col 2286 $ PetscScalar val 2287 $ ... 2288 $ call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr) 2289 You must instead have 2290 $ PetscInt row(1),col(1) 2291 $ PetscScalar val(1) 2292 $ ... 2293 $ call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr) 2294 2295 2296 See the example src/vec/vec/examples/tutorials/ex20f90.F90 for an example that can use all four approaches 2297 2298 Developer Notes: The finclude/petscXXXdef.h contain all the #defines (would be typedefs in C code) these 2299 automatically include their predecessors; for example finclude/petscvecdef.h includes finclude/petscisdef.h 2300 2301 The finclude/petscXXXX.h contain all the parameter statements for that package. These automatically include 2302 their finclude/petscXXXdef.h file but DO NOT automatically include their predecessors; for example 2303 finclude/petscvec.h does NOT automatically include finclude/petscis.h 2304 2305 The finclude/ftn-custom/petscXXXdef.h90 are not intended to be used directly in code, they define the 2306 Fortran data type type(XXX) (for example type(Vec)) when PETSc is ./configure with the --with-fortran-datatypes option. 2307 2308 The finclude/ftn-custom/petscXXX.h90 (not included directly by code) contain interface definitions for 2309 the PETSc Fortran stubs that have different bindings then their C version (for example VecGetArrayF90). 2310 2311 The finclude/ftn-auto/petscXXX.h90 (not included directly by code) contain interface definitions generated 2312 automatically by "make allfortranstubs". 2313 2314 The finclude/petscXXX.h90 includes the custom finclude/ftn-custom/petscXXX.h90 and if ./configure 2315 was run with --with-fortran-interfaces it also includes the finclude/ftn-auto/petscXXX.h90 These DO NOT automatically 2316 include their predecessors 2317 2318 Level: beginner 2319 2320 M*/ 2321 2322 PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t); 2323 PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t); 2324 PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t); 2325 PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t); 2326 PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]); 2327 PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t); 2328 PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t); 2329 2330 PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]); 2331 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]); 2332 PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*); 2333 PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]); 2334 PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]); 2335 PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]); 2336 PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt*,PetscInt*,PetscInt*); 2337 PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]); 2338 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]); 2339 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]); 2340 PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]); 2341 PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]); 2342 PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]); 2343 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]); 2344 PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]); 2345 PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]); 2346 PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**); 2347 PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt*,const PetscInt*,PetscInt,const PetscInt*,const PetscInt*,PetscInt*,PetscInt**,PetscInt**); 2348 2349 PETSC_EXTERN PetscErrorCode PetscSetDisplay(void); 2350 PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t); 2351 2352 /*J 2353 PetscRandomType - String with the name of a PETSc randomizer 2354 2355 Level: beginner 2356 2357 Notes: to use the SPRNG you must have ./configure PETSc 2358 with the option --download-sprng 2359 2360 .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate() 2361 J*/ 2362 typedef const char* PetscRandomType; 2363 #define PETSCRAND "rand" 2364 #define PETSCRAND48 "rand48" 2365 #define PETSCSPRNG "sprng" 2366 2367 /* Logging support */ 2368 PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID; 2369 2370 PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void); 2371 2372 /*S 2373 PetscRandom - Abstract PETSc object that manages generating random numbers 2374 2375 Level: intermediate 2376 2377 Concepts: random numbers 2378 2379 .seealso: PetscRandomCreate(), PetscRandomGetValue(), PetscRandomType 2380 S*/ 2381 typedef struct _p_PetscRandom* PetscRandom; 2382 2383 /* Dynamic creation and loading functions */ 2384 PETSC_EXTERN PetscFunctionList PetscRandomList; 2385 PETSC_EXTERN PetscBool PetscRandomRegisterAllCalled; 2386 2387 PETSC_EXTERN PetscErrorCode PetscRandomRegisterAll(void); 2388 PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom)); 2389 PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType); 2390 PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom); 2391 PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType*); 2392 PETSC_STATIC_INLINE PetscErrorCode PetscRandomViewFromOptions(PetscRandom A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 2393 PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer); 2394 2395 PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*); 2396 PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*); 2397 PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*); 2398 PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*); 2399 PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar); 2400 PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long); 2401 PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *); 2402 PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom); 2403 PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*); 2404 2405 PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t); 2406 PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t); 2407 PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t); 2408 PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]); 2409 PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t); 2410 PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *); 2411 PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *); 2412 2413 PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscDataType); 2414 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscDataType); 2415 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,void*,PetscInt,PetscDataType,PetscBool ); 2416 PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,void*,PetscInt,PetscDataType,PetscBool ); 2417 PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *); 2418 PETSC_EXTERN PetscErrorCode PetscBinaryClose(int); 2419 PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool *); 2420 PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool *); 2421 PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t); 2422 PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *); 2423 PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *); 2424 PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*); 2425 PETSC_EXTERN PetscErrorCode PetscWebServe(MPI_Comm,int); 2426 2427 /* 2428 In binary files variables are stored using the following lengths, 2429 regardless of how they are stored in memory on any one particular 2430 machine. Use these rather then sizeof() in computing sizes for 2431 PetscBinarySeek(). 2432 */ 2433 #define PETSC_BINARY_INT_SIZE (32/8) 2434 #define PETSC_BINARY_FLOAT_SIZE (32/8) 2435 #define PETSC_BINARY_CHAR_SIZE (8/8) 2436 #define PETSC_BINARY_SHORT_SIZE (16/8) 2437 #define PETSC_BINARY_DOUBLE_SIZE (64/8) 2438 #define PETSC_BINARY_SCALAR_SIZE sizeof(PetscScalar) 2439 2440 /*E 2441 PetscBinarySeekType - argument to PetscBinarySeek() 2442 2443 Level: advanced 2444 2445 .seealso: PetscBinarySeek(), PetscBinarySynchronizedSeek() 2446 E*/ 2447 typedef enum {PETSC_BINARY_SEEK_SET = 0,PETSC_BINARY_SEEK_CUR = 1,PETSC_BINARY_SEEK_END = 2} PetscBinarySeekType; 2448 PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*); 2449 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*); 2450 PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt); 2451 2452 PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]); 2453 PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool ); 2454 PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void); 2455 PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*); 2456 PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void); 2457 PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void); 2458 2459 PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*); 2460 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**); 2461 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**); 2462 PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**); 2463 PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**); 2464 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscInt,const PetscMPIInt*,const void*,PetscInt*,PetscMPIInt**,void*) PetscAttrMPIPointerWithType(6,3); 2465 2466 /*E 2467 PetscBuildTwoSidedType - algorithm for setting up two-sided communication 2468 2469 $ PETSC_BUILDTWOSIDED_ALLREDUCE - classical algorithm using an MPI_Allreduce with 2470 $ a buffer of length equal to the communicator size. Not memory-scalable due to 2471 $ the large reduction size. Requires only MPI-1. 2472 $ PETSC_BUILDTWOSIDED_IBARRIER - nonblocking algorithm based on MPI_Issend and MPI_Ibarrier. 2473 $ Proved communication-optimal in Hoefler, Siebert, and Lumsdaine (2010). Requires MPI-3. 2474 2475 Level: developer 2476 2477 .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedSetType(), PetscCommBuildTwoSidedGetType() 2478 E*/ 2479 typedef enum { 2480 PETSC_BUILDTWOSIDED_NOTSET = -1, 2481 PETSC_BUILDTWOSIDED_ALLREDUCE = 0, 2482 PETSC_BUILDTWOSIDED_IBARRIER = 1 2483 /* Updates here must be accompanied by updates in finclude/petscsys.h and the string array in mpits.c */ 2484 } PetscBuildTwoSidedType; 2485 PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[]; 2486 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType); 2487 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*); 2488 2489 PETSC_EXTERN PetscErrorCode PetscSSEIsEnabled(MPI_Comm,PetscBool *,PetscBool *); 2490 2491 /*E 2492 InsertMode - Whether entries are inserted or added into vectors or matrices 2493 2494 Level: beginner 2495 2496 .seealso: VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(), 2497 VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), 2498 MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd() 2499 E*/ 2500 typedef enum {NOT_SET_VALUES, INSERT_VALUES, ADD_VALUES, MAX_VALUES, INSERT_ALL_VALUES, ADD_ALL_VALUES, INSERT_BC_VALUES, ADD_BC_VALUES} InsertMode; 2501 2502 /*MC 2503 INSERT_VALUES - Put a value into a vector or matrix, overwrites any previous value 2504 2505 Level: beginner 2506 2507 .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(), 2508 VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), ADD_VALUES, 2509 MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES 2510 2511 M*/ 2512 2513 /*MC 2514 ADD_VALUES - Adds a value into a vector or matrix, if there previously was no value, just puts the 2515 value into that location 2516 2517 Level: beginner 2518 2519 .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(), 2520 VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), INSERT_VALUES, 2521 MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES 2522 2523 M*/ 2524 2525 /*MC 2526 MAX_VALUES - Puts the maximum of the scattered/gathered value and the current value into each location 2527 2528 Level: beginner 2529 2530 .seealso: InsertMode, VecScatterBegin(), VecScatterEnd(), ADD_VALUES, INSERT_VALUES 2531 2532 M*/ 2533 2534 PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject); 2535 2536 typedef enum {PETSC_SUBCOMM_GENERAL=0,PETSC_SUBCOMM_CONTIGUOUS=1,PETSC_SUBCOMM_INTERLACED=2} PetscSubcommType; 2537 PETSC_EXTERN const char *const PetscSubcommTypes[]; 2538 2539 /*S 2540 PetscSubcomm - Context of MPI subcommunicators, used by PCREDUNDANT 2541 2542 Level: advanced 2543 2544 Concepts: communicator, create 2545 S*/ 2546 typedef struct _n_PetscSubcomm* PetscSubcomm; 2547 2548 struct _n_PetscSubcomm { 2549 MPI_Comm parent; /* parent communicator */ 2550 MPI_Comm dupparent; /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */ 2551 MPI_Comm comm; /* this communicator */ 2552 PetscMPIInt n; /* num of subcommunicators under the parent communicator */ 2553 PetscMPIInt color; /* color of processors belong to this communicator */ 2554 PetscMPIInt *subsize; /* size of subcommunicator[color] */ 2555 PetscSubcommType type; 2556 }; 2557 2558 PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*); 2559 PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*); 2560 PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt); 2561 PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType); 2562 PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt); 2563 PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer); 2564 PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm); 2565 2566 /*S 2567 PetscSegBuffer - a segmented extendable buffer 2568 2569 Level: developer 2570 2571 .seealso: PetscSegBufferCreate(), PetscSegBufferGet(), PetscSegBufferExtract(), PetscSegBufferDestroy() 2572 S*/ 2573 typedef struct _n_PetscSegBuffer *PetscSegBuffer; 2574 PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*); 2575 PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*); 2576 PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*); 2577 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*); 2578 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*); 2579 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*); 2580 PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*); 2581 PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t); 2582 2583 /* Type-safe wrapper to encourage use of PETSC_RESTRICT. Does not use PetscFunctionBegin because the error handling 2584 * prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as 2585 * possible. */ 2586 PETSC_STATIC_INLINE PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg,PetscInt count,PetscInt *PETSC_RESTRICT *slot) {return PetscSegBufferGet(seg,count,(void**)slot);} 2587 2588 extern PetscSegBuffer PetscCitationsList; 2589 #undef __FUNCT__ 2590 #define __FUNCT__ "PetscCitationsRegister" 2591 /*@C 2592 PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code. 2593 2594 Not Collective - only what is registered on rank 0 of PETSC_COMM_WORLD will be printed 2595 2596 Input Parameters: 2597 + cite - the bibtex item, formated to displayed on multiple lines nicely 2598 - set - a boolean variable initially set to PETSC_FALSE; this is used to insure only a single registration of the citation 2599 2600 Level: intermediate 2601 2602 Options Database: 2603 . -citations [filenmae] - print out the bibtex entries for the given computation 2604 @*/ 2605 PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set) 2606 { 2607 size_t len; 2608 char *vstring; 2609 PetscErrorCode ierr; 2610 2611 PetscFunctionBegin; 2612 if (set && *set) PetscFunctionReturn(0); 2613 ierr = PetscStrlen(cit,&len);CHKERRQ(ierr); 2614 ierr = PetscSegBufferGet(PetscCitationsList,(PetscInt)len,&vstring);CHKERRQ(ierr); 2615 ierr = PetscMemcpy(vstring,cit,len);CHKERRQ(ierr); 2616 if (set) *set = PETSC_TRUE; 2617 PetscFunctionReturn(0); 2618 } 2619 2620 2621 /* Reset __FUNCT__ in case the user does not define it themselves */ 2622 #undef __FUNCT__ 2623 #define __FUNCT__ "User provided function" 2624 2625 #endif 2626