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