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