Actual source code: petscsys.h

petsc-3.9.0 2018-04-07
Report Typos and Errors
  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: */
  7: /* ========================================================================== */
  8: /*
  9:    petscconf.h is contained in ${PETSC_ARCH}/include/petscconf.h it is
 10:    found automatically by the compiler due to the -I${PETSC_DIR}/${PETSC_ARCH}/include.
 11:    For --prefix installs the ${PETSC_ARCH}/ does not exist and petscconf.h is in the same
 12:    directory as the other PETSc include files.
 13: */
 14: #include <petscconf.h>
 15: #include <petscfix.h>

 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

 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

 46: /* ========================================================================== */
 47: /*
 48:    Since PETSc manages its own extern "C" handling users should never include PETSc include
 49:    files within extern "C". This will generate a compiler error if a user does put the include
 50:    file within an extern "C".
 51: */
 52: #if defined(__cplusplus)
 53: void assert_never_put_petsc_headers_inside_an_extern_c(int); void assert_never_put_petsc_headers_inside_an_extern_c(double);
 54: #endif

 56: #if defined(__cplusplus)
 57: #  define PETSC_RESTRICT PETSC_CXX_RESTRICT
 58: #else
 59: #  define PETSC_RESTRICT PETSC_C_RESTRICT
 60: #endif

 62: #if defined(__cplusplus)
 63: #  define PETSC_STATIC_INLINE PETSC_CXX_STATIC_INLINE
 64: #else
 65: #  define PETSC_STATIC_INLINE PETSC_C_STATIC_INLINE
 66: #endif

 68: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES) /* For Win32 shared libraries */
 69: #  define  __declspec(dllexport)
 70: #  define PETSC_DLLIMPORT __declspec(dllimport)
 71: #  define PETSC_VISIBILITY_INTERNAL
 72: #elif defined(PETSC_USE_VISIBILITY_CXX) && defined(__cplusplus)
 73: #  define  __attribute__((visibility ("default")))
 74: #  define PETSC_DLLIMPORT __attribute__((visibility ("default")))
 75: #  define PETSC_VISIBILITY_INTERNAL __attribute__((visibility ("hidden")))
 76: #elif defined(PETSC_USE_VISIBILITY_C) && !defined(__cplusplus)
 77: #  define  __attribute__((visibility ("default")))
 78: #  define PETSC_DLLIMPORT __attribute__((visibility ("default")))
 79: #  define PETSC_VISIBILITY_INTERNAL __attribute__((visibility ("hidden")))
 80: #else
 81: #  define 
 82: #  define PETSC_DLLIMPORT
 83: #  define PETSC_VISIBILITY_INTERNAL
 84: #endif

 86: #if defined(petsc_EXPORTS)      /* CMake defines this when building the shared library */
 87: #  define PETSC_VISIBILITY_PUBLIC 
 88: #else  /* Win32 users need this to import symbols from petsc.dll */
 89: #  define PETSC_VISIBILITY_PUBLIC PETSC_DLLIMPORT
 90: #endif

 92: /*
 93:     Functions tagged with PETSC_EXTERN in the header files are
 94:   always defined as extern "C" when compiled with C++ so they may be
 95:   used from C and are always visible in the shared libraries
 96: */
 97: #if defined(__cplusplus)
 98: #define PETSC_EXTERN extern "C" PETSC_VISIBILITY_PUBLIC
 99: #define PETSC_EXTERN_TYPEDEF extern "C"
100: #define PETSC_INTERN extern "C" PETSC_VISIBILITY_INTERNAL
101: #else
102: #define PETSC_EXTERN extern PETSC_VISIBILITY_PUBLIC
103: #define PETSC_EXTERN_TYPEDEF
104: #define PETSC_INTERN extern PETSC_VISIBILITY_INTERNAL
105: #endif

107:  #include <petscversion.h>
108: #define PETSC_AUTHOR_INFO  "       The PETSc Team\n    petsc-maint@mcs.anl.gov\n http://www.mcs.anl.gov/petsc/\n"

110: /* ========================================================================== */

112: /*
113:     Defines the interface to MPI allowing the use of all MPI functions.

115:     PETSc does not use the C++ binding of MPI at ALL. The following flag
116:     makes sure the C++ bindings are not included. The C++ bindings REQUIRE
117:     putting mpi.h before ANY C++ include files, we cannot control this
118:     with all PETSc users. Users who want to use the MPI C++ bindings can include
119:     mpicxx.h directly in their code
120: */
121: #if !defined(MPICH_SKIP_MPICXX)
122: #  define MPICH_SKIP_MPICXX 1
123: #endif
124: #if !defined(OMPI_SKIP_MPICXX)
125: #  define OMPI_SKIP_MPICXX 1
126: #endif
127: #include <mpi.h>

129: /*
130:    Perform various sanity checks that the correct mpi.h is being included at compile time.
131:    This usually happens because
132:       * either an unexpected mpi.h is in the default compiler path (i.e. in /usr/include) or
133:       * an extra include path -I/something (which contains the unexpected mpi.h) is being passed to the compiler
134: */
135: #if defined(PETSC_HAVE_MPIUNI)
136: #  if !defined(__MPIUNI_H)
137: #    error "PETSc was configured with --with-mpi=0 but now appears to be compiling using a different mpi.h"
138: #  endif
139: #elif defined(PETSC_HAVE_I_MPI_NUMVERSION)
140: #  if !defined(I_MPI_NUMVERSION)
141: #    error "PETSc was configured with I_MPI but now appears to be compiling using a non-I_MPI mpi.h"
142: #  elif I_MPI_NUMVERSION != PETSC_HAVE_I_MPI_NUMVERSION
143: #    error "PETSc was configured with one I_MPI mpi.h version but now appears to be compiling using a different I_MPI mpi.h version"
144: #  endif
145: #elif defined(PETSC_HAVE_MVAPICH2_NUMVERSION)
146: #  if !defined(MVAPICH2_NUMVERSION)
147: #    error "PETSc was configured with MVAPICH2 but now appears to be compiling using a non-MVAPICH2 mpi.h"
148: #  elif MVAPICH2_NUMVERSION != PETSC_HAVE_MVAPICH2_NUMVERSION
149: #    error "PETSc was configured with one MVAPICH2 mpi.h version but now appears to be compiling using a different MVAPICH2 mpi.h version"
150: #  endif
151: #elif defined(PETSC_HAVE_MPICH_NUMVERSION)
152: #  if !defined(MPICH_NUMVERSION) || defined(MVAPICH2_NUMVERSION) || defined(I_MPI_NUMVERSION)
153: #    error "PETSc was configured with MPICH but now appears to be compiling using a non-MPICH mpi.h"
154: #  elif (MPICH_NUMVERSION/100000 != PETSC_HAVE_MPICH_NUMVERSION/100000) || (MPICH_NUMVERSION%100000/1000 < PETSC_HAVE_MPICH_NUMVERSION%100000/1000)
155: #    error "PETSc was configured with one MPICH mpi.h version but now appears to be compiling using a different MPICH mpi.h version"
156: #  endif
157: #elif defined(PETSC_HAVE_OMPI_MAJOR_VERSION)
158: #  if !defined(OMPI_MAJOR_VERSION)
159: #    error "PETSc was configured with OpenMPI but now appears to be compiling using a non-OpenMPI mpi.h"
160: #  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)
161: #    error "PETSc was configured with one OpenMPI mpi.h version but now appears to be compiling using a different OpenMPI mpi.h version"
162: #  endif
163: #elif defined(OMPI_MAJOR_VERSION) || defined(MPICH_NUMVERSION)
164: #  error "PETSc was configured with undetermined MPI - but now appears to be compiling using either of OpenMPI or a MPICH variant"
165: #endif

167: /*
168:     Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler
169:     see the top of mpicxx.h in the MPICH2 distribution.
170: */
171: #include <stdio.h>

173: /* MSMPI on 32bit windows requires this yukky hack - that breaks MPI standard compliance */
174: #if !defined(MPIAPI)
175: #define MPIAPI
176: #endif

178: /*
179:    Support for Clang (>=3.2) matching type tag arguments with void* buffer types.
180:    This allows the compiler to detect cases where the MPI datatype argument passed to a MPI routine
181:    does not match the actual type of the argument being passed in
182: */
183: #if defined(__has_attribute) && defined(works_with_const_which_is_not_true)
184: #  if __has_attribute(argument_with_type_tag) && __has_attribute(pointer_with_type_tag) && __has_attribute(type_tag_for_datatype)
185: #    define PetscAttrMPIPointerWithType(bufno,typeno) __attribute__((pointer_with_type_tag(MPI,bufno,typeno)))
186: #    define PetscAttrMPITypeTag(type)                 __attribute__((type_tag_for_datatype(MPI,type)))
187: #    define PetscAttrMPITypeTagLayoutCompatible(type) __attribute__((type_tag_for_datatype(MPI,type,layout_compatible)))
188: #  endif
189: #endif
190: #if !defined(PetscAttrMPIPointerWithType)
191: #  define PetscAttrMPIPointerWithType(bufno,typeno)
192: #  define PetscAttrMPITypeTag(type)
193: #  define PetscAttrMPITypeTagLayoutCompatible(type)
194: #endif

196: /*MC
197:     PetscErrorCode - datatype used for return error code from almost all PETSc functions

199:     Level: beginner

201: .seealso: CHKERRQ, SETERRQ
202: M*/
203: typedef int PetscErrorCode;

205: /*MC

207:     PetscClassId - A unique id used to identify each PETSc class.

209:     Notes:
210:     Use PetscClassIdRegister() to obtain a new value for a new class being created. Usually
211:          XXXInitializePackage() calls it for each class it defines.

213:     Developer Notes:
214:     Internal integer stored in the _p_PetscObject data structure.
215:          These are all computed by an offset from the lowest one, PETSC_SMALLEST_CLASSID.

217:     Level: developer

219: .seealso: PetscClassIdRegister(), PetscLogEventRegister(), PetscHeaderCreate()
220: M*/
221: typedef int PetscClassId;


224: /*MC
225:     PetscMPIInt - datatype used to represent 'int' parameters to MPI functions.

227:     Level: intermediate

229:     Notes:
230:     usually this is the same as PetscInt, but if PETSc was built with --with-64-bit-indices but
231:            standard C/Fortran integers are 32 bit then this is NOT the same as PetscInt; it remains 32 bit.

233:     PetscMPIIntCast(a,&b) checks if the given PetscInt a will fit in a PetscMPIInt, if not it
234:       generates a PETSC_ERR_ARG_OUTOFRANGE error.

236: .seealso: PetscBLASInt, PetscInt, PetscMPIIntCast()

238: M*/
239: typedef int PetscMPIInt;

241: /*MC
242:     PetscEnum - datatype used to pass enum types within PETSc functions.

244:     Level: intermediate

246: .seealso: PetscOptionsGetEnum(), PetscOptionsEnum(), PetscBagRegisterEnum()
247: M*/
248: typedef enum { ENUM_DUMMY } PetscEnum;
249: PETSC_EXTERN MPI_Datatype MPIU_ENUM PetscAttrMPITypeTag(PetscEnum);

251: #if defined(PETSC_HAVE_STDINT_H)
252: #include <stdint.h>
253: #endif
254: #if defined (PETSC_HAVE_INTTYPES_H)
256: #include <inttypes.h>
257: # if !defined(PRId64)
258: # define PRId64 "ld"
259: # endif
260: #endif

262: typedef short PetscShort;
263: typedef char PetscChar;
264: typedef float PetscFloat;

266: /*MC
267:   PetscInt - PETSc type that represents an integer, used primarily to
268:       represent size of arrays and indexing into arrays. Its size can be configured with the option --with-64-bit-indices to be either 32-bit (default) or 64-bit.

270:   Notes:
271:   For MPI calls that require datatypes, use MPIU_INT as the datatype for PetscInt. It will automatically work correctly regardless of the size of PetscInt.

273:   Level: beginner

275: .seealso: PetscBLASInt, PetscMPIInt, PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_SCALAR, MPIU_COMPLEX, MPIU_INT
276: M*/

