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