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