278: /*MC
279:    MPIU_INT - MPI datatype corresponding to PetscInt

281:    Notes:
282:    In MPI calls that require an MPI datatype that matches a PetscInt or array of PetscInt values, pass this value.

284:    Level: beginner

286: .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_SCALAR, MPIU_COMPLEX
287: M*/

289: #if defined(PETSC_HAVE_STDINT_H) && defined(PETSC_HAVE_INTTYPES_H) && defined(PETSC_HAVE_MPI_INT64_T) /* MPI_INT64_T is not guaranteed to be a macro */
290: typedef int64_t PetscInt64;
291: # define MPIU_INT64 MPI_INT64_T
292: # define PetscInt64_FMT PRId64
293: #elif (PETSC_SIZEOF_LONG_LONG == 8)
294: typedef long long PetscInt64;
295: # define MPIU_INT64 MPI_LONG_LONG_INT
296: # define PetscInt64_FMT "lld"
297: #elif defined(PETSC_HAVE___INT64)
298: typedef __int64 PetscInt64;
299: # define MPIU_INT64 MPI_INT64_T
300: # define PetscInt64_FMT "ld"
301: #else
302: #error "cannot determine PetscInt64 type"
303: #endif

305: #if PETSC_SIZEOF_VOID_P == 4
306: #define MPIU_FORTRANADDR MPI_INT
307: #else
308: #define MPIU_FORTRANADDR MPIU_INT64
309: #endif

311: #if defined(PETSC_USE_64BIT_INDICES)
312: typedef PetscInt64 PetscInt;
313: #define MPIU_INT MPIU_INT64
314: #define PetscInt_FMT PetscInt64_FMT
315: #else
316: typedef int PetscInt;
317: #define MPIU_INT MPI_INT
318: #define PetscInt_FMT "d"
319: #endif

321: /*MC
322:    PetscBLASInt - datatype used to represent 'int' parameters to BLAS/LAPACK functions.

324:    Notes:
325:     Usually this is the same as PetscInt, but if PETSc was built with --with-64-bit-indices but
326:            standard C/Fortran integers are 32 bit then this is NOT the same as PetscInt it remains 32 bit
327:            (except on very rare BLAS/LAPACK implementations that support 64 bit integers see the note below).

329:     PetscErrorCode PetscBLASIntCast(a,&b) checks if the given PetscInt a will fit in a PetscBLASInt, if not it
330:       generates a PETSC_ERR_ARG_OUTOFRANGE error

332:    Installation Notes:
333:     The 64bit versions of MATLAB ship with BLAS and LAPACK that use 64 bit integers for sizes etc,
334:      if you run ./configure with the option
335:      --with-blaslapack-lib=[/Applications/MATLAB_R2010b.app/bin/maci64/libmwblas.dylib,/Applications/MATLAB_R2010b.app/bin/maci64/libmwlapack.dylib]
336:      but you need to also use --known-64-bit-blas-indices.

338:         MKL also ships with 64 bit integer versions of the BLAS and LAPACK, if you select those you must also ./configure with
339:         --known-64-bit-blas-indices

341:         OpenBLAS can also be built to use 64 bit integers. The ./configure options --download-openblas -download-openblas-64-bit-blas-indices
342:         will build a 64 bit integer version

344:     Developer Notes:
345:      Eventually ./configure should automatically determine the size of the integers used by BLAS/LAPACK.

347:      External packages such as hypre, ML, SuperLU etc do not provide any support for passing 64 bit integers to BLAS/LAPACK so cannot
348:      be used with PETSc if you have set PetscBLASInt to long int.

350:    Level: intermediate

352: .seealso: PetscMPIInt, PetscInt, PetscBLASIntCast()

354: M*/
355: #if defined(PETSC_HAVE_64BIT_BLAS_INDICES)
356: typedef PetscInt64 PetscBLASInt;
357: #else
358: typedef int PetscBLASInt;
359: #endif

361: /*
362:     For the rare cases when one needs to send a size_t object with MPI
363: */
364: #if (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_INT)
365: #define MPIU_SIZE_T MPI_UNSIGNED
366: #elif  (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_LONG)
367: #define MPIU_SIZE_T MPI_UNSIGNED_LONG
368: #elif  (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_LONG_LONG)
369: #define MPIU_SIZE_T MPI_UNSIGNED_LONG_LONG
370: #else
371: #error "Unknown size for size_t! Send us a bugreport at petsc-maint@mcs.anl.gov"
372: #endif

374: /*
375:       You can use PETSC_STDOUT as a replacement of stdout. You can also change
376:     the value of PETSC_STDOUT to redirect all standard output elsewhere
377: */
378: PETSC_EXTERN FILE* PETSC_STDOUT;

380: /*
381:       You can use PETSC_STDERR as a replacement of stderr. You can also change
382:     the value of PETSC_STDERR to redirect all standard error elsewhere
383: */
384: PETSC_EXTERN FILE* PETSC_STDERR;

386: /*MC
387:     PetscUnlikely - hints the compiler that the given condition is usually FALSE

389:     Synopsis:
390:     #include <petscsys.h>
391:     PetscBool  PetscUnlikely(PetscBool  cond)

393:     Not Collective

395:     Input Parameters:
396: .   cond - condition or expression

398:     Notes:
399:     This returns the same truth value, it is only a hint to compilers that the resulting
400:     branch is unlikely.

402:     Level: advanced

404: .seealso: PetscLikely(), CHKERRQ
405: M*/

407: /*MC
408:     PetscLikely - hints the compiler that the given condition is usually TRUE

410:     Synopsis:
411:     #include <petscsys.h>
412:     PetscBool  PetscLikely(PetscBool  cond)

414:     Not Collective

416:     Input Parameters:
417: .   cond - condition or expression

419:     Notes:
420:     This returns the same truth value, it is only a hint to compilers that the resulting
421:     branch is likely.

423:     Level: advanced

425: .seealso: PetscUnlikely()
426: M*/
427: #if defined(PETSC_HAVE_BUILTIN_EXPECT)
428: #  define PetscUnlikely(cond)   __builtin_expect(!!(cond),0)
429: #  define PetscLikely(cond)     __builtin_expect(!!(cond),1)
430: #else
431: #  define PetscUnlikely(cond)   (cond)
432: #  define PetscLikely(cond)     (cond)
433: #endif

435: /*
436:     Declare extern C stuff after including external header files
437: */


440: /*
441:        Basic PETSc constants
442: */

444: /*E
445:     PetscBool  - Logical variable. Actually an int in C and a logical in Fortran.

447:    Level: beginner

449:    Developer Note:
450:    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
451:       boolean values. It is not easy to have a simple macro that that will work properly in all circumstances with all three mechanisms.

453: .seealso: PETSC_TRUE, PETSC_FALSE, PetscNot()
454: E*/
455: typedef enum { PETSC_FALSE,PETSC_TRUE } PetscBool;
456: PETSC_EXTERN const char *const PetscBools[];
457: PETSC_EXTERN MPI_Datatype MPIU_BOOL PetscAttrMPITypeTag(PetscBool);

459: /*
460:     Defines elementary mathematics functions and constants.
461: */
462:  #include <petscmath.h>

464: /*E
465:     PetscCopyMode  - Determines how an array passed to certain functions is copied or retained

467:    Level: beginner

469: $   PETSC_COPY_VALUES - the array values are copied into new space, the user is free to reuse or delete the passed in array
470: $   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
471: $                       delete the array. The array MUST have been obtained with PetscMalloc(). Hence this mode cannot be used in Fortran.
472: $   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
473:                         the array but the user must delete the array after the object is destroyed.

475: E*/
476: typedef enum { PETSC_COPY_VALUES, PETSC_OWN_POINTER, PETSC_USE_POINTER} PetscCopyMode;
477: PETSC_EXTERN const char *const PetscCopyModes[];

479: /*MC
480:     PETSC_FALSE - False value of PetscBool

482:     Level: beginner

484:     Note:
485:     Zero integer

487: .seealso: PetscBool, PETSC_TRUE
488: M*/

490: /*MC
491:     PETSC_TRUE - True value of PetscBool

493:     Level: beginner

495:     Note:
496:     Nonzero integer

498: .seealso: PetscBool, PETSC_FALSE
499: M*/

501: /*MC
502:     PETSC_IGNORE - same as NULL, means PETSc will ignore this argument

504:    Level: beginner

506:    Note:
507:    Accepted by many PETSc functions to not set a parameter and instead use
508:           some default

510:    Fortran Notes:
511:    This macro does not exist in Fortran; you must use PETSC_NULL_INTEGER,
512:           PETSC_NULL_DOUBLE_PRECISION etc

514: .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_DETERMINE

516: M*/
517: #define PETSC_IGNORE         NULL

519: /* This is deprecated */
520: #define PETSC_NULL NULL

522: /*MC
523:     PETSC_DECIDE - standard way of passing in integer or floating point parameter
524:        where you wish PETSc to use the default.

526:    Level: beginner

528: .seealso: PETSC_DEFAULT, PETSC_IGNORE, PETSC_DETERMINE

530: M*/
531: #define PETSC_DECIDE  -1

533: /*MC
534:     PETSC_DETERMINE - standard way of passing in integer or floating point parameter
535:        where you wish PETSc to compute the required value.

537:    Level: beginner

539:    Developer Note:
540:    I would like to use const PetscInt PETSC_DETERMINE = PETSC_DECIDE; but for
541:      some reason this is not allowed by the standard even though PETSC_DECIDE is a constant value.

543: .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_IGNORE, VecSetSizes()

545: M*/
546: #define PETSC_DETERMINE PETSC_DECIDE

548: /*MC
549:     PETSC_DEFAULT - standard way of passing in integer or floating point parameter
550:        where you wish PETSc to use the default.

552:    Level: beginner

554:    Fortran Notes:
555:    You need to use PETSC_DEFAULT_INTEGER or PETSC_DEFAULT_REAL.

557: .seealso: PETSC_DECIDE, PETSC_IGNORE, PETSC_DETERMINE

559: M*/
560: #define PETSC_DEFAULT  -2

562: /*MC
563:     PETSC_COMM_WORLD - the equivalent of the MPI_COMM_WORLD communicator which represents
564:            all the processs that PETSc knows about.

566:    Level: beginner

568:    Notes:
569:    By default PETSC_COMM_WORLD and MPI_COMM_WORLD are identical unless you wish to
570:           run PETSc on ONLY a subset of MPI_COMM_WORLD. In that case create your new (smaller)
571:           communicator, call it, say comm, and set PETSC_COMM_WORLD = comm BEFORE calling
572:           PetscInitialize(), but after MPI_Init() has been called.

574:           The value of PETSC_COMM_WORLD should never be USED/accessed before PetscInitialize()
575:           is called because it may not have a valid value yet.

577: .seealso: PETSC_COMM_SELF

579: M*/
580: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;

582: /*MC
583:     PETSC_COMM_SELF - This is always MPI_COMM_SELF

585:    Level: beginner

587:    Notes:
588:    Do not USE/access or set this variable before PetscInitialize() has been called.

590: .seealso: PETSC_COMM_WORLD

592: M*/
593: #define PETSC_COMM_SELF MPI_COMM_SELF

595: PETSC_EXTERN PetscBool PetscBeganMPI;
596: PETSC_EXTERN PetscBool PetscInitializeCalled;
597: PETSC_EXTERN PetscBool PetscFinalizeCalled;
598: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;
599: PETSC_EXTERN PetscBool PetscCUDASynchronize;

601: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm));
602: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*);
603: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm*);

605: /*MC
606:    PetscMalloc - Allocates memory, One should use PetscNew(), PetscMalloc1() or PetscCalloc1() usually instead of this

608:    Synopsis:
609:     #include <petscsys.h>
610:    PetscErrorCode PetscMalloc(size_t m,void **result)

612:    Not Collective

614:    Input Parameter:
615: .  m - number of bytes to allocate

617:    Output Parameter:
618: .  result - memory allocated

620:    Level: beginner

622:    Notes:
623:    Memory is always allocated at least double aligned

625:    It is safe to allocate size 0 and pass the resulting pointer (which may or may not be NULL) to PetscFree().

627: .seealso: PetscFree(), PetscNew()

629:   Concepts: memory allocation

631: M*/
632: #define PetscMalloc(a,b)  ((*PetscTrMalloc)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__,(void**)(b)))

