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