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