634: /*MC
635:    PetscRealloc - Rellocates memory

637:    Synopsis:
638:     #include <petscsys.h>
639:    PetscErrorCode PetscRealloc(size_t m,void **result)

641:    Not Collective

643:    Input Parameters:
644: +  m - number of bytes to allocate
645: -  result - initial memory allocated

647:    Output Parameter:
648: .  result - new memory allocated

650:    Level: developer

652:    Notes:
653:    Memory is always allocated at least double aligned

655: .seealso: PetscMalloc(), PetscFree(), PetscNew()

657:   Concepts: memory allocation

659: M*/
660: #define PetscRealloc(a,b)  ((*PetscTrRealloc)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__,(void**)(b)))

662: /*MC
663:    PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment

665:    Synopsis:
666:     #include <petscsys.h>
667:    void *PetscAddrAlign(void *addr)

669:    Not Collective

671:    Input Parameters:
672: .  addr - address to align (any pointer type)

674:    Level: developer

676: .seealso: PetscMallocAlign()

678:   Concepts: memory allocation
679: M*/
680: #define PetscAddrAlign(a) (void*)((((PETSC_UINTPTR_T)(a))+(PETSC_MEMALIGN-1)) & ~(PETSC_MEMALIGN-1))

682: /*MC
683:    PetscMalloc1 - Allocates an array of memory aligned to PETSC_MEMALIGN

685:    Synopsis:
686:     #include <petscsys.h>
687:    PetscErrorCode PetscMalloc1(size_t m1,type **r1)

689:    Not Collective

691:    Input Parameter:
692: .  m1 - number of elements to allocate  (may be zero)

694:    Output Parameter:
695: .  r1 - memory allocated in first chunk

697:    Note:
698:    This uses the sizeof() of the memory type requested to determine the total memory to be allocated, therefore you should not
699:          multiply the number of elements requested by the sizeof() the type. For example use
700: $  PetscInt *id;
701: $  PetscMalloc1(10,&id);
702:         not
703: $  PetscInt *id;
704: $  PetscMalloc1(10*sizeof(PetscInt),&id);

706:         Does not zero the memory allocatd, used PetscCalloc1() to obtain memory that has been zeroed.

708:    Level: beginner

710: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc1(), PetscMalloc2()

712:   Concepts: memory allocation

714: M*/
715: #define PetscMalloc1(m1,r1) PetscMallocA(1,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1))

717: /*MC
718:    PetscCalloc1 - Allocates a cleared (zeroed) array of memory aligned to PETSC_MEMALIGN

720:    Synopsis:
721:     #include <petscsys.h>
722:    PetscErrorCode PetscCalloc1(size_t m1,type **r1)

724:    Not Collective

726:    Input Parameter:
727: .  m1 - number of elements to allocate in 1st chunk  (may be zero)

729:    Output Parameter:
730: .  r1 - memory allocated in first chunk

732:    Notes:
733:    See PetsMalloc1() for more details on usage.

735:    Level: beginner

737: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc1(), PetscCalloc2()

739:   Concepts: memory allocation

741: M*/
742: #define PetscCalloc1(m1,r1) PetscMallocA(1,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1))

744: /*MC
745:    PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN

747:    Synopsis:
748:     #include <petscsys.h>
749:    PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)

751:    Not Collective

753:    Input Parameter:
754: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
755: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

757:    Output Parameter:
758: +  r1 - memory allocated in first chunk
759: -  r2 - memory allocated in second chunk

761:    Level: developer

763: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc1(), PetscCalloc2()

765:   Concepts: memory allocation

767: M*/
768: #define PetscMalloc2(m1,r1,m2,r2) PetscMallocA(2,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2))

770: /*MC
771:    PetscCalloc2 - Allocates 2 cleared (zeroed) arrays of memory both aligned to PETSC_MEMALIGN

773:    Synopsis:
774:     #include <petscsys.h>
775:    PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)

777:    Not Collective

779:    Input Parameter:
780: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
781: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

783:    Output Parameter:
784: +  r1 - memory allocated in first chunk
785: -  r2 - memory allocated in second chunk

787:    Level: developer

789: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc1(), PetscMalloc2()

791:   Concepts: memory allocation
792: M*/
793: #define PetscCalloc2(m1,r1,m2,r2) PetscMallocA(2,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2))

795: /*MC
796:    PetscMalloc3 - Allocates 3 arrays of memory, all aligned to PETSC_MEMALIGN

798:    Synopsis:
799:     #include <petscsys.h>
800:    PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

802:    Not Collective

804:    Input Parameter:
805: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
806: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
807: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

809:    Output Parameter:
810: +  r1 - memory allocated in first chunk
811: .  r2 - memory allocated in second chunk
812: -  r3 - memory allocated in third chunk

814:    Level: developer

816: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc3(), PetscFree3()

818:   Concepts: memory allocation

820: M*/
821: #define PetscMalloc3(m1,r1,m2,r2,m3,r3) PetscMallocA(3,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3))

823: /*MC
824:    PetscCalloc3 - Allocates 3 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN

826:    Synopsis:
827:     #include <petscsys.h>
828:    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

830:    Not Collective

832:    Input Parameter:
833: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
834: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
835: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

837:    Output Parameter:
838: +  r1 - memory allocated in first chunk
839: .  r2 - memory allocated in second chunk
840: -  r3 - memory allocated in third chunk

842:    Level: developer

844: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc2(), PetscMalloc3(), PetscFree3()

846:   Concepts: memory allocation
847: M*/
848: #define PetscCalloc3(m1,r1,m2,r2,m3,r3) PetscMallocA(3,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3))

850: /*MC
851:    PetscMalloc4 - Allocates 4 arrays of memory, all aligned to PETSC_MEMALIGN

853:    Synopsis:
854:     #include <petscsys.h>
855:    PetscErrorCode PetscMalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)

857:    Not Collective

859:    Input Parameter:
860: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
861: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
862: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
863: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

865:    Output Parameter:
866: +  r1 - memory allocated in first chunk
867: .  r2 - memory allocated in second chunk
868: .  r3 - memory allocated in third chunk
869: -  r4 - memory allocated in fourth chunk

871:    Level: developer

873: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc4(), PetscFree4()

875:   Concepts: memory allocation

877: M*/
878: #define PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4) PetscMallocA(4,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4))

880: /*MC
881:    PetscCalloc4 - Allocates 4 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN

883:    Synopsis:
884:     #include <petscsys.h>
885:    PetscErrorCode PetscCalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)

887:    Not Collective

889:    Input Parameters:
890: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
891: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
892: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
893: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

895:    Output Parameters:
896: +  r1 - memory allocated in first chunk
897: .  r2 - memory allocated in second chunk
898: .  r3 - memory allocated in third chunk
899: -  r4 - memory allocated in fourth chunk

901:    Level: developer

903: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc4(), PetscFree4()

905:   Concepts: memory allocation

907: M*/
908: #define PetscCalloc4(m1,r1,m2,r2,m3,r3,m4,r4) PetscMallocA(4,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4))

910: /*MC
911:    PetscMalloc5 - Allocates 5 arrays of memory, all aligned to PETSC_MEMALIGN

913:    Synopsis:
914:     #include <petscsys.h>
915:    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)

917:    Not Collective

919:    Input Parameters:
920: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
921: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
922: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
923: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
924: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

926:    Output Parameters:
927: +  r1 - memory allocated in first chunk
928: .  r2 - memory allocated in second chunk
929: .  r3 - memory allocated in third chunk
930: .  r4 - memory allocated in fourth chunk
931: -  r5 - memory allocated in fifth chunk

933:    Level: developer

935: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc5(), PetscFree5()

937:   Concepts: memory allocation

939: M*/
940: #define PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5) PetscMallocA(5,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4),(size_t)(m5)*sizeof(**(r5)),(r5))

942: /*MC
943:    PetscCalloc5 - Allocates 5 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN

945:    Synopsis:
946:     #include <petscsys.h>
947:    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)

949:    Not Collective

951:    Input Parameters:
952: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
953: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
954: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
955: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
956: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

958:    Output Parameters:
959: +  r1 - memory allocated in first chunk
960: .  r2 - memory allocated in second chunk
961: .  r3 - memory allocated in third chunk
962: .  r4 - memory allocated in fourth chunk
963: -  r5 - memory allocated in fifth chunk

965:    Level: developer

967: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc5(), PetscFree5()

969:   Concepts: memory allocation

971: M*/
972: #define PetscCalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5) PetscMallocA(5,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4),(size_t)(m5)*sizeof(**(r5)),(r5))

974: /*MC
975:    PetscMalloc6 - Allocates 6 arrays of memory, all aligned to PETSC_MEMALIGN

977:    Synopsis:
978:     #include <petscsys.h>
979:    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)

981:    Not Collective

983:    Input Parameters:
984: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
985: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
986: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
987: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
988: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
989: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

991:    Output Parameteasr:
992: +  r1 - memory allocated in first chunk
993: .  r2 - memory allocated in second chunk
994: .  r3 - memory allocated in third chunk
995: .  r4 - memory allocated in fourth chunk
996: .  r5 - memory allocated in fifth chunk
997: -  r6 - memory allocated in sixth chunk

999:    Level: developer

1001: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc6(), PetscFree3(), PetscFree4(), PetscFree5(), PetscFree6()

1003:   Concepts: memory allocation

1005: M*/
1006: #define PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) PetscMallocA(6,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4),(size_t)(m5)*sizeof(**(r5)),(r5),(size_t)(m6)*sizeof(**(r6)),(r6))

1008: /*MC
1009:    PetscCalloc6 - Allocates 6 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN

1011:    Synopsis:
1012:     #include <petscsys.h>
1013:    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)

1015:    Not Collective

1017:    Input Parameters:
1018: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
1019: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
1020: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
1021: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
1022: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
1023: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

1025:    Output Parameters:
1026: +  r1 - memory allocated in first chunk
1027: .  r2 - memory allocated in second chunk
1028: .  r3 - memory allocated in third chunk
1029: .  r4 - memory allocated in fourth chunk
1030: .  r5 - memory allocated in fifth chunk
1031: -  r6 - memory allocated in sixth chunk

1033:    Level: developer

1035: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscMalloc6(), PetscFree6()

1037:   Concepts: memory allocation
1038: M*/
1039: #define PetscCalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) PetscMallocA(6,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4),(size_t)(m5)*sizeof(**(r5)),(r5),(size_t)(m6)*sizeof(**(r6)),(r6))

1041: /*MC
1042:    PetscMalloc7 - Allocates 7 arrays of memory, all aligned to PETSC_MEMALIGN

1044:    Synopsis:
1045:     #include <petscsys.h>
1046:    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)

1048:    Not Collective

1050:    Input Parameters:
1051: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
1052: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
1053: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
1054: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
1055: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
1056: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
1057: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

1059:    Output Parameters:
1060: +  r1 - memory allocated in first chunk
1061: .  r2 - memory allocated in second chunk
1062: .  r3 - memory allocated in third chunk
1063: .  r4 - memory allocated in fourth chunk
1064: .  r5 - memory allocated in fifth chunk
1065: .  r6 - memory allocated in sixth chunk
1066: -  r7 - memory allocated in seventh chunk

1068:    Level: developer

1070: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc7(), PetscFree7()

1072:   Concepts: memory allocation

1074: M*/
1075: #define PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) PetscMallocA(7,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4),(size_t)(m5)*sizeof(**(r5)),(r5),(size_t)(m6)*sizeof(**(r6)),(r6),(size_t)(m7)*sizeof(**(r7)),(r7))

