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