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