1077: /*MC
1078:    PetscCalloc7 - Allocates 7 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN

1080:    Synopsis:
1081:     #include <petscsys.h>
1082:    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)

1084:    Not Collective

1086:    Input Parameters:
1087: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
1088: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
1089: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
1090: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
1091: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
1092: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
1093: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

1095:    Output Parameters:
1096: +  r1 - memory allocated in first chunk
1097: .  r2 - memory allocated in second chunk
1098: .  r3 - memory allocated in third chunk
1099: .  r4 - memory allocated in fourth chunk
1100: .  r5 - memory allocated in fifth chunk
1101: .  r6 - memory allocated in sixth chunk
1102: -  r7 - memory allocated in seventh chunk

1104:    Level: developer

1106: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscMalloc7(), PetscFree7()

1108:   Concepts: memory allocation
1109: M*/
1110: #define PetscCalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) PetscMallocA(7,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4),(size_t)(m5)*sizeof(**(r5)),(r5),(size_t)(m6)*sizeof(**(r6)),(r6),(size_t)(m7)*sizeof(**(r7)),(r7))

1112: /*MC
1113:    PetscNew - Allocates memory of a particular type, zeros the memory! Aligned to PETSC_MEMALIGN

1115:    Synopsis:
1116:     #include <petscsys.h>
1117:    PetscErrorCode PetscNew(type **result)

1119:    Not Collective

1121:    Output Parameter:
1122: .  result - memory allocated, sized to match pointer type

1124:    Level: beginner

1126: .seealso: PetscFree(), PetscMalloc(), PetscNewLog(), PetscCalloc1(), PetscMalloc1()

1128:   Concepts: memory allocation

1130: M*/
1131: #define PetscNew(b)      PetscCalloc1(1,(b))

1133: /*MC
1134:    PetscNewLog - Allocates memory of a type matching pointer, zeros the memory! Aligned to PETSC_MEMALIGN. Associates the memory allocated
1135:          with the given object using PetscLogObjectMemory().

1137:    Synopsis:
1138:     #include <petscsys.h>
1139:    PetscErrorCode PetscNewLog(PetscObject obj,type **result)

1141:    Not Collective

1143:    Input Parameter:
1144: .  obj - object memory is logged to

1146:    Output Parameter:
1147: .  result - memory allocated, sized to match pointer type

1149:    Level: developer

1151: .seealso: PetscFree(), PetscMalloc(), PetscNew(), PetscLogObjectMemory(), PetscCalloc1(), PetscMalloc1()

1153:   Concepts: memory allocation

1155: M*/
1156: #define PetscNewLog(o,b) (PetscNew((b)) || PetscLogObjectMemory((PetscObject)o,sizeof(**(b))))

1158: /*MC
1159:    PetscFree - Frees memory

1161:    Synopsis:
1162:     #include <petscsys.h>
1163:    PetscErrorCode PetscFree(void *memory)

1165:    Not Collective

1167:    Input Parameter:
1168: .   memory - memory to free (the pointer is ALWAYS set to NULL upon sucess)

1170:    Level: beginner

1172:    Notes:
1173:    Do not free memory obtained with PetscMalloc2(), PetscCalloc2() etc, they must be freed with PetscFree2() etc.

1175:    It is safe to call PetscFree() on a NULL pointer.

1177: .seealso: PetscNew(), PetscMalloc(), PetscNewLog(), PetscMalloc1(), PetscCalloc1()

1179:   Concepts: memory allocation

1181: M*/
1182: #define PetscFree(a)   ((*PetscTrFree)((void*)(a),__LINE__,PETSC_FUNCTION_NAME,__FILE__) || ((a) = 0,0))

1184: /*MC
1185:    PetscFree2 - Frees 2 chunks of memory obtained with PetscMalloc2()

1187:    Synopsis:
1188:     #include <petscsys.h>
1189:    PetscErrorCode PetscFree2(void *memory1,void *memory2)

1191:    Not Collective

1193:    Input Parameters:
1194: +   memory1 - memory to free
1195: -   memory2 - 2nd memory to free

1197:    Level: developer

1199:    Notes: Memory must have been obtained with PetscMalloc2()

1201: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree()

1203:   Concepts: memory allocation

1205: M*/
1206: #define PetscFree2(m1,m2)   PetscFreeA(2,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2))

1208: /*MC
1209:    PetscFree3 - Frees 3 chunks of memory obtained with PetscMalloc3()

1211:    Synopsis:
1212:     #include <petscsys.h>
1213:    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)

1215:    Not Collective

1217:    Input Parameters:
1218: +   memory1 - memory to free
1219: .   memory2 - 2nd memory to free
1220: -   memory3 - 3rd memory to free

1222:    Level: developer

1224:    Notes: Memory must have been obtained with PetscMalloc3()

1226: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3()

1228:   Concepts: memory allocation

1230: M*/
1231: #define PetscFree3(m1,m2,m3)   PetscFreeA(3,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3))

1233: /*MC
1234:    PetscFree4 - Frees 4 chunks of memory obtained with PetscMalloc4()

1236:    Synopsis:
1237:     #include <petscsys.h>
1238:    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)

1240:    Not Collective

1242:    Input Parameters:
1243: +   m1 - memory to free
1244: .   m2 - 2nd memory to free
1245: .   m3 - 3rd memory to free
1246: -   m4 - 4th memory to free

1248:    Level: developer

1250:    Notes: Memory must have been obtained with PetscMalloc4()

1252: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4()

1254:   Concepts: memory allocation

1256: M*/
1257: #define PetscFree4(m1,m2,m3,m4)   PetscFreeA(4,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4))

1259: /*MC
1260:    PetscFree5 - Frees 5 chunks of memory obtained with PetscMalloc5()

1262:    Synopsis:
1263:     #include <petscsys.h>
1264:    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)

1266:    Not Collective

1268:    Input Parameters:
1269: +   m1 - memory to free
1270: .   m2 - 2nd memory to free
1271: .   m3 - 3rd memory to free
1272: .   m4 - 4th memory to free
1273: -   m5 - 5th memory to free

1275:    Level: developer

1277:    Notes: Memory must have been obtained with PetscMalloc5()

1279: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5()

1281:   Concepts: memory allocation

1283: M*/
1284: #define PetscFree5(m1,m2,m3,m4,m5)   PetscFreeA(5,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5))

1286: /*MC
1287:    PetscFree6 - Frees 6 chunks of memory obtained with PetscMalloc6()

1289:    Synopsis:
1290:     #include <petscsys.h>
1291:    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)

1293:    Not Collective

1295:    Input Parameters:
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


1304:    Level: developer

1306:    Notes: Memory must have been obtained with PetscMalloc6()

1308: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6()

1310:   Concepts: memory allocation

1312: M*/
1313: #define PetscFree6(m1,m2,m3,m4,m5,m6)   PetscFreeA(6,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5),&(m6))

1315: /*MC
1316:    PetscFree7 - Frees 7 chunks of memory obtained with PetscMalloc7()

1318:    Synopsis:
1319:     #include <petscsys.h>
1320:    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)

1322:    Not Collective

1324:    Input Parameters:
1325: +   m1 - memory to free
1326: .   m2 - 2nd memory to free
1327: .   m3 - 3rd memory to free
1328: .   m4 - 4th memory to free
1329: .   m5 - 5th memory to free
1330: .   m6 - 6th memory to free
1331: -   m7 - 7th memory to free


1334:    Level: developer

1336:    Notes: Memory must have been obtained with PetscMalloc7()

1338: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1339:           PetscMalloc7()

1341:   Concepts: memory allocation

1343: M*/
1344: #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   PetscFreeA(7,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5),&(m6),&(m7))

1346: PETSC_EXTERN PetscErrorCode PetscMallocA(int,PetscBool,int,const char *,const char *,size_t,void *,...);
1347: PETSC_EXTERN PetscErrorCode PetscFreeA(int,int,const char *,const char *,void *,...);
1348: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,int,const char[],const char[],void**);
1349: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1350: PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t,int,const char[],const char[],void**);
1351: PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1352: PETSC_EXTERN PetscErrorCode PetscMallocSet(PetscErrorCode (*)(size_t,int,const char[],const char[],void**),PetscErrorCode (*)(void*,int,const char[],const char[]));
1353: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);

1355: /*
1356:   Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1357: */
1358: PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1359: PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);

1361: /*
1362:     PetscLogDouble variables are used to contain double precision numbers
1363:   that are not used in the numerical computations, but rather in logging,
1364:   timing etc.
1365: */
1366: typedef double PetscLogDouble;
1367: #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE

1369: /*
1370:    Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1371: */
1372: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1373: PETSC_EXTERN PetscErrorCode PetscMallocDumpLog(FILE *);
1374: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1375: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1376: PETSC_EXTERN PetscErrorCode PetscMallocDebug(PetscBool);
1377: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*);
1378: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1379: PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLog(void);
1380: PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLogThreshold(PetscLogDouble);
1381: PETSC_EXTERN PetscErrorCode PetscMallocGetDumpLog(PetscBool*);

1383: /*E
1384:     PetscDataType - Used for handling different basic data types.

1386:    Level: beginner

1388:    Notes:
1389:    Use of this should be avoided if one can directly use MPI_Datatype instead.

1391:    Developer comment:
1392:    It would be nice if we could always just use MPI Datatypes, why can we not?

1394: .seealso: PetscBinaryRead(), PetscBinaryWrite(), PetscDataTypeToMPIDataType(),
1395:           PetscDataTypeGetSize()

1397: E*/
1398: typedef enum {PETSC_INT = 0,PETSC_DOUBLE = 1,PETSC_COMPLEX = 2, PETSC_LONG = 3 ,PETSC_SHORT = 4,PETSC_FLOAT = 5,
1399:               PETSC_CHAR = 6,PETSC_BIT_LOGICAL = 7,PETSC_ENUM = 8,PETSC_BOOL=9, PETSC___FLOAT128 = 10,PETSC_OBJECT = 11, PETSC_FUNCTION = 12, PETSC_STRING = 13, PETSC___FP16 = 14,PETSC_STRUCT, PETSC_DATATYPE_UNKNOWN} PetscDataType;
1400: PETSC_EXTERN const char *const PetscDataTypes[];

1402: #if defined(PETSC_USE_COMPLEX)
1403: #define  PETSC_SCALAR  PETSC_COMPLEX
1404: #else
1405: #if defined(PETSC_USE_REAL_SINGLE)
1406: #define  PETSC_SCALAR  PETSC_FLOAT
1407: #elif defined(PETSC_USE_REAL___FLOAT128)
1408: #define  PETSC_SCALAR  PETSC___FLOAT128
1409: #elif defined(PETSC_USE_REAL___FP16)
1410: #define  PETSC_SCALAR  PETSC___FP16
1411: #else
1412: #define  PETSC_SCALAR  PETSC_DOUBLE
1413: #endif
1414: #endif
1415: #if defined(PETSC_USE_REAL_SINGLE)
1416: #define  PETSC_REAL  PETSC_FLOAT
1417: #elif defined(PETSC_USE_REAL___FLOAT128)
1418: #define  PETSC_REAL  PETSC___FLOAT128
1419: #elif defined(PETSC_USE_REAL___FP16)
1420: #define  PETSC_REAL  PETSC___FP16
1421: #else
1422: #define  PETSC_REAL  PETSC_DOUBLE
1423: #endif
1424: #define  PETSC_FORTRANADDR  PETSC_LONG

1426: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1427: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1428: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1429: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);

1431: /*
1432:     Basic memory and string operations. These are usually simple wrappers
1433:    around the basic Unix system calls, but a few of them have additional
1434:    functionality and/or error checking.
1435: */
1436: PETSC_EXTERN PetscErrorCode PetscBitMemcpy(void*,PetscInt,const void*,PetscInt,PetscInt,PetscDataType);
1437: PETSC_EXTERN PetscErrorCode PetscMemmove(void*,void *,size_t);
1438: PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1439: PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1440: PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***);
1441: PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1442: PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool  *);
1443: PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool  *);
1444: PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1445: PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1446: PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1447: PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1448: PETSC_EXTERN PetscErrorCode PetscStrlcat(char[],const char[],size_t);
1449: PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1450: PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1451: PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1452: PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1453: PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1454: PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1455: PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1456: PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1457: PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*);
1458: PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1459: PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1460: PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1461: PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1462: PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt,const char *const*,char***);
1463: PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt,char***);
1464: PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);

1466: PETSC_EXTERN void PetscStrcmpNoError(const char[],const char[],PetscBool  *);

