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