1468: /*S
1469:     PetscToken - 'Token' used for managing tokenizing strings

1471:   Level: intermediate

1473: .seealso: PetscTokenCreate(), PetscTokenFind(), PetscTokenDestroy()
1474: S*/
1475: typedef struct _p_PetscToken* PetscToken;

1477: PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1478: PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1479: PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);
1480: PETSC_EXTERN PetscErrorCode PetscStrInList(const char[],const char[],char,PetscBool*);

1482: PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1483: PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);

1485: /*
1486:    These are MPI operations for MPI_Allreduce() etc
1487: */
1488: PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1489: #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1490: PETSC_EXTERN MPI_Op MPIU_SUM;
1491: #else
1492: #define MPIU_SUM MPI_SUM
1493: #endif
1494: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1495: PETSC_EXTERN MPI_Op MPIU_MAX;
1496: PETSC_EXTERN MPI_Op MPIU_MIN;
1497: #else
1498: #define MPIU_MAX MPI_MAX
1499: #define MPIU_MIN MPI_MIN
1500: #endif
1501: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);

1503: PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1504: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);

1506: /*S
1507:      PetscObject - any PETSc object, PetscViewer, Mat, Vec, KSP etc

1509:    Level: beginner

1511:    Note:
1512:    This is the base class from which all PETSc objects are derived from.

1514: .seealso:  PetscObjectDestroy(), PetscObjectView(), PetscObjectGetName(), PetscObjectSetName(), PetscObjectReference(), PetscObjectDereference()
1515: S*/
1516: typedef struct _p_PetscObject* PetscObject;

1518: /*MC
1519:     PetscObjectId - unique integer Id for a PetscObject

1521:     Level: developer

1523:     Notes:
1524:     Unlike pointer values, object ids are never reused.

1526: .seealso: PetscObjectState, PetscObjectGetId()
1527: M*/
1528: #if defined(PETSC_USING_F90) && !defined(PETSC_USE_FORTRANKIND) /* compaq F90 */
1529: typedef int PetscObjectId;
1530: #else
1531: typedef PetscInt64 PetscObjectId;
1532: #endif

1534: /*MC
1535:     PetscObjectState - integer state for a PetscObject

1537:     Level: developer

1539:     Notes:
1540:     Object state is always-increasing and (for objects that track state) can be used to determine if an object has
1541:     changed since the last time you interacted with it.  It is 64-bit so that it will not overflow for a very long time.

1543: .seealso: PetscObjectId, PetscObjectStateGet(), PetscObjectStateIncrease(), PetscObjectStateSet()
1544: M*/
1545: #if defined(PETSC_USING_F90) && !defined(PETSC_USE_FORTRANKIND) /* compaq F90 */
1546: typedef int PetscObjectState;
1547: #else
1548: typedef PetscInt64 PetscObjectState;
1549: #endif

1551: /*S
1552:      PetscFunctionList - Linked list of functions, possibly stored in dynamic libraries, accessed
1553:       by string name

1555:    Level: advanced

1557: .seealso:  PetscFunctionListAdd(), PetscFunctionListDestroy(), PetscOpFlist
1558: S*/
1559: typedef struct _n_PetscFunctionList *PetscFunctionList;

1561: /*E
1562:   PetscFileMode - Access mode for a file.

1564:   Level: beginner

1566: $  FILE_MODE_READ - open a file at its beginning for reading
1567: $  FILE_MODE_WRITE - open a file at its beginning for writing (will create if the file does not exist)
1568: $  FILE_MODE_APPEND - open a file at end for writing
1569: $  FILE_MODE_UPDATE - open a file for updating, meaning for reading and writing
1570: $  FILE_MODE_APPEND_UPDATE - open a file for updating, meaning for reading and writing, at the end

1572: .seealso: PetscViewerFileSetMode()
1573: E*/
1574: typedef enum {FILE_MODE_READ, FILE_MODE_WRITE, FILE_MODE_APPEND, FILE_MODE_UPDATE, FILE_MODE_APPEND_UPDATE} PetscFileMode;
1575: extern const char *const PetscFileModes[];

1577: /*
1578:     Defines PETSc error handling.
1579: */
1580:  #include <petscerror.h>

1582: #define PETSC_SMALLEST_CLASSID  1211211
1583: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1584: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1585: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);

1587: /*
1588:    Routines that get memory usage information from the OS
1589: */
1590: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1591: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1592: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1593: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);

1595: PETSC_EXTERN PetscErrorCode PetscInfoAllow(PetscBool ,const char []);
1596: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);

1598: /*
1599:    Initialization of PETSc
1600: */
1601: PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1602: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1603: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1604: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1605: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1606: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1607: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1608: PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1609: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1610: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);

1612: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1613: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);

1615: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1616: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1617: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1618: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);

1620: /*
1621:      These are so that in extern C code we can caste function pointers to non-extern C
1622:    function pointers. Since the regular C++ code expects its function pointers to be C++
1623: */
1624: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1625: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1626: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);

1628: /*
1629:     Functions that can act on any PETSc object.
1630: */
1631: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1632: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1633: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1634: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1635: PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1636: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1637: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1638: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1639: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1640: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1641: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1642: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1643: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1644: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1645: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt *);
1646: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1647: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1648: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject *);
1649: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void));
1650: #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d))
1651: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1652: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1653: PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1654: PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject,PetscObject);
1655: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt *);

1657:  #include <petscviewertypes.h>
1658:  #include <petscoptions.h>

1660: PETSC_EXTERN PetscErrorCode PetscObjectsListGetGlobalNumbering(MPI_Comm,PetscInt,PetscObject*,PetscInt*,PetscInt*);

1662: PETSC_EXTERN PetscErrorCode PetscMemoryShowUsage(PetscViewer,const char[]);
1663: PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer,const char[]);
1664: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1665: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1666: #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1667: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1668: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1669: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1670: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1671: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1672: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1673: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1674: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1675: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,PetscObject,const char[]);
1676: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1677: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1678: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject,const char[],PetscBool *);
1679: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1680: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1681: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);

1683: #if defined(PETSC_HAVE_SAWS)
1684: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1685: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1686: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1687: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1688: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1689: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1690: PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1691: PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1692: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1693: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);

1695: #else
1696: #define PetscSAWsBlock()                        0
1697: #define PetscObjectSAWsViewOff(obj)             0
1698: #define PetscObjectSAWsSetBlock(obj,flg)        0
1699: #define PetscObjectSAWsBlock(obj)               0
1700: #define PetscObjectSAWsGrantAccess(obj)         0
1701: #define PetscObjectSAWsTakeAccess(obj)          0
1702: #define PetscStackViewSAWs()                    0
1703: #define PetscStackSAWsViewOff()                 0
1704: #define PetscStackSAWsTakeAccess()
1705: #define PetscStackSAWsGrantAccess()

1707: #endif

1709: typedef void* PetscDLHandle;
1710: typedef enum {PETSC_DL_DECIDE=0,PETSC_DL_NOW=1,PETSC_DL_LOCAL=2} PetscDLMode;
1711: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1712: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1713: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);

1715: #if defined(PETSC_USE_DEBUG)
1716: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1717: #endif
1718: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);

1720: /*S
1721:      PetscObjectList - Linked list of PETSc objects, each accessable by string name

1723:    Level: developer

1725:    Notes:
1726:    Used by PetscObjectCompose() and PetscObjectQuery()

1728: .seealso:  PetscObjectListAdd(), PetscObjectListDestroy(), PetscObjectListFind(), PetscObjectCompose(), PetscObjectQuery(), PetscFunctionList
1729: S*/
1730: typedef struct _n_PetscObjectList *PetscObjectList;

1732: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1733: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1734: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1735: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1736: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1737: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);

1739: /*
1740:     Dynamic library lists. Lists of names of routines in objects or in dynamic
1741:   link libraries that will be loaded as needed.
1742: */

1744: #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1745: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1746: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1747: #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1748: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1749: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[]);
1750: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1751: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1752: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);

1754: /*S
1755:      PetscDLLibrary - Linked list of dynamics libraries to search for functions

1757:    Level: advanced

1759: .seealso:  PetscDLLibraryOpen()
1760: S*/
1761: typedef struct _n_PetscDLLibrary *PetscDLLibrary;
1762: PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1763: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1764: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1765: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1766: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1767: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1768: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1769: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);

1771: /*
1772:      Useful utility routines
1773: */
1774: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1775: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1776: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1777: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1778: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1779: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);
1780: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm,PetscInt[],PetscInt[]);
1781: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm,PetscReal[],PetscReal[]);

1783: /*
1784:     PetscNot - negates a logical type value and returns result as a PetscBool

1786:     Notes:
1787:     This is useful in cases like
1788: $     int        *a;
1789: $     PetscBool  flag = PetscNot(a)
1790:      where !a would not return a PetscBool because we cannot provide a cast from int to PetscBool in C.
1791: */
1792: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)

1794: /*MC
1795:     PetscHelpPrintf - Prints help messages.

1797:    Synopsis:
1798:     #include <petscsys.h>
1799:      PetscErrorCode (*PetscHelpPrintf)(const char format[],...);

1801:     Not Collective

1803:     Input Parameters:
1804: .   format - the usual printf() format string

1806:    Level: developer

1808:     Fortran Note:
1809:     This routine is not supported in Fortran.

1811:     Concepts: help messages^printing
1812:     Concepts: printing^help messages

1814: .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1815: M*/
1816: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);

1818: /*
1819:      Defines PETSc profiling.
1820: */
1821:  #include <petsclog.h>

1823: /*
1824:       Simple PETSc parallel IO for ASCII printing
1825: */
1826: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1827: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1828: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1829: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1830: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1831: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1832: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);

1834: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1835: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1836: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);

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*);
1841: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1842: #endif

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*[]);

1852: /*S
1853:      PetscContainer - Simple PETSc object that contains a pointer to any required data

1855:    Level: advanced

1857: .seealso:  PetscObject, PetscContainerCreate()
1858: S*/
1859: PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1860: typedef struct _p_PetscContainer*  PetscContainer;
1861: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void **);
1862: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void *);
1863: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1864: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer *);
1865: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void*));
1866: PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void*);

1868: /*
1869:    For use in debuggers
1870: */
1871: PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1872: PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1873: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1874: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1875: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);

1877: #include <stddef.h>
1878: #include <string.h>             /* for memcpy, memset */
1879: #if defined(PETSC_HAVE_STDLIB_H)
1880: #include <stdlib.h>
1881: #endif

1883: #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1884: #include <xmmintrin.h>
1885: #endif

1887: /*@C
1888:    PetscMemcpy - Copies n bytes, beginning at location b, to the space
1889:    beginning at location a. The two memory regions CANNOT overlap, use
1890:    PetscMemmove() in that case.

1892:    Not Collective

1894:    Input Parameters:
1895: +  b - pointer to initial memory space
1896: -  n - length (in bytes) of space to copy

1898:    Output Parameter:
1899: .  a - pointer to copy space

1901:    Level: intermediate

1903:    Compile Option:
1904:     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1905:                                   for memory copies on double precision values.
1906:     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1907:                                   for memory copies on double precision values.
1908:     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1909:                                   for memory copies on double precision values.

1911:    Note:
1912:    This routine is analogous to memcpy().

1914:    Not available from Fortran

1916:    Developer Note: this is inlined for fastest performance

1918:   Concepts: memory^copying
1919:   Concepts: copying^memory

1921: .seealso: PetscMemmove(), PetscStrallocpy()

1923: @*/
1924: PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n)
1925: {
1926: #if defined(PETSC_USE_DEBUG)
1927:   size_t al = (size_t) a,bl = (size_t) b;
1928:   size_t nl = (size_t) n;
1930:   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1931:   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1932: #else
1934: #endif
1935:   if (a != b && n > 0) {
1936: #if defined(PETSC_USE_DEBUG)
1937:     if ((al > bl && (al - bl) < nl) || (bl - al) < nl)  SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1938:               or make sure your copy regions and lengths are correct. \n\
1939:               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1940: #endif
1941: #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1942:    if (!(a % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1943:       size_t len = n/sizeof(PetscScalar);
1944: #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1945:       PetscBLASInt   one = 1,blen;
1947:       PetscBLASIntCast(len,&blen);
1948:       PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one));
1949: #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1950:       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1951: #else
1952:       size_t      i;
1953:       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1954:       for (i=0; i<len; i++) y[i] = x[i];
1955: #endif
1956:     } else {
1957:       memcpy((char*)(a),(char*)(b),n);
1958:     }
1959: #else
1960:     memcpy((char*)(a),(char*)(b),n);
1961: #endif
1962:   }
1963:   return(0);
1964: }

1966: /*@C
1967:    PetscMemzero - Zeros the specified memory.

1969:    Not Collective

1971:    Input Parameters:
1972: +  a - pointer to beginning memory location
1973: -  n - length (in bytes) of memory to initialize

1975:    Level: intermediate

1977:    Compile Option:
1978:    PETSC_PREFER_BZERO - on certain machines (the IBM RS6000) the bzero() routine happens
1979:   to be faster than the memset() routine. This flag causes the bzero() routine to be used.

1981:    Not available from Fortran

1983:    Developer Note: this is inlined for fastest performance

1985:    Concepts: memory^zeroing
1986:    Concepts: zeroing^memory

1988: .seealso: PetscMemcpy()
1989: @*/
1990: PETSC_STATIC_INLINE PetscErrorCode  PetscMemzero(void *a,size_t n)
1991: {
1992:   if (n > 0) {
1993: #if defined(PETSC_USE_DEBUG)
1994:     if (!a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer");
1995: #endif
1996: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1997:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1998:       size_t      i,len = n/sizeof(PetscScalar);
1999:       PetscScalar *x = (PetscScalar*)a;
2000:       for (i=0; i<len; i++) x[i] = 0.0;
2001:     } else {
2002: #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
2003:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
2004:       PetscInt len = n/sizeof(PetscScalar);
2005:       fortranzero_(&len,(PetscScalar*)a);
2006:     } else {
2007: #endif
2008: #if defined(PETSC_PREFER_BZERO)
2009:       bzero((char *)a,n);
2010: #else
2011:       memset((char*)a,0,n);
2012: #endif
2013: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
2014:     }
2015: #endif
2016:   }
2017:   return 0;
2018: }

2020: /*MC
2021:    PetscPrefetchBlock - Prefetches a block of memory

2023:    Synopsis:
2024:     #include <petscsys.h>
2025:     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)

2027:    Not Collective

2029:    Input Parameters:
2030: +  a - pointer to first element to fetch (any type but usually PetscInt or PetscScalar)
2031: .  n - number of elements to fetch
2032: .  rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors)
2033: -  t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note

2035:    Level: developer

2037:    Notes:
2038:    The last two arguments (rw and t) must be compile-time constants.

2040:    Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality.  Not all architectures offer
2041:    equivalent locality hints, but the following macros are always defined to their closest analogue.
2042: +  PETSC_PREFETCH_HINT_NTA - Non-temporal.  Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched).
2043: .  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.
2044: .  PETSC_PREFETCH_HINT_T1 - Fetch to level 2 and higher (not L1).
2045: -  PETSC_PREFETCH_HINT_T2 - Fetch to high-level cache only.  (On many systems, T0 and T1 are equivalent.)

2047:    This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid
2048:    address).

2050:    Concepts: memory
2051: M*/
2052: #define PetscPrefetchBlock(a,n,rw,t) do {                               \
2053:     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
2054:     for ( ; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
2055:   } while (0)

2057: /*
2058:       Determine if some of the kernel computation routines use
2059:    Fortran (rather than C) for the numerical calculations. On some machines
2060:    and compilers (like complex numbers) the Fortran version of the routines
2061:    is faster than the C/C++ versions. The flag --with-fortran-kernels
2062:    should be used with ./configure to turn these on.
2063: */
2064: #if defined(PETSC_USE_FORTRAN_KERNELS)

2066: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
2067: #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
2068: #endif

2070: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
2071: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
2072: #endif

2074: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
2075: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
2076: #endif

2078: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
2079: #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
2080: #endif

2082: #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
2083: #define PETSC_USE_FORTRAN_KERNEL_NORM
2084: #endif

2086: #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
2087: #define PETSC_USE_FORTRAN_KERNEL_MAXPY
2088: #endif

2090: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
2091: #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
2092: #endif

2094: #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
2095: #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
2096: #endif

2098: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
2099: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
2100: #endif

2102: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
2103: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
2104: #endif

2106: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
2107: #define PETSC_USE_FORTRAN_KERNEL_MDOT
2108: #endif

2110: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
2111: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
2112: #endif

2114: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
2115: #define PETSC_USE_FORTRAN_KERNEL_AYPX
2116: #endif

2118: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
2119: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
2120: #endif

2122: #endif

2124: /*
2125:     Macros for indicating code that should be compiled with a C interface,
2126:    rather than a C++ interface. Any routines that are dynamically loaded
2127:    (such as the PCCreate_XXX() routines) must be wrapped so that the name
2128:    mangler does not change the functions symbol name. This just hides the
2129:    ugly extern "C" {} wrappers.
2130: */
2131: #if defined(__cplusplus)
2132: #define EXTERN_C_BEGIN extern "C" {
2133: #define EXTERN_C_END }
2134: #else
2135: #define EXTERN_C_BEGIN
2136: #define EXTERN_C_END
2137: #endif

2139: /* --------------------------------------------------------------------*/

2141: /*MC
2142:     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
2143:         communication

2145:    Level: beginner

2147:    Note: This manual page is a place-holder because MPICH does not have a manual page for MPI_Comm

2149: .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
2150: M*/

2152: /*MC
2153:    PetscScalar - PETSc type that represents either a double precision real number, a double precision
2154:        complex number, a single precision real number, a __float128 real or complex or a __fp16 real - if the code is configured
2155:        with --with-scalar-type=real,complex --with-precision=single,double,__float128,__fp16

2157:    Notes:
2158:    For MPI calls that require datatypes, use MPIU_SCALAR as the datatype for PetscScalar and MPIU_SUM, MPIU_MAX etc for operations. They will automatically work correctly regardless of the size of PetscScalar.

2160:    Level: beginner

2162: .seealso: PetscReal, PetscComplex, PetscInt, MPIU_REAL, MPIU_SCALAR, MPIU_COMPLEX, MPIU_INT
2163: M*/

2165: /*MC
2166:    PetscComplex - PETSc type that represents a complex number with precision matching that of PetscReal.

2168:    Synopsis:
2169:    #include <petscsys.h>
2170:    PetscComplex number = 1. + 2.*PETSC_i;

2172:    Notes:
2173:    For MPI calls that require datatypes, use MPIU_COMPLEX as the datatype for PetscComplex and MPIU_SUM etc for operations.
2174:           They will automatically work correctly regardless of the size of PetscComplex.

2176:           See PetscScalar for details on how to ./configure the size of PetscReal

2178:           Complex numbers are automatically available if PETSc was able to find a working complex implementation

2180:    Level: beginner

2182: .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_SCALAR, MPIU_COMPLEX, MPIU_INT, PETSC_i
2183: M*/

2185: /*MC
2186:    PetscReal - PETSc type that represents a real number version of PetscScalar


2189:    Notes:
2190:    For MPI calls that require datatypes, use MPIU_REAL as the datatype for PetscScalar and MPIU_SUM, MPIU_MAX, etc. for operations.
2191:           They will automatically work correctly regardless of the size of PetscReal.

2193:           See PetscScalar for details on how to ./configure the size of PetscReal.

2195:    Level: beginner

2197: .seealso: PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_SCALAR, MPIU_COMPLEX, MPIU_INT
2198: M*/

2200: /*MC
2201:    MPIU_SCALAR - MPI datatype corresponding to PetscScalar

2203:    Notes:
2204:    In MPI calls that require an MPI datatype that matches a PetscScalar or array of PetscScalar values, pass this value.

2206:    Level: beginner

2208: .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_COMPLEX, MPIU_INT
2209: M*/

2211: /*MC
2212:    MPIU_COMPLEX - MPI datatype corresponding to PetscComplex

2214:    Notes:
2215:    In MPI calls that require an MPI datatype that matches a PetscComplex or array of PetscComplex values, pass this value.

2217:    Level: beginner

2219: .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_SCALAR, MPIU_COMPLEX, MPIU_INT, PETSC_i
2220: M*/

2222: /*MC
2223:    MPIU_REAL - MPI datatype corresponding to PetscReal

2225:    Notes:
2226:    In MPI calls that require an MPI datatype that matches a PetscReal or array of PetscReal values, pass this value.

2228:    Level: beginner

2230: .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_SCALAR, MPIU_COMPLEX, MPIU_INT
2231: M*/

2233: #if defined(PETSC_HAVE_MPIIO)
2234: #if !defined(PETSC_WORDS_BIGENDIAN)
2235: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2236: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2237: #else
2238: #define MPIU_File_write_all(a,b,c,d,e) MPI_File_write_all(a,b,c,d,e)
2239: #define MPIU_File_read_all(a,b,c,d,e) MPI_File_read_all(a,b,c,d,e)
2240: #endif
2241: #endif

2243: /* the following petsc_static_inline require petscerror.h */

2245: /* Limit MPI to 32-bits */
2246: #define PETSC_MPI_INT_MAX  2147483647
2247: #define PETSC_MPI_INT_MIN -2147483647
2248: /* Limit BLAS to 32-bits */
2249: #define PETSC_BLAS_INT_MAX  2147483647
2250: #define PETSC_BLAS_INT_MIN -2147483647

2252: /*@C
2253:     PetscBLASIntCast - casts a PetscInt (which may be 64 bits in size) to a PetscBLASInt (which may be 32 bits in size), generates an
2254:          error if the PetscBLASInt is not large enough to hold the number.

2256:    Not Collective

2258:    Input Parameter:
2259: .     a - the PetscInt value

2261:    Output Parameter:
2262: .     b - the resulting PetscBLASInt value

2264:    Level: advanced

2266:    Not available from Fortran

2268: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast()
2269: @*/
2270: PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
2271: {
2273:   *b =  (PetscBLASInt)(a);
2274: #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
2275:   if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK");
2276: #endif
2277:   return(0);
2278: }

2280: /*@C
2281:     PetscMPIIntCast - casts a PetscInt (which may be 64 bits in size) to a PetscMPIInt (which may be 32 bits in size), generates an
2282:          error if the PetscMPIInt is not large enough to hold the number.

2284:    Not Collective

2286:    Input Parameter:
2287: .     a - the PetscInt value

2289:    Output Parameter:
2290: .     b - the resulting PetscMPIInt value

2292:    Level: advanced

2294:    Not available from Fortran

2296: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast()
2297: @*/
2298: PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
2299: {
2301:   *b =  (PetscMPIInt)(a);
2302: #if defined(PETSC_USE_64BIT_INDICES)
2303:   if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for MPI");
2304: #endif
2305:   return(0);
2306: }

2308: #define PetscInt64Mult(a,b)   ((PetscInt64)(a))*((PetscInt64)(b))

2310: /*@C

2312:    PetscRealIntMultTruncate - Computes the product of a positive PetscReal and a positive PetscInt and truncates the value to slightly less than the maximal possible value

2314:    Not Collective

2316:    Input Parameter:
2317: +     a - the PetscReal value
2318: -     b - the second value

2320:    Returns:
2321:       the result as a PetscInt value

2323:    Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2324:    Use PetscIntMultTruncate() to compute the product of two positive PetscInt and truncate to fit a PetscInt
2325:    Use PetscIntMultError() to compute the product of two PetscInt if you wish to generate an error if the result will not fit in a PetscInt

2327:    Developers Note:
2328:    We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks.

2330:    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.

2332:    Not available from Fortran

2334:    Level: advanced

2336: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2337: @*/
2338: PETSC_STATIC_INLINE PetscInt PetscRealIntMultTruncate(PetscReal a,PetscInt b)
2339: {
2340:   PetscInt64 r;

2342:   r  =  (PetscInt64) (a*(PetscReal)b);
2343:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2344:   return (PetscInt) r;
2345: }

2347: /*@C

2349:    PetscIntMultTruncate - Computes the product of two positive PetscInt and truncates the value to slightly less than the maximal possible value

2351:    Not Collective

2353:    Input Parameter:
2354: +     a - the PetscInt value
2355: -     b - the second value

2357:    Returns:
2358:       the result as a PetscInt value

2360:    Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2361:    Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2362:    Use PetscIntMultError() to compute the product of two PetscInt if you wish to generate an error if the result will not fit in a PetscInt

2364:    Not available from Fortran

2366:    Developers Note: We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks.

2368:    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.

2370:    Level: advanced

2372: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2373: @*/
2374: PETSC_STATIC_INLINE PetscInt PetscIntMultTruncate(PetscInt a,PetscInt b)
2375: {
2376:   PetscInt64 r;

2378:   r  =  PetscInt64Mult(a,b);
2379:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2380:   return (PetscInt) r;
2381: }

2383: /*@C

2385:    PetscIntSumTruncate - Computes the sum of two positive PetscInt and truncates the value to slightly less than the maximal possible value

2387:    Not Collective

2389:    Input Parameter:
2390: +     a - the PetscInt value
2391: -     b - the second value

2393:    Returns:
2394:      the result as a PetscInt value

2396:    Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2397:    Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2398:    Use PetscIntMultError() to compute the product of two PetscInt if you wish to generate an error if the result will not fit in a PetscInt

2400:    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.

2402:    Not available from Fortran

2404:    Level: advanced

2406: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2407: @*/
2408: PETSC_STATIC_INLINE PetscInt PetscIntSumTruncate(PetscInt a,PetscInt b)
2409: {
2410:   PetscInt64 r;

2412:   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2413:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2414:   return (PetscInt) r;
2415: }

2417: /*@C

2419:    PetscIntMultError - Computes the product of two positive PetscInt and generates an error with overflow.

2421:    Not Collective

2423:    Input Parameter:
2424: +     a - the PetscInt value
2425: -     b - the second value

2427:    Output Parameter:ma
2428: .     result - the result as a PetscInt value, or NULL if you do not want the result, you just want to check if it overflows

2430:    Use PetscInt64Mult() to compute the product of two 32 bit PetscInt and store in a PetscInt64
2431:    Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt

2433:    Not available from Fortran

2435:    Developers Note: We currently assume that PetscInt addition does not overflow, this is obviously wrong but requires many more checks.

2437:    Level: advanced

2439: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64(), PetscIntSumError()
2440: @*/
2441: PETSC_STATIC_INLINE PetscErrorCode PetscIntMultError(PetscInt a,PetscInt b,PetscInt *result)
2442: {
2443:   PetscInt64 r;

2446:   r  =  PetscInt64Mult(a,b);
2447: #if !defined(PETSC_USE_64BIT_INDICES)
2448:   if (r > PETSC_MAX_INT) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Product of two integer %d %d overflow, you must ./configure PETSc with --with-64-bit-indices for the case you are running",a,b);
2449: #endif
2450:   if (result) *result = (PetscInt) r;
2451:   return(0);
2452: }

2454: /*@C

2456:    PetscIntSumError - Computes the sum of two positive PetscInt and generates an error with overflow.

2458:    Not Collective

2460:    Input Parameter:
2461: +     a - the PetscInt value
2462: -     b - the second value

2464:    Output Parameter:ma
2465: .     c - the result as a PetscInt value,  or NULL if you do not want the result, you just want to check if it overflows

2467:    Use PetscInt64Mult() to compute the product of two 32 bit PetscInt and store in a PetscInt64
2468:    Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt

2470:    Not available from Fortran

2472:    Level: advanced

2474: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2475: @*/
2476: PETSC_STATIC_INLINE PetscErrorCode PetscIntSumError(PetscInt a,PetscInt b,PetscInt *result)
2477: {
2478:   PetscInt64 r;

2481:   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2482: #if !defined(PETSC_USE_64BIT_INDICES)
2483:   if (r > PETSC_MAX_INT) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Sum of two integer %d %d overflow, you must ./configure PETSc with --with-64-bit-indices for the case you are running",a,b);
2484: #endif
2485:   if (result) *result = (PetscInt) r;
2486:   return(0);
2487: }

2489: /*
2490:      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2491: */
2492: #if defined(hz)
2493: #undef hz
2494: #endif

2496: /*  For arrays that contain filenames or paths */


2499: #if defined(PETSC_HAVE_LIMITS_H)
2500: #include <limits.h>
2501: #endif
2502: #if defined(PETSC_HAVE_SYS_PARAM_H)
2503: #include <sys/param.h>
2504: #endif
2505: #if defined(PETSC_HAVE_SYS_TYPES_H)
2506: #include <sys/types.h>
2507: #endif
2508: #if defined(MAXPATHLEN)
2509: #  define PETSC_MAX_PATH_LEN     MAXPATHLEN
2510: #elif defined(MAX_PATH)
2511: #  define PETSC_MAX_PATH_LEN     MAX_PATH
2512: #elif defined(_MAX_PATH)
2513: #  define PETSC_MAX_PATH_LEN     _MAX_PATH
2514: #else
2515: #  define PETSC_MAX_PATH_LEN     4096
2516: #endif

2518: /*MC

2520:     PETSC_VERSION - This manual page provides information about how PETSc documents and uses its version information. This information is available to both C/C++
2521:                     and Fortran compilers when petscsys.h is included.


2524:     The current PETSc version and the API for accessing it are defined in petscversion.h

2526:     The complete version number is given as the triple  PETSC_VERSION_MAJOR.PETSC_VERSION_MINOR.PETSC_VERSION_SUBMINOR (in short hand x.y.z)

2528:     A change in the minor version number (y) indicates possible/likely changes in the PETSc API. Note this is different than with the semantic versioning convention
2529:     where only a change in the major version number (x) indicates a change in the API.

2531:     A subminor greater than zero indicates a patch release. Version x.y.z maintains source and binary compatibility with version x.y.w for all z and w

2533:     Use the macros PETSC_VERSION_EQ(x,y,z), PETSC_VERSION_LT(x,y,z), PETSC_VERSION_LE(x,y,z), PETSC_VERSION_GT(x,y,z),
2534:     PETSC_VERSION_GE(x,y,z) to determine if the current version is equal to, less than, less than or equal to, greater than or greater than or equal to a given
2535:     version number (x.y.z).

2537:     PETSC_RELEASE_DATE is the date the x.y version was released (i.e. the version before any patch releases)

2539:     PETSC_VERSION_DATE is the date the x.y.z version was released

2541:     PETSC_VERSION_GIT is the last git commit to the repository given in the form vx.y.z-wwwww

2543:     PETSC_VERSION_DATE_GIT is the date of the last git commit to the repository

2545:     Level: intermediate

2547:     PETSC_VERSION_() and PETSC_VERSION_PATCH are deprecated and will eventually be removed. For several releases PETSC_VERSION_PATCH is always 0

2549: M*/

2551: /*MC

2553:     UsingFortran - To use PETSc with Fortran you must use both PETSc include files and modules. At the beginning
2554:       of every function and module definition you need something like

2556: $
2557: $#include "petsc/finclude/petscXXX.h"
2558: $         use petscXXX

2560:      You can declare PETSc variables using either of the following.

2562: $    XXX variablename
2563: $    type(tXXX) variablename

2565:     For example,

2567: $#include "petsc/finclude/petscvec.h"
2568: $         use petscvec
2569: $
2570: $    Vec b
2571: $    type(tVec) x

2573:     Level: beginner

2575: M*/

2577: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2578: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2579: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2580: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2581: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2582: PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2583: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2584: PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt*,PetscInt*,PetscInt*,PetscInt*);

2586: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2587: PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt*,PetscInt[]);
2588: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2589: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2590: PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt*);
2591: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2592: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2593: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2594: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt[],PetscInt[],PetscInt[]);
2595: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2596: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2597: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2598: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2599: PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt,PetscInt[],void*,size_t,void*);
2600: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2601: PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt,PetscReal[],PetscInt[]);
2602: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2603: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2604: PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal,PetscInt,const PetscReal[],PetscReal,PetscInt*);
2605: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2606: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2607: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2608: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt[],const PetscInt[],PetscInt,const PetscInt[],const PetscInt[],PetscInt*,PetscInt**,PetscInt**);
2609: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt*,PetscInt**);
2610: PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt,const PetscMPIInt[],PetscInt,const PetscMPIInt[],PetscInt*,PetscMPIInt**);

2612: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2613: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);

2615: /*J
2616:     PetscRandomType - String with the name of a PETSc randomizer

2618:    Level: beginner

2620:    Notes:
2621:    To use SPRNG or RANDOM123 you must have ./configure PETSc
2622:    with the option --download-sprng or --download-random123

2624: .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2625: J*/
2626: typedef const char* PetscRandomType;
2627: #define PETSCRAND       "rand"
2628: #define PETSCRAND48     "rand48"
2629: #define PETSCSPRNG      "sprng"
2630: #define PETSCRANDER48   "rander48"
2631: #define PETSCRANDOM123  "random123"

2633: /* Logging support */
2634: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;

2636: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);

2638: /*S
2639:      PetscRandom - Abstract PETSc object that manages generating random numbers

2641:    Level: intermediate

2643:   Concepts: random numbers

2645: .seealso:  PetscRandomCreate(), PetscRandomGetValue(), PetscRandomType
2646: S*/
2647: typedef struct _p_PetscRandom*   PetscRandom;

2649: /* Dynamic creation and loading functions */
2650: PETSC_EXTERN PetscFunctionList PetscRandomList;

2652: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2653: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2654: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2655: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType*);
2656: PETSC_STATIC_INLINE PetscErrorCode PetscRandomViewFromOptions(PetscRandom A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
2657: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);

2659: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2660: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2661: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2662: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2663: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2664: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2665: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2666: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2667: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);

2669: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2670: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2671: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2672: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2673: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2674: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2675: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);
2676: PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2677: PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);

2679: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscDataType);
2680: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscDataType);
2681: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,void*,PetscInt,PetscDataType,PetscBool );
2682: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,void*,PetscInt,PetscDataType,PetscBool );
2683: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2684: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2685: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2686: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2687: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2688: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2689: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2690: #if defined(PETSC_USE_SOCKET_VIEWER)
2691: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);
2692: #endif

2694: /*
2695:    In binary files variables are stored using the following lengths,
2696:   regardless of how they are stored in memory on any one particular
2697:   machine. Use these rather then sizeof() in computing sizes for
2698:   PetscBinarySeek().
2699: */
2700: #define PETSC_BINARY_INT_SIZE   (32/8)
2701: #define PETSC_BINARY_FLOAT_SIZE  (32/8)
2702: #define PETSC_BINARY_CHAR_SIZE  (8/8)
2703: #define PETSC_BINARY_SHORT_SIZE  (16/8)
2704: #define PETSC_BINARY_DOUBLE_SIZE  (64/8)
2705: #define PETSC_BINARY_SCALAR_SIZE  sizeof(PetscScalar)

2707: /*E
2708:   PetscBinarySeekType - argument to PetscBinarySeek()

2710:   Level: advanced

2712: .seealso: PetscBinarySeek(), PetscBinarySynchronizedSeek()
2713: E*/
2714: typedef enum {PETSC_BINARY_SEEK_SET = 0,PETSC_BINARY_SEEK_CUR = 1,PETSC_BINARY_SEEK_END = 2} PetscBinarySeekType;
2715: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2716: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2717: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);

2719: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2720: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool );
2721: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2722: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2723: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2724: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);

2726: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2727: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2728: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2729: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2730: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2731: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt*,const void*,PetscMPIInt*,PetscMPIInt**,void*)
2732:   PetscAttrMPIPointerWithType(6,3);
2733: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2734:                                                     PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2735:                                                     PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2736:   PetscAttrMPIPointerWithType(6,3);
2737: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2738:                                                        MPI_Request**,MPI_Request**,
2739:                                                        PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2740:                                                        PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2741:   PetscAttrMPIPointerWithType(6,3);

2743: /*E
2744:     PetscBuildTwoSidedType - algorithm for setting up two-sided communication

2746: $  PETSC_BUILDTWOSIDED_ALLREDUCE - classical algorithm using an MPI_Allreduce with
2747: $      a buffer of length equal to the communicator size. Not memory-scalable due to
2748: $      the large reduction size. Requires only MPI-1.
2749: $  PETSC_BUILDTWOSIDED_IBARRIER - nonblocking algorithm based on MPI_Issend and MPI_Ibarrier.
2750: $      Proved communication-optimal in Hoefler, Siebert, and Lumsdaine (2010). Requires MPI-3.
2751: $  PETSC_BUILDTWOSIDED_REDSCATTER - similar to above, but use more optimized function
2752: $      that only communicates the part of the reduction that is necessary.  Requires MPI-2.

2754:    Level: developer

2756: .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedSetType(), PetscCommBuildTwoSidedGetType()
2757: E*/
2758: typedef enum {
2759:   PETSC_BUILDTWOSIDED_NOTSET = -1,
2760:   PETSC_BUILDTWOSIDED_ALLREDUCE = 0,
2761:   PETSC_BUILDTWOSIDED_IBARRIER = 1,
2762:   PETSC_BUILDTWOSIDED_REDSCATTER = 2
2763:   /* Updates here must be accompanied by updates in finclude/petscsys.h and the string array in mpits.c */
2764: } PetscBuildTwoSidedType;
2765: PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2766: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2767: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);

2769: PETSC_EXTERN PetscErrorCode PetscSSEIsEnabled(MPI_Comm,PetscBool  *,PetscBool  *);

2771: /*E
2772:   InsertMode - Whether entries are inserted or added into vectors or matrices

2774:   Level: beginner

2776: .seealso: VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2777:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(),
2778:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd()
2779: E*/
2780:  typedef enum {NOT_SET_VALUES, INSERT_VALUES, ADD_VALUES, MAX_VALUES, INSERT_ALL_VALUES, ADD_ALL_VALUES, INSERT_BC_VALUES, ADD_BC_VALUES} InsertMode;

2782: /*MC
2783:     INSERT_VALUES - Put a value into a vector or matrix, overwrites any previous value

2785:     Level: beginner

2787: .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2788:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), ADD_VALUES,
2789:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES

2791: M*/

2793: /*MC
2794:     ADD_VALUES - Adds a value into a vector or matrix, if there previously was no value, just puts the
2795:                 value into that location

2797:     Level: beginner

2799: .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2800:           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), INSERT_VALUES,
2801:           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES

2803: M*/

2805: /*MC
2806:     MAX_VALUES - Puts the maximum of the scattered/gathered value and the current value into each location

2808:     Level: beginner

2810: .seealso: InsertMode, VecScatterBegin(), VecScatterEnd(), ADD_VALUES, INSERT_VALUES

2812: M*/

2814: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);

2816: typedef enum {PETSC_SUBCOMM_GENERAL=0,PETSC_SUBCOMM_CONTIGUOUS=1,PETSC_SUBCOMM_INTERLACED=2} PetscSubcommType;
2817: PETSC_EXTERN const char *const PetscSubcommTypes[];

2819: /*S
2820:    PetscSubcomm - A decomposition of an MPI communicator into subcommunicators

2822:    Notes:
2823:    After a call to PetscSubcommSetType(), PetscSubcommSetTypeGeneral(), or PetscSubcommSetFromOptions() one may call
2824: $     PetscSubcommChild() returns the associated subcommunicator on this process
2825: $     PetscSubcommContiguousParent() returns a parent communitor but with all child of the same subcommunicator having contiguous rank

2827:    Sample Usage:
2828:        PetscSubcommCreate()
2829:        PetscSubcommSetNumber()
2830:        PetscSubcommSetType(PETSC_SUBCOMM_INTERLACED);
2831:        ccomm = PetscSubcommChild()
2832:        PetscSubcommDestroy()

2834:    Level: advanced

2836:    Concepts: communicator, create

2838:    Notes:
2839: $   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
2840: $   PETSC_SUBCOMM_CONTIGUOUS - each new communicator contains a set of process with contiguous ranks in the original MPI communicator
2841: $   PETSC_SUBCOMM_INTERLACED - each new communictor contains a set of processes equally far apart in rank from the others in that new communicator

2843:    Example: Consider a communicator with six processes split into 3 subcommunicators.
2844: $     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
2845: $     PETSC_SUBCOMM_INTERLACED - the first communicator contains rank 0,3, the second 1,4 and the third 2,5

2847:    Developer Notes:
2848:    This is used in objects such as PCREDUNDANT to manage the subcommunicators on which the redundant computations
2849:       are performed.


2852: .seealso: PetscSubcommCreate(), PetscSubcommSetNumber(), PetscSubcommSetType(), PetscSubcommView(), PetscSubcommSetFromOptions()

2854: S*/
2855: typedef struct _n_PetscSubcomm* PetscSubcomm;

2857: struct _n_PetscSubcomm {
2858:   MPI_Comm         parent;           /* parent communicator */
2859:   MPI_Comm         dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2860:   MPI_Comm         child;            /* the sub-communicator */
2861:   PetscMPIInt      n;                /* num of subcommunicators under the parent communicator */
2862:   PetscMPIInt      color;            /* color of processors belong to this communicator */
2863:   PetscMPIInt      *subsize;         /* size of subcommunicator[color] */
2864:   PetscSubcommType type;
2865:   char             *subcommprefix;
2866: };

2868: PETSC_STATIC_INLINE MPI_Comm PetscSubcommParent(PetscSubcomm scomm) {return scomm->parent;}
2869: PETSC_STATIC_INLINE MPI_Comm PetscSubcommChild(PetscSubcomm scomm) {return scomm->child;}
2870: PETSC_STATIC_INLINE MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm) {return scomm->dupparent;}
2871: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2872: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2873: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2874: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2875: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2876: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2877: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2878: PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm,const char[]);

2880: /*S
2881:      PetscHeap - A simple class for managing heaps

2883:    Level: intermediate

2885:   Concepts: random numbers

2887: .seealso:  PetscHeapCreate(), PetscHeapAdd(), PetscHeapPop(), PetscHeapPeek(), PetscHeapStash(), PetscHeapUnstash(), PetscHeapView(), PetscHeapDestroy()
2888: S*/
2889: typedef struct _PetscHeap *PetscHeap;

2891: PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt,PetscHeap*);
2892: PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap,PetscInt,PetscInt);
2893: PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap,PetscInt*,PetscInt*);
2894: PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap,PetscInt*,PetscInt*);
2895: PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap,PetscInt,PetscInt);
2896: PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2897: PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap*);
2898: PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap,PetscViewer);

2900: PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2901: typedef struct _n_PetscCommShared* PetscCommShared;
2902: PETSC_EXTERN PetscErrorCode PetscCommSharedGet(MPI_Comm,PetscCommShared*);
2903: PETSC_EXTERN PetscErrorCode PetscCommSharedGlobalToLocal(PetscCommShared,PetscMPIInt,PetscMPIInt*);
2904: PETSC_EXTERN PetscErrorCode PetscCommSharedGetComm(PetscCommShared,MPI_Comm*);

2906: /*S
2907:    PetscSegBuffer - a segmented extendable buffer

2909:    Level: developer

2911: .seealso: PetscSegBufferCreate(), PetscSegBufferGet(), PetscSegBufferExtract(), PetscSegBufferDestroy()
2912: S*/
2913: typedef struct _n_PetscSegBuffer *PetscSegBuffer;
2914: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2915: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2916: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2917: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2918: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2919: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2920: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2921: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);


2925:  * prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as
2926:  * possible. */
2927: PETSC_STATIC_INLINE PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg,PetscInt count,PetscInt *PETSC_RESTRICT *slot) {return PetscSegBufferGet(seg,(size_t)count,(void**)slot);}

2929: typedef struct _n_PetscOptionsHelpPrinted *PetscOptionsHelpPrinted;
2930: extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2931: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted*);
2932: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted*);
2933: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted,const char*,const char*,PetscBool*);

2935: PETSC_EXTERN PetscSegBuffer PetscCitationsList;

2937: /*@C
2938:       PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code.

2940:      Not Collective - only what is registered on rank 0 of PETSC_COMM_WORLD will be printed

2942:      Input Parameters:
2943: +      cite - the bibtex item, formated to displayed on multiple lines nicely
2944: -      set - a boolean variable initially set to PETSC_FALSE; this is used to insure only a single registration of the citation

2946:    Level: intermediate

2948:    Not available from Fortran

2950:      Options Database:
2951: .     -citations [filename]   - print out the bibtex entries for the given computation
2952: @*/
2953: PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2954: {
2955:   size_t         len;
2956:   char           *vstring;

2960:   if (set && *set) return(0);
2961:   PetscStrlen(cit,&len);
2962:   PetscSegBufferGet(PetscCitationsList,len,&vstring);
2963:   PetscMemcpy(vstring,cit,len);
2964:   if (set) *set = PETSC_TRUE;
2965:   return(0);
2966: }

2968: PETSC_EXTERN PetscErrorCode PetscURLShorten(const char[],char[],size_t);
2969: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm,char[],char[],size_t);
2970: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm,const char[],char[],size_t);
2971: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm,const char[],const char []);

2973: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm,char[],char[],size_t);
2974: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm,const char[],char[],char[],size_t);

2976: PETSC_EXTERN PetscErrorCode PetscGlobusGetTransfers(MPI_Comm,const char[],char[],size_t);

2978: PETSC_EXTERN PetscErrorCode PetscTextBelt(MPI_Comm,const char[],const char[],PetscBool*);
2979: PETSC_EXTERN PetscErrorCode PetscTellMyCell(MPI_Comm,const char[],const char[],PetscBool*);

2981: PETSC_EXTERN PetscErrorCode PetscPullJSONValue(const char[],const char[],char[],size_t,PetscBool*);
2982: PETSC_EXTERN PetscErrorCode PetscPushJSONValue(char[],const char[],const char[],size_t);


2985: #if defined(PETSC_USE_DEBUG)
2986: /*
2987:    Verify that all processes in the communicator have called this from the same line of code
2988:  */
2989: PETSC_EXTERN PetscErrorCode PetscAllreduceBarrierCheck(MPI_Comm,PetscMPIInt,int,const char*,const char *);

2991: /*MC
2992:    MPIU_Allreduce - a PETSc replacement for MPI_Allreduce() that tries to determine if the call from all the MPI processes occur from the
2993:                     same place in the PETSc code. This helps to detect bugs where different MPI processes follow different code paths
2994:                     resulting in inconsistent and incorrect calls to MPI_Allreduce().

2996:    Collective on MPI_Comm

2998:    Synopsis:
2999:      #include <petscsys.h>
3000:      PetscErrorCode MPIU_Allreduce(void *indata,void *outdata,PetscMPIInt count,MPI_Datatype datatype, MPI_Op op, MPI_Comm comm);

3002:    Input Parameters:
3003: +  indata - pointer to the input data to be reduced
3004: .  count - the number of MPI data items in indata and outdata
3005: .  datatype - the MPI datatype, for example MPI_INT
3006: .  op - the MPI operation, for example MPI_SUM
3007: -  comm - the MPI communicator on which the operation occurs

3009:    Output Parameter:
3010: .  outdata - the reduced values

3012:    Notes:
3013:    In optimized mode this directly calls MPI_Allreduce()

3015:    Level: developer

3017: .seealso: MPI_Allreduce()
3018: M*/
3019: #define MPIU_Allreduce(a,b,c,d,e,fcomm) (PetscAllreduceBarrierCheck(fcomm,c,__LINE__,PETSC_FUNCTION_NAME,__FILE__) || MPI_Allreduce(a,b,c,d,e,fcomm))
3020: #else
3021: #define MPIU_Allreduce(a,b,c,d,e,fcomm) MPI_Allreduce(a,b,c,d,e,fcomm)
3022: #endif

3024: #if defined(PETSC_HAVE_MPI_WIN_CREATE_FEATURE)
3025: PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint,PetscMPIInt,MPI_Info,MPI_Comm,void*,MPI_Win*);
3026: PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win,PetscMPIInt,MPI_Aint*,PetscMPIInt*,void*);
3027: #endif

3029: /*
3030:     Returned from PETSc functions that are called from MPI, such as related to attributes
3031: */
3032: PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CLASS,PETSC_MPI_ERROR_CODE;

3034: #endif