Actual source code: petscsys.h

petsc-3.13.6 2020-09-29
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: */
  5: #if !defined(PETSCSYS_H)
  6: #define PETSCSYS_H
  7: /* ========================================================================== */
  8: /*
  9:    petscconf.h is contained in ${PETSC_ARCH}/include/petscconf.h it is
 10:    found automatically by the compiler due to the -I${PETSC_DIR}/${PETSC_ARCH}/include.
 11:    For --prefix installs the ${PETSC_ARCH}/ does not exist and petscconf.h is in the same
 12:    directory as the other PETSc include files.
 13: */
 14: #include <petscconf.h>
 15: #include <petscfix.h>

 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:  #include <petscsystypes.h>

 38: /* ========================================================================== */
 39: /*
 40:    This facilitates using the C version of PETSc from C++ and the C++ version from C.
 41: */
 42: #if defined(__cplusplus)
 43: #  define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_CXX
 44: #else
 45: #  define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_C
 46: #endif

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

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

 64: #if defined(__cplusplus)
 65: #  define PETSC_INLINE PETSC_CXX_INLINE
 66: #else
 67: #  define PETSC_INLINE PETSC_C_INLINE
 68: #endif

 70: #define PETSC_STATIC_INLINE static PETSC_INLINE

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

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

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

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

114: /* ========================================================================== */

116: /*
117:     Defines the interface to MPI allowing the use of all MPI functions.

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

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

181: /*
182:     Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler
183:     see the top of mpicxx.h in the MPICH2 distribution.
184: */
185: #include <stdio.h>

187: /* MSMPI on 32bit windows requires this yukky hack - that breaks MPI standard compliance */
188: #if !defined(MPIAPI)
189: #define MPIAPI
190: #endif

192: /*
193:    Support for Clang (>=3.2) matching type tag arguments with void* buffer types.
194:    This allows the compiler to detect cases where the MPI datatype argument passed to a MPI routine
195:    does not match the actual type of the argument being passed in
196: */
197: #if defined(__has_attribute) && defined(works_with_const_which_is_not_true)
198: #  if __has_attribute(argument_with_type_tag) && __has_attribute(pointer_with_type_tag) && __has_attribute(type_tag_for_datatype)
199: #    define PetscAttrMPIPointerWithType(bufno,typeno) __attribute__((pointer_with_type_tag(MPI,bufno,typeno)))
200: #    define PetscAttrMPITypeTag(type)                 __attribute__((type_tag_for_datatype(MPI,type)))
201: #    define PetscAttrMPITypeTagLayoutCompatible(type) __attribute__((type_tag_for_datatype(MPI,type,layout_compatible)))
202: #  endif
203: #endif
204: #if !defined(PetscAttrMPIPointerWithType)
205: #  define PetscAttrMPIPointerWithType(bufno,typeno)
206: #  define PetscAttrMPITypeTag(type)
207: #  define PetscAttrMPITypeTagLayoutCompatible(type)
208: #endif

210: PETSC_EXTERN MPI_Datatype MPIU_ENUM PetscAttrMPITypeTag(PetscEnum);
211: PETSC_EXTERN MPI_Datatype MPIU_BOOL PetscAttrMPITypeTag(PetscBool);

213: /*MC
214:    MPIU_INT - MPI datatype corresponding to PetscInt

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

219:    Level: beginner

221: .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_SCALAR, MPIU_COMPLEX
222: M*/

224: #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 */
225: #  define MPIU_INT64 MPI_INT64_T
226: #  define PetscInt64_FMT PRId64
227: #elif (PETSC_SIZEOF_LONG_LONG == 8)
228: #  define MPIU_INT64 MPI_LONG_LONG_INT
229: #  define PetscInt64_FMT "lld"
230: #elif defined(PETSC_HAVE___INT64)
231: #  define MPIU_INT64 MPI_INT64_T
232: #  define PetscInt64_FMT "ld"
233: #else
234: #  error "cannot determine PetscInt64 type"
235: #endif

237: PETSC_EXTERN MPI_Datatype MPIU_FORTRANADDR;

239: #if defined(PETSC_USE_64BIT_INDICES)
240: #  define MPIU_INT MPIU_INT64
241: #  define PetscInt_FMT PetscInt64_FMT
242: #else
243: #  define MPIU_INT MPI_INT
244: #  define PetscInt_FMT "d"
245: #endif

247: /*
248:     For the rare cases when one needs to send a size_t object with MPI
249: */
250: PETSC_EXTERN MPI_Datatype MPIU_SIZE_T;

252: /*
253:       You can use PETSC_STDOUT as a replacement of stdout. You can also change
254:     the value of PETSC_STDOUT to redirect all standard output elsewhere
255: */
256: PETSC_EXTERN FILE* PETSC_STDOUT;

258: /*
259:       You can use PETSC_STDERR as a replacement of stderr. You can also change
260:     the value of PETSC_STDERR to redirect all standard error elsewhere
261: */
262: PETSC_EXTERN FILE* PETSC_STDERR;

264: /*MC
265:     PetscUnlikely - hints the compiler that the given condition is usually FALSE

267:     Synopsis:
268:  #include <petscsys.h>
269:     PetscBool  PetscUnlikely(PetscBool  cond)

271:     Not Collective

273:     Input Parameters:
274: .   cond - condition or expression

276:     Notes:
277:     This returns the same truth value, it is only a hint to compilers that the resulting
278:     branch is unlikely.

280:     Level: advanced

282: .seealso: PetscLikely(), CHKERRQ
283: M*/

285: /*MC
286:     PetscLikely - hints the compiler that the given condition is usually TRUE

288:     Synopsis:
289:  #include <petscsys.h>
290:     PetscBool  PetscLikely(PetscBool  cond)

292:     Not Collective

294:     Input Parameters:
295: .   cond - condition or expression

297:     Notes:
298:     This returns the same truth value, it is only a hint to compilers that the resulting
299:     branch is likely.

301:     Level: advanced

303: .seealso: PetscUnlikely()
304: M*/
305: #if defined(PETSC_HAVE_BUILTIN_EXPECT)
306: #  define PetscUnlikely(cond)   __builtin_expect(!!(cond),0)
307: #  define PetscLikely(cond)     __builtin_expect(!!(cond),1)
308: #else
309: #  define PetscUnlikely(cond)   (cond)
310: #  define PetscLikely(cond)     (cond)
311: #endif

313: /*
314:     Declare extern C stuff after including external header files
315: */

317: PETSC_EXTERN const char *const PetscBools[];

319: /*
320:     Defines elementary mathematics functions and constants.
321: */
322:  #include <petscmath.h>

324: PETSC_EXTERN const char *const PetscCopyModes[];

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

329:    Level: beginner

331:    Note:
332:    Accepted by many PETSc functions to not set a parameter and instead use
333:           some default

335:    Fortran Notes:
336:    This macro does not exist in Fortran; you must use PETSC_NULL_INTEGER,
337:           PETSC_NULL_DOUBLE_PRECISION etc

339: .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_DETERMINE

341: M*/
342: #define PETSC_IGNORE NULL

344: /* This is deprecated */
345: #define PETSC_NULL NULL

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

351:    Level: beginner

353: .seealso: PETSC_DEFAULT, PETSC_IGNORE, PETSC_DETERMINE

355: M*/
356: #define PETSC_DECIDE -1

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

362:    Level: beginner

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

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

370: M*/
371: #define PETSC_DETERMINE PETSC_DECIDE

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

377:    Level: beginner

379:    Fortran Notes:
380:    You need to use PETSC_DEFAULT_INTEGER or PETSC_DEFAULT_REAL.

382: .seealso: PETSC_DECIDE, PETSC_IGNORE, PETSC_DETERMINE

384: M*/
385: #define PETSC_DEFAULT -2

387: /*MC
388:     PETSC_COMM_WORLD - the equivalent of the MPI_COMM_WORLD communicator which represents
389:            all the processes that PETSc knows about.

391:    Level: beginner

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

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

402: .seealso: PETSC_COMM_SELF

404: M*/
405: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;

407: /*MC
408:     PETSC_COMM_SELF - This is always MPI_COMM_SELF

410:    Level: beginner

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

415: .seealso: PETSC_COMM_WORLD

417: M*/
418: #define PETSC_COMM_SELF MPI_COMM_SELF

420: PETSC_EXTERN PetscBool PetscBeganMPI;
421: PETSC_EXTERN PetscBool PetscInitializeCalled;
422: PETSC_EXTERN PetscBool PetscFinalizeCalled;
423: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;
424: PETSC_EXTERN PetscBool PetscCUDASynchronize;

426: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm));
427: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*);
428: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm*);

430: #if defined(PETSC_HAVE_CUDA)
431: PETSC_EXTERN PetscErrorCode PetscCUDAInitialize(MPI_Comm);
432: #endif

434: #if defined(PETSC_HAVE_ELEMENTAL)
435: PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
436: PETSC_EXTERN PetscErrorCode PetscElementalInitialized(PetscBool*);
437: PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
438: #endif

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

443:    Synopsis:
444:  #include <petscsys.h>
445:    PetscErrorCode PetscMalloc(size_t m,void **result)

447:    Not Collective

449:    Input Parameter:
450: .  m - number of bytes to allocate

452:    Output Parameter:
453: .  result - memory allocated

455:    Level: beginner

457:    Notes:
458:    Memory is always allocated at least double aligned

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

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

464: M*/
465: #define PetscMalloc(a,b)  ((*PetscTrMalloc)((a),PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(void**)(b)))

467: /*MC
468:    PetscRealloc - Rellocates memory

470:    Synopsis:
471:  #include <petscsys.h>
472:    PetscErrorCode PetscRealloc(size_t m,void **result)

474:    Not Collective

476:    Input Parameters:
477: +  m - number of bytes to allocate
478: -  result - previous memory

480:    Output Parameter:
481: .  result - new memory allocated

483:    Level: developer

485:    Notes:
486:    Memory is always allocated at least double aligned

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

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

493: /*MC
494:    PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment

496:    Synopsis:
497:  #include <petscsys.h>
498:    void *PetscAddrAlign(void *addr)

500:    Not Collective

502:    Input Parameters:
503: .  addr - address to align (any pointer type)

505:    Level: developer

507: .seealso: PetscMallocAlign()

509: M*/
510: #define PetscAddrAlign(a) (void*)((((PETSC_UINTPTR_T)(a))+(PETSC_MEMALIGN-1)) & ~(PETSC_MEMALIGN-1))

512: /*MC
513:    PetscMalloc1 - Allocates an array of memory aligned to PETSC_MEMALIGN

515:    Synopsis:
516:  #include <petscsys.h>
517:    PetscErrorCode PetscMalloc1(size_t m1,type **r1)

519:    Not Collective

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

524:    Output Parameter:
525: .  r1 - memory allocated

527:    Note:
528:    This uses the sizeof() of the memory type requested to determine the total memory to be allocated, therefore you should not
529:          multiply the number of elements requested by the sizeof() the type. For example use
530: $  PetscInt *id;
531: $  PetscMalloc1(10,&id);
532:         not
533: $  PetscInt *id;
534: $  PetscMalloc1(10*sizeof(PetscInt),&id);

536:         Does not zero the memory allocated, use PetscCalloc1() to obtain memory that has been zeroed.

538:    Level: beginner

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

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

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

548:    Synopsis:
549:  #include <petscsys.h>
550:    PetscErrorCode PetscCalloc1(size_t m1,type **r1)

552:    Not Collective

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

557:    Output Parameter:
558: .  r1 - memory allocated

560:    Notes:
561:    See PetsMalloc1() for more details on usage.

563:    Level: beginner

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

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

570: /*MC
571:    PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN

573:    Synopsis:
574:  #include <petscsys.h>
575:    PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)

577:    Not Collective

579:    Input Parameter:
580: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
581: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

583:    Output Parameter:
584: +  r1 - memory allocated in first chunk
585: -  r2 - memory allocated in second chunk

587:    Level: developer

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

591: M*/
592: #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))

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

597:    Synopsis:
598:  #include <petscsys.h>
599:    PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)

601:    Not Collective

603:    Input Parameter:
604: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
605: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

607:    Output Parameter:
608: +  r1 - memory allocated in first chunk
609: -  r2 - memory allocated in second chunk

611:    Level: developer

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

615: M*/
616: #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))

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

621:    Synopsis:
622:  #include <petscsys.h>
623:    PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

625:    Not Collective

627:    Input Parameter:
628: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
629: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
630: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

632:    Output Parameter:
633: +  r1 - memory allocated in first chunk
634: .  r2 - memory allocated in second chunk
635: -  r3 - memory allocated in third chunk

637:    Level: developer

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

641: M*/
642: #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))

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

647:    Synopsis:
648:  #include <petscsys.h>
649:    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

651:    Not Collective

653:    Input Parameter:
654: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
655: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
656: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

658:    Output Parameter:
659: +  r1 - memory allocated in first chunk
660: .  r2 - memory allocated in second chunk
661: -  r3 - memory allocated in third chunk

663:    Level: developer

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

667: M*/
668: #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))

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

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

677:    Not Collective

679:    Input Parameter:
680: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
681: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
682: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
683: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

685:    Output Parameter:
686: +  r1 - memory allocated in first chunk
687: .  r2 - memory allocated in second chunk
688: .  r3 - memory allocated in third chunk
689: -  r4 - memory allocated in fourth chunk

691:    Level: developer

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

695: M*/
696: #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))

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

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

705:    Not Collective

707:    Input Parameters:
708: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
709: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
710: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
711: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

713:    Output Parameters:
714: +  r1 - memory allocated in first chunk
715: .  r2 - memory allocated in second chunk
716: .  r3 - memory allocated in third chunk
717: -  r4 - memory allocated in fourth chunk

719:    Level: developer

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

723: M*/
724: #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))

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

729:    Synopsis:
730:  #include <petscsys.h>
731:    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)

733:    Not Collective

735:    Input Parameters:
736: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
737: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
738: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
739: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
740: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

742:    Output Parameters:
743: +  r1 - memory allocated in first chunk
744: .  r2 - memory allocated in second chunk
745: .  r3 - memory allocated in third chunk
746: .  r4 - memory allocated in fourth chunk
747: -  r5 - memory allocated in fifth chunk

749:    Level: developer

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

753: M*/
754: #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))

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

759:    Synopsis:
760:  #include <petscsys.h>
761:    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)

763:    Not Collective

765:    Input Parameters:
766: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
767: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
768: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
769: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
770: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

772:    Output Parameters:
773: +  r1 - memory allocated in first chunk
774: .  r2 - memory allocated in second chunk
775: .  r3 - memory allocated in third chunk
776: .  r4 - memory allocated in fourth chunk
777: -  r5 - memory allocated in fifth chunk

779:    Level: developer

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

783: M*/
784: #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))

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

789:    Synopsis:
790:  #include <petscsys.h>
791:    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)

793:    Not Collective

795:    Input Parameters:
796: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
797: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
798: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
799: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
800: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
801: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

803:    Output Parameteasr:
804: +  r1 - memory allocated in first chunk
805: .  r2 - memory allocated in second chunk
806: .  r3 - memory allocated in third chunk
807: .  r4 - memory allocated in fourth chunk
808: .  r5 - memory allocated in fifth chunk
809: -  r6 - memory allocated in sixth chunk

811:    Level: developer

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

815: M*/
816: #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))

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

821:    Synopsis:
822:  #include <petscsys.h>
823:    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)

825:    Not Collective

827:    Input Parameters:
828: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
829: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
830: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
831: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
832: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
833: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

835:    Output Parameters:
836: +  r1 - memory allocated in first chunk
837: .  r2 - memory allocated in second chunk
838: .  r3 - memory allocated in third chunk
839: .  r4 - memory allocated in fourth chunk
840: .  r5 - memory allocated in fifth chunk
841: -  r6 - memory allocated in sixth chunk

843:    Level: developer

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

847: M*/
848: #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))

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

853:    Synopsis:
854:  #include <petscsys.h>
855:    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)

857:    Not Collective

859:    Input Parameters:
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)
864: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
865: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
866: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

868:    Output Parameters:
869: +  r1 - memory allocated in first chunk
870: .  r2 - memory allocated in second chunk
871: .  r3 - memory allocated in third chunk
872: .  r4 - memory allocated in fourth chunk
873: .  r5 - memory allocated in fifth chunk
874: .  r6 - memory allocated in sixth chunk
875: -  r7 - memory allocated in seventh chunk

877:    Level: developer

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

881: M*/
882: #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))

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

887:    Synopsis:
888:  #include <petscsys.h>
889:    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)

891:    Not Collective

893:    Input Parameters:
894: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
895: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
896: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
897: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
898: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
899: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
900: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

902:    Output Parameters:
903: +  r1 - memory allocated in first chunk
904: .  r2 - memory allocated in second chunk
905: .  r3 - memory allocated in third chunk
906: .  r4 - memory allocated in fourth chunk
907: .  r5 - memory allocated in fifth chunk
908: .  r6 - memory allocated in sixth chunk
909: -  r7 - memory allocated in seventh chunk

911:    Level: developer

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

915: M*/
916: #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))

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

921:    Synopsis:
922:  #include <petscsys.h>
923:    PetscErrorCode PetscNew(type **result)

925:    Not Collective

927:    Output Parameter:
928: .  result - memory allocated, sized to match pointer type

930:    Level: beginner

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

934: M*/
935: #define PetscNew(b)      PetscCalloc1(1,(b))

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

941:    Synopsis:
942:  #include <petscsys.h>
943:    PetscErrorCode PetscNewLog(PetscObject obj,type **result)

945:    Not Collective

947:    Input Parameter:
948: .  obj - object memory is logged to

950:    Output Parameter:
951: .  result - memory allocated, sized to match pointer type

953:    Level: developer

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

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

960: /*MC
961:    PetscFree - Frees memory

963:    Synopsis:
964:  #include <petscsys.h>
965:    PetscErrorCode PetscFree(void *memory)

967:    Not Collective

969:    Input Parameter:
970: .   memory - memory to free (the pointer is ALWAYS set to NULL upon success)

972:    Level: beginner

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

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

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

981: M*/
982: #define PetscFree(a)   ((*PetscTrFree)((void*)(a),__LINE__,PETSC_FUNCTION_NAME,__FILE__) || ((a) = NULL,0))

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

987:    Synopsis:
988:  #include <petscsys.h>
989:    PetscErrorCode PetscFree2(void *memory1,void *memory2)

991:    Not Collective

993:    Input Parameters:
994: +   memory1 - memory to free
995: -   memory2 - 2nd memory to free

997:    Level: developer

999:    Notes:
1000:     Memory must have been obtained with PetscMalloc2()

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

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

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

1010:    Synopsis:
1011:  #include <petscsys.h>
1012:    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)

1014:    Not Collective

1016:    Input Parameters:
1017: +   memory1 - memory to free
1018: .   memory2 - 2nd memory to free
1019: -   memory3 - 3rd memory to free

1021:    Level: developer

1023:    Notes:
1024:     Memory must have been obtained with PetscMalloc3()

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

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

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

1034:    Synopsis:
1035:  #include <petscsys.h>
1036:    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)

1038:    Not Collective

1040:    Input Parameters:
1041: +   m1 - memory to free
1042: .   m2 - 2nd memory to free
1043: .   m3 - 3rd memory to free
1044: -   m4 - 4th memory to free

1046:    Level: developer

1048:    Notes:
1049:     Memory must have been obtained with PetscMalloc4()

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

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

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

1059:    Synopsis:
1060:  #include <petscsys.h>
1061:    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)

1063:    Not Collective

1065:    Input Parameters:
1066: +   m1 - memory to free
1067: .   m2 - 2nd memory to free
1068: .   m3 - 3rd memory to free
1069: .   m4 - 4th memory to free
1070: -   m5 - 5th memory to free

1072:    Level: developer

1074:    Notes:
1075:     Memory must have been obtained with PetscMalloc5()

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

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

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

1085:    Synopsis:
1086:  #include <petscsys.h>
1087:    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)

1089:    Not Collective

1091:    Input Parameters:
1092: +   m1 - memory to free
1093: .   m2 - 2nd memory to free
1094: .   m3 - 3rd memory to free
1095: .   m4 - 4th memory to free
1096: .   m5 - 5th memory to free
1097: -   m6 - 6th memory to free


1100:    Level: developer

1102:    Notes:
1103:     Memory must have been obtained with PetscMalloc6()

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

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

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

1113:    Synopsis:
1114:  #include <petscsys.h>
1115:    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)

1117:    Not Collective

1119:    Input Parameters:
1120: +   m1 - memory to free
1121: .   m2 - 2nd memory to free
1122: .   m3 - 3rd memory to free
1123: .   m4 - 4th memory to free
1124: .   m5 - 5th memory to free
1125: .   m6 - 6th memory to free
1126: -   m7 - 7th memory to free


1129:    Level: developer

1131:    Notes:
1132:     Memory must have been obtained with PetscMalloc7()

1134: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1135:           PetscMalloc7()

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

1140: PETSC_EXTERN PetscErrorCode PetscMallocA(int,PetscBool,int,const char *,const char *,size_t,void *,...);
1141: PETSC_EXTERN PetscErrorCode PetscFreeA(int,int,const char *,const char *,void *,...);
1142: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,PetscBool,int,const char[],const char[],void**);
1143: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1144: PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t,int,const char[],const char[],void**);
1145: PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1146: PETSC_EXTERN PetscErrorCode PetscMallocSet(PetscErrorCode (*)(size_t,PetscBool,int,const char[],const char[],void**),PetscErrorCode (*)(void*,int,const char[],const char[]),PetscErrorCode (*)(size_t,int,const char[],const char[], void **));
1147: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);

1149: /*
1150:   Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1151: */
1152: PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1153: PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);
1154: #if defined(PETSC_HAVE_CUDA)
1155: PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void);
1156: PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void);
1157: #endif

1159: #define MPIU_PETSCLOGDOUBLE  MPI_DOUBLE
1160: #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION

1162: /*
1163:    Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1164: */
1165: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1166: PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1167: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1168: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1169: PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1170: PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int,PetscLogDouble*);
1171: PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool,PetscBool);
1172: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*,PetscBool*,PetscBool*);
1173: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1174: PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1175: PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool*);

1177: PETSC_EXTERN const char *const PetscDataTypes[];
1178: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1179: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1180: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1181: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);

1183: /*
1184:     Basic memory and string operations. These are usually simple wrappers
1185:    around the basic Unix system calls, but a few of them have additional
1186:    functionality and/or error checking.
1187: */
1188: PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1189: PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1190: PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***);
1191: PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1192: PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool  *);
1193: PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool  *);
1194: PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1195: PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1196: PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1197: PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1198: PETSC_EXTERN PetscErrorCode PetscStrlcat(char[],const char[],size_t);
1199: PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1200: PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1201: PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1202: PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1203: PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1204: PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1205: PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1206: PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1207: PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*);
1208: PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1209: PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1210: PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1211: PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1212: PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt,const char *const*,char***);
1213: PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt,char***);
1214: PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);

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

1218: PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1219: PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1220: PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);

1222: PETSC_EXTERN PetscErrorCode PetscStrInList(const char[],const char[],char,PetscBool*);
1223: PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1224: PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);

1226: /*
1227:    These are MPI operations for MPI_Allreduce() etc
1228: */
1229: PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1230: #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1231: PETSC_EXTERN MPI_Op MPIU_SUM;
1232: #else
1233: #define MPIU_SUM MPI_SUM
1234: #endif
1235: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1236: PETSC_EXTERN MPI_Op MPIU_MAX;
1237: PETSC_EXTERN MPI_Op MPIU_MIN;
1238: #else
1239: #define MPIU_MAX MPI_MAX
1240: #define MPIU_MIN MPI_MIN
1241: #endif
1242: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);

1244: PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1245: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);

1247: PETSC_EXTERN const char *const PetscFileModes[];

1249: /*
1250:     Defines PETSc error handling.
1251: */
1252:  #include <petscerror.h>

1254: #define PETSC_SMALLEST_CLASSID  1211211
1255: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1256: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1257: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);
1258: PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject,PetscObjectId*);
1259: PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject,PetscObjectId,PetscBool*);


1262: /*
1263:    Routines that get memory usage information from the OS
1264: */
1265: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1266: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1267: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1268: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);

1270: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);

1272: /*
1273:    Initialization of PETSc
1274: */
1275: PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1276: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1277: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1278: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1279: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1280: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1281: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1282: PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1283: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1284: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);

1286: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1287: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);

1289: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1290: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1291: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1292: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);

1294: PETSC_EXTERN PetscErrorCode PetscMonitorCompare(PetscErrorCode (*)(void),void *,PetscErrorCode (*)(void**),PetscErrorCode (*)(void),void *,PetscErrorCode (*)(void**),PetscBool *);


1297: /*
1298:      These are so that in extern C code we can caste function pointers to non-extern C
1299:    function pointers. Since the regular C++ code expects its function pointers to be C++
1300: */
1301: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1302: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1303: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);

1305: /*
1306:     Functions that can act on any PETSc object.
1307: */
1308: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1309: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1310: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1311: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1312: PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1313: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1314: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1315: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1316: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1317: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1318: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1319: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1320: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1321: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1322: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt *);
1323: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1324: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1325: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject *);
1326: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void));
1327: #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d))
1328: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1329: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1330: PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1331: PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject,PetscObject);
1332: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt *);

1334:  #include <petscviewertypes.h>
1335:  #include <petscoptions.h>

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

1339: PETSC_EXTERN PetscErrorCode PetscMemoryShowUsage(PetscViewer,const char[]);
1340: PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer,const char[]);
1341: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1342: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1343: #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1344: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1345: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1346: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1347: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1348: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1349: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1350: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1351: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1352: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,PetscObject,const char[]);
1353: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1354: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1355: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject,const char[],PetscBool *);
1356: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1357: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1358: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1359: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);

1361: #if defined(PETSC_HAVE_SAWS)
1362: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1363: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1364: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1365: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1366: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1367: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1368: PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1369: PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1370: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1371: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);

1373: #else
1374: #define PetscSAWsBlock()                        0
1375: #define PetscObjectSAWsViewOff(obj)             0
1376: #define PetscObjectSAWsSetBlock(obj,flg)        0
1377: #define PetscObjectSAWsBlock(obj)               0
1378: #define PetscObjectSAWsGrantAccess(obj)         0
1379: #define PetscObjectSAWsTakeAccess(obj)          0
1380: #define PetscStackViewSAWs()                    0
1381: #define PetscStackSAWsViewOff()                 0
1382: #define PetscStackSAWsTakeAccess()
1383: #define PetscStackSAWsGrantAccess()

1385: #endif

1387: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1388: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1389: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);

1391: #if defined(PETSC_USE_DEBUG)
1392: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1393: #endif
1394: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);

1396: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1397: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1398: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1399: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1400: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1401: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);

1403: /*
1404:     Dynamic library lists. Lists of names of routines in objects or in dynamic
1405:   link libraries that will be loaded as needed.
1406: */

1408: #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1409: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1410: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1411: #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1412: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1413: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[],const char[]);
1414: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1415: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1416: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);

1418: PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1419: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1420: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1421: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1422: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1423: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1424: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1425: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);

1427: /*
1428:      Useful utility routines
1429: */
1430: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1431: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1432: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1433: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1434: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1435: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);
1436: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm,PetscInt[],PetscInt[]);
1437: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm,PetscReal[],PetscReal[]);

1439: /*MC
1440:     PetscNot - negates a logical type value and returns result as a PetscBool

1442:     Notes:
1443:     This is useful in cases like
1444: $     int        *a;
1445: $     PetscBool  flag = PetscNot(a)
1446:      where !a would not return a PetscBool because we cannot provide a cast from int to PetscBool in C.

1448:     Level: beginner

1450:     .seealso : PetscBool, PETSC_TRUE, PETSC_FALSE
1451: M*/
1452: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)

1454: /*MC
1455:     PetscHelpPrintf - Prints help messages.

1457:    Synopsis:
1458:  #include <petscsys.h>
1459:      PetscErrorCode (*PetscHelpPrintf)(const char format[],...);

1461:     Not Collective

1463:     Input Parameters:
1464: .   format - the usual printf() format string

1466:    Level: developer

1468:     Fortran Note:
1469:     This routine is not supported in Fortran.


1472: .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1473: M*/
1474: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);

1476: /*
1477:      Defines PETSc profiling.
1478: */
1479:  #include <petsclog.h>

1481: /*
1482:       Simple PETSc parallel IO for ASCII printing
1483: */
1484: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1485: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1486: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1487: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1488: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1489: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1490: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);
1491: PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[],size_t,const char*,PetscInt,const PetscReal[]);

1493: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1494: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1495: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);

1497: PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char*,size_t*);
1498: PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char*,char *);

1500: #if defined(PETSC_HAVE_POPEN)
1501: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1502: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*);
1503: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1504: #endif

1506: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1507: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1508: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*);
1509: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1510: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1511: PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1512: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);

1514: PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1515: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void **);
1516: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void *);
1517: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1518: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer *);
1519: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void*));
1520: PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void*);

1522: /*
1523:    For use in debuggers
1524: */
1525: PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1526: PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1527: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1528: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1529: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);

1531: #include <stddef.h>
1532: #include <string.h>             /* for memcpy, memset */
1533: #include <stdlib.h>

1535: #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1536: #include <xmmintrin.h>
1537: #endif

1539: /*@C
1540:    PetscMemmove - Copies n bytes, beginning at location b, to the space
1541:    beginning at location a. Copying  between regions that overlap will
1542:    take place correctly. Use PetscMemcpy() if the locations do not overlap

1544:    Not Collective

1546:    Input Parameters:
1547: +  b - pointer to initial memory space
1548: .  a - pointer to copy space
1549: -  n - length (in bytes) of space to copy

1551:    Level: intermediate

1553:    Note:
1554:    PetscArraymove() is preferred
1555:    This routine is analogous to memmove().

1557:    Developers Note: This is inlined for performance

1559: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscStrallocpy(),
1560:           PetscArraymove()
1561: @*/
1562: PETSC_STATIC_INLINE PetscErrorCode PetscMemmove(void *a,const void *b,size_t n)
1563: {
1565:   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to null pointer");
1566:   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1567: #if !defined(PETSC_HAVE_MEMMOVE)
1568:   if (a < b) {
1569:     if (a <= b - n) memcpy(a,b,n);
1570:     else {
1571:       memcpy(a,b,(int)(b - a));
1572:       PetscMemmove(b,b + (int)(b - a),n - (int)(b - a));
1573:     }
1574:   } else {
1575:     if (b <= a - n) memcpy(a,b,n);
1576:     else {
1577:       memcpy(b + n,b + (n - (int)(a - b)),(int)(a - b));
1578:       PetscMemmove(a,b,n - (int)(a - b));
1579:     }
1580:   }
1581: #else
1582:   memmove((char*)(a),(char*)(b),n);
1583: #endif
1584:   return(0);
1585: }

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

1592:    Not Collective

1594:    Input Parameters:
1595: +  b - pointer to initial memory space
1596: -  n - length (in bytes) of space to copy

1598:    Output Parameter:
1599: .  a - pointer to copy space

1601:    Level: intermediate

1603:    Compile Option:
1604:     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1605:                                   for memory copies on double precision values.
1606:     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1607:                                   for memory copies on double precision values.
1608:     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1609:                                   for memory copies on double precision values.

1611:    Note:
1612:    Prefer PetscArraycpy()
1613:    This routine is analogous to memcpy().
1614:    Not available from Fortran

1616:    Developer Note: this is inlined for fastest performance

1618: .seealso: PetscMemzero(), PetscMemcmp(), PetscArrayzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy()

1620: @*/
1621: PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n)
1622: {
1623: #if defined(PETSC_USE_DEBUG)
1624:   size_t al = (size_t) a,bl = (size_t) b;
1625:   size_t nl = (size_t) n;
1627:   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1628:   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1629: #else
1631: #endif
1632:   if (a != b && n > 0) {
1633: #if defined(PETSC_USE_DEBUG)
1634:     if ((al > bl && (al - bl) < nl) || (bl - al) < nl)  SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1635:               or make sure your copy regions and lengths are correct. \n\
1636:               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1637: #endif
1638: #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1639:    if (!(a % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1640:       size_t len = n/sizeof(PetscScalar);
1641: #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1642:       PetscBLASInt   one = 1,blen;
1644:       PetscBLASIntCast(len,&blen);
1645:       PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one));
1646: #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1647:       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1648: #else
1649:       size_t      i;
1650:       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1651:       for (i=0; i<len; i++) y[i] = x[i];
1652: #endif
1653:     } else {
1654:       memcpy((char*)(a),(char*)(b),n);
1655:     }
1656: #else
1657:     memcpy((char*)(a),(char*)(b),n);
1658: #endif
1659:   }
1660:   return(0);
1661: }

1663: /*@C
1664:    PetscMemzero - Zeros the specified memory.

1666:    Not Collective

1668:    Input Parameters:
1669: +  a - pointer to beginning memory location
1670: -  n - length (in bytes) of memory to initialize

1672:    Level: intermediate

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

1678:    Not available from Fortran
1679:    Prefer PetscArrayzero()

1681:    Developer Note: this is inlined for fastest performance

1683: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy()
1684: @*/
1685: PETSC_STATIC_INLINE PetscErrorCode  PetscMemzero(void *a,size_t n)
1686: {
1687:   if (n > 0) {
1688: #if defined(PETSC_USE_DEBUG)
1689:     if (!a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer");
1690: #endif
1691: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1692:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1693:       size_t      i,len = n/sizeof(PetscScalar);
1694:       PetscScalar *x = (PetscScalar*)a;
1695:       for (i=0; i<len; i++) x[i] = 0.0;
1696:     } else {
1697: #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1698:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1699:       PetscInt len = n/sizeof(PetscScalar);
1700:       fortranzero_(&len,(PetscScalar*)a);
1701:     } else {
1702: #endif
1703: #if defined(PETSC_PREFER_BZERO)
1704:       bzero((char *)a,n);
1705: #else
1706:       memset((char*)a,0,n);
1707: #endif
1708: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1709:     }
1710: #endif
1711:   }
1712:   return 0;
1713: }

1715: /*MC
1716:    PetscArraycmp - Compares two arrays in memory.

1718:    Synopsis:
1719:  #include <petscsys.h>
1720:     PetscErrorCode PetscArraycmp(const anytype *str1,const anytype *str2,size_t cnt,PetscBool *e)

1722:    Not Collective

1724:    Input Parameters:
1725: +  str1 - First array
1726: .  str2 - Second array
1727: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1729:    Output Parameters:
1730: .   e - PETSC_TRUE if equal else PETSC_FALSE.

1732:    Level: intermediate

1734:    Note:
1735:    This routine is a preferred replacement to PetscMemcmp()
1736:    The arrays must be of the same type

1738: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy(),
1739:           PetscArraymove()
1740: M*/
1741: #define  PetscArraycmp(str1,str2,cnt,e) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcmp(str1,str2,(size_t)(cnt)*sizeof(*(str1)),e))

1743: /*MC
1744:    PetscArraymove - Copies from one array in memory to another, the arrays may overlap. Use PetscArraycpy() when the arrays
1745:                     do not overlap

1747:    Synopsis:
1748:  #include <petscsys.h>
1749:     PetscErrorCode PetscArraymove(anytype *str1,const anytype *str2,size_t cnt)

1751:    Not Collective

1753:    Input Parameters:
1754: +  str1 - First array
1755: .  str2 - Second array
1756: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1758:    Level: intermediate

1760:    Note:
1761:    This routine is a preferred replacement to PetscMemmove()
1762:    The arrays must be of the same type

1764: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycpy(), PetscMemmove(), PetscArraycmp(), PetscStrallocpy()
1765: M*/
1766: #define  PetscArraymove(str1,str2,cnt) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemmove(str1,str2,(size_t)(cnt)*sizeof(*(str1))))

1768: /*MC
1769:    PetscArraycpy - Copies from one array in memory to another

1771:    Synopsis:
1772:  #include <petscsys.h>
1773:     PetscErrorCode PetscArraycpy(anytype *str1,const anytype *str2,size_t cnt)

1775:    Not Collective

1777:    Input Parameters:
1778: +  str1 - First array (destination)
1779: .  str2 - Second array (source)
1780: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1782:    Level: intermediate

1784:    Note:
1785:    This routine is a preferred replacement to PetscMemcpy()
1786:    The arrays must be of the same type

1788: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraymove(), PetscMemmove(), PetscArraycmp(), PetscStrallocpy()
1789: M*/
1790: #define  PetscArraycpy(str1,str2,cnt) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcpy(str1,str2,(size_t)(cnt)*sizeof(*(str1))))

1792: /*MC
1793:    PetscArrayzero - Zeros an array in memory.

1795:    Synopsis:
1796:  #include <petscsys.h>
1797:     PetscErrorCode PetscArrayzero(anytype *str1,size_t cnt)

1799:    Not Collective

1801:    Input Parameters:
1802: +  str1 - array
1803: -  cnt  - Count of the array, not in bytes, but number of entries in the array

1805:    Level: intermediate

1807:    Note:
1808:    This routine is a preferred replacement to PetscMemzero()

1810: .seealso: PetscMemcpy(), PetscMemcmp(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy(), PetscArraymove()
1811: M*/
1812: #define  PetscArrayzero(str1,cnt) PetscMemzero(str1,(size_t)(cnt)*sizeof(*(str1)))

1814: /*MC
1815:    PetscPrefetchBlock - Prefetches a block of memory

1817:    Synopsis:
1818:  #include <petscsys.h>
1819:     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)

1821:    Not Collective

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

1829:    Level: developer

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

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

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

1844: M*/
1845: #define PetscPrefetchBlock(a,n,rw,t) do {                               \
1846:     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
1847:     for ( ; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
1848:   } while (0)

1850: /*
1851:       Determine if some of the kernel computation routines use
1852:    Fortran (rather than C) for the numerical calculations. On some machines
1853:    and compilers (like complex numbers) the Fortran version of the routines
1854:    is faster than the C/C++ versions. The flag --with-fortran-kernels
1855:    should be used with ./configure to turn these on.
1856: */
1857: #if defined(PETSC_USE_FORTRAN_KERNELS)

1859: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1860: #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1861: #endif

1863: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
1864: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
1865: #endif

1867: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1868: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
1869: #endif

1871: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1872: #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
1873: #endif

1875: #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
1876: #define PETSC_USE_FORTRAN_KERNEL_NORM
1877: #endif

1879: #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
1880: #define PETSC_USE_FORTRAN_KERNEL_MAXPY
1881: #endif

1883: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1884: #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
1885: #endif

1887: #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1888: #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
1889: #endif

1891: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
1892: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
1893: #endif

1895: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1896: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
1897: #endif

1899: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
1900: #define PETSC_USE_FORTRAN_KERNEL_MDOT
1901: #endif

1903: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
1904: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
1905: #endif

1907: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
1908: #define PETSC_USE_FORTRAN_KERNEL_AYPX
1909: #endif

1911: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
1912: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
1913: #endif

1915: #endif

1917: /*
1918:     Macros for indicating code that should be compiled with a C interface,
1919:    rather than a C++ interface. Any routines that are dynamically loaded
1920:    (such as the PCCreate_XXX() routines) must be wrapped so that the name
1921:    mangler does not change the functions symbol name. This just hides the
1922:    ugly extern "C" {} wrappers.
1923: */
1924: #if defined(__cplusplus)
1925: #  define EXTERN_C_BEGIN extern "C" {
1926: #  define EXTERN_C_END }
1927: #else
1928: #  define EXTERN_C_BEGIN
1929: #  define EXTERN_C_END
1930: #endif

1932: /* --------------------------------------------------------------------*/

1934: /*MC
1935:     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
1936:         communication

1938:    Level: beginner

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

1942: .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
1943: M*/

1945: #if defined(PETSC_HAVE_MPIIO)
1946: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1947: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1948: PETSC_EXTERN PetscErrorCode MPIU_File_write_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1949: PETSC_EXTERN PetscErrorCode MPIU_File_read_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1950: PETSC_EXTERN PetscErrorCode MPIU_File_write_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1951: PETSC_EXTERN PetscErrorCode MPIU_File_read_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1952: #endif

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

1956: /* Limit MPI to 32-bits */
1957: #define PETSC_MPI_INT_MAX  2147483647
1958: #define PETSC_MPI_INT_MIN -2147483647
1959: /* Limit BLAS to 32-bits */
1960: #define PETSC_BLAS_INT_MAX  2147483647
1961: #define PETSC_BLAS_INT_MIN -2147483647

1963: /*@C
1964:     PetscIntCast - casts a PetscInt64 (which is 64 bits in size) to a PetscInt (which may be 32 bits in size), generates an
1965:          error if the PetscInt is not large enough to hold the number.

1967:    Not Collective

1969:    Input Parameter:
1970: .     a - the PetscInt64 value

1972:    Output Parameter:
1973: .     b - the resulting PetscInt value

1975:    Level: advanced

1977:    Notes: If integers needed for the applications are too large to fit in 32 bit ints you can ./configure using --with-64-bit-indices to make PetscInt use 64 bit ints

1979:    Not available from Fortran

1981: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscBLASIntCast()
1982: @*/
1983: PETSC_STATIC_INLINE PetscErrorCode PetscIntCast(PetscInt64 a,PetscInt *b)
1984: {
1986: #if !defined(PETSC_USE_64BIT_INDICES)
1987:   if ((a) > PETSC_MAX_INT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Integer too long for PetscInt, you may need to ./configure using --with-64-bit-indices");
1988: #endif
1989:   *b =  (PetscInt)(a);
1990:   return(0);
1991: }

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

1997:    Not Collective

1999:    Input Parameter:
2000: .     a - the PetscInt value

2002:    Output Parameter:
2003: .     b - the resulting PetscBLASInt value

2005:    Level: advanced

2007:    Not available from Fortran

2009: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscIntCast()
2010: @*/
2011: PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
2012: {
2014: #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
2015:   *b = 0;
2016:   if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK");
2017: #endif
2018:   *b =  (PetscBLASInt)(a);
2019:   return(0);
2020: }

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

2026:    Not Collective

2028:    Input Parameter:
2029: .     a - the PetscInt value

2031:    Output Parameter:
2032: .     b - the resulting PetscMPIInt value

2034:    Level: advanced

2036:    Not available from Fortran

2038: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntCast()
2039: @*/
2040: PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
2041: {
2043: #if defined(PETSC_USE_64BIT_INDICES)
2044:   *b = 0;
2045:   if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for MPI");
2046: #endif
2047:   *b =  (PetscMPIInt)(a);
2048:   return(0);
2049: }

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

2053: /*@C

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

2057:    Not Collective

2059:    Input Parameter:
2060: +     a - the PetscReal value
2061: -     b - the second value

2063:    Returns:
2064:       the result as a PetscInt value

2066:    Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2067:    Use PetscIntMultTruncate() to compute the product of two positive PetscInt and truncate to fit a PetscInt
2068:    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

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

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

2075:    Not available from Fortran

2077:    Level: advanced

2079: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2080: @*/
2081: PETSC_STATIC_INLINE PetscInt PetscRealIntMultTruncate(PetscReal a,PetscInt b)
2082: {
2083:   PetscInt64 r;

2085:   r  =  (PetscInt64) (a*(PetscReal)b);
2086:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2087:   return (PetscInt) r;
2088: }

2090: /*@C

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

2094:    Not Collective

2096:    Input Parameter:
2097: +     a - the PetscInt value
2098: -     b - the second value

2100:    Returns:
2101:       the result as a PetscInt value

2103:    Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2104:    Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2105:    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

2107:    Not available from Fortran

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

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

2113:    Level: advanced

2115: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2116: @*/
2117: PETSC_STATIC_INLINE PetscInt PetscIntMultTruncate(PetscInt a,PetscInt b)
2118: {
2119:   PetscInt64 r;

2121:   r  =  PetscInt64Mult(a,b);
2122:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2123:   return (PetscInt) r;
2124: }

2126: /*@C

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

2130:    Not Collective

2132:    Input Parameter:
2133: +     a - the PetscInt value
2134: -     b - the second value

2136:    Returns:
2137:      the result as a PetscInt value

2139:    Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2140:    Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2141:    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

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

2145:    Not available from Fortran

2147:    Level: advanced

2149: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2150: @*/
2151: PETSC_STATIC_INLINE PetscInt PetscIntSumTruncate(PetscInt a,PetscInt b)
2152: {
2153:   PetscInt64 r;

2155:   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2156:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2157:   return (PetscInt) r;
2158: }

2160: /*@C

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

2164:    Not Collective

2166:    Input Parameter:
2167: +     a - the PetscInt value
2168: -     b - the second value

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

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

2176:    Not available from Fortran

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

2180:    Level: advanced

2182: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64(), PetscIntSumError()
2183: @*/
2184: PETSC_STATIC_INLINE PetscErrorCode PetscIntMultError(PetscInt a,PetscInt b,PetscInt *result)
2185: {
2186:   PetscInt64 r;

2189:   r  =  PetscInt64Mult(a,b);
2190: #if !defined(PETSC_USE_64BIT_INDICES)
2191:   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);
2192: #endif
2193:   if (result) *result = (PetscInt) r;
2194:   return(0);
2195: }

2197: /*@C

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

2201:    Not Collective

2203:    Input Parameter:
2204: +     a - the PetscInt value
2205: -     b - the second value

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

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

2213:    Not available from Fortran

2215:    Level: advanced

2217: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2218: @*/
2219: PETSC_STATIC_INLINE PetscErrorCode PetscIntSumError(PetscInt a,PetscInt b,PetscInt *result)
2220: {
2221:   PetscInt64 r;

2224:   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2225: #if !defined(PETSC_USE_64BIT_INDICES)
2226:   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);
2227: #endif
2228:   if (result) *result = (PetscInt) r;
2229:   return(0);
2230: }

2232: /*
2233:      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2234: */
2235: #if defined(hz)
2236: #  undef hz
2237: #endif

2239: #include <limits.h>

2241: /* The number of bits in a byte */

2243: #define PETSC_BITS_PER_BYTE CHAR_BIT

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

2247: #if defined(PETSC_HAVE_SYS_PARAM_H)
2248: #  include <sys/param.h>
2249: #endif
2250: #if defined(PETSC_HAVE_SYS_TYPES_H)
2251: #  include <sys/types.h>
2252: #endif
2253: #if defined(MAXPATHLEN)
2254: #  define PETSC_MAX_PATH_LEN MAXPATHLEN
2255: #elif defined(MAX_PATH)
2256: #  define PETSC_MAX_PATH_LEN MAX_PATH
2257: #elif defined(_MAX_PATH)
2258: #  define PETSC_MAX_PATH_LEN _MAX_PATH
2259: #else
2260: #  define PETSC_MAX_PATH_LEN 4096
2261: #endif

2263: /*MC

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


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

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

2273:     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
2274:     where only a change in the major version number (x) indicates a change in the API.

2276:     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

2278:     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),
2279:     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
2280:     version number (x.y.z).

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

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

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

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

2290:     Level: intermediate

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

2294: M*/

2296: /*MC

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

2301: $
2302: $#include "petsc/finclude/petscXXX.h"
2303: $         use petscXXX

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

2307: $    XXX variablename
2308: $    type(tXXX) variablename

2310:     For example,

2312: $#include "petsc/finclude/petscvec.h"
2313: $         use petscvec
2314: $
2315: $    Vec b
2316: $    type(tVec) x

2318:     Level: beginner

2320: M*/

2322: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2323: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2324: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2325: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2326: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2327: PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2328: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2329: PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt*,PetscInt*,PetscInt*,PetscInt*);

2331: PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt,const PetscInt[],PetscBool*);
2332: PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt,const PetscMPIInt[],PetscBool*);
2333: PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt,const PetscReal[],PetscBool*);
2334: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2335: PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt,PetscInt[]);
2336: PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt*,PetscInt[]);
2337: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2339: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2340: PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt*);
2341: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2342: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2343: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2344: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt[],PetscInt[],PetscInt[]);
2345: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2346: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2347: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2348: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt,PetscMPIInt[],PetscInt[]);
2349: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2350: PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt,PetscInt[],void*,size_t,void*);
2351: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2352: PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt,PetscReal[],PetscInt[]);
2353: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2354: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2355: PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal,PetscInt,const PetscReal[],PetscReal,PetscInt*);
2356: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2357: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2358: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2359: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt[],const PetscInt[],PetscInt,const PetscInt[],const PetscInt[],PetscInt*,PetscInt**,PetscInt**);
2360: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt*,PetscInt**);
2361: PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt,const PetscMPIInt[],PetscInt,const PetscMPIInt[],PetscInt*,PetscMPIInt**);

2363: PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *);

2365: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2366: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);

2368: /*J
2369:     PetscRandomType - String with the name of a PETSc randomizer

2371:    Level: beginner

2373:    Notes:
2374:    To use SPRNG or RANDOM123 you must have ./configure PETSc
2375:    with the option --download-sprng or --download-random123

2377: .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2378: J*/
2379: typedef const char* PetscRandomType;
2380: #define PETSCRAND       "rand"
2381: #define PETSCRAND48     "rand48"
2382: #define PETSCSPRNG      "sprng"
2383: #define PETSCRANDER48   "rander48"
2384: #define PETSCRANDOM123  "random123"

2386: /* Logging support */
2387: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;

2389: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);

2391: /* Dynamic creation and loading functions */
2392: PETSC_EXTERN PetscFunctionList PetscRandomList;

2394: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2395: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2396: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2397: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType*);
2398: PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom,PetscObject,const char[]);
2399: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);

2401: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2402: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2403: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2404: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2405: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2406: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2407: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2408: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2409: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);

2411: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2412: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2413: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2414: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2415: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2416: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2417: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);
2418: PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2419: PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2420: PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);

2422: PETSC_STATIC_INLINE PetscBool PetscBinaryBigEndian(void) {long _petsc_v = 1; return ((char*)&_petsc_v)[0] ? PETSC_FALSE : PETSC_TRUE;}

2424: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscInt*,PetscDataType);
2425: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscInt*,PetscDataType);
2426: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,const void*,PetscInt,PetscDataType);
2427: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,const void*,PetscInt,PetscDataType);
2428: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2429: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2430: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2431: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2432: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2433: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2434: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2435: #if defined(PETSC_USE_SOCKET_VIEWER)
2436: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);
2437: #endif

2439: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2440: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2441: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);

2443: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2444: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool );
2445: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2446: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2447: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2448: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);

2450: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2451: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2452: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2453: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2454: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2455: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt*,const void*,PetscMPIInt*,PetscMPIInt**,void*)
2456:   PetscAttrMPIPointerWithType(6,3);
2457: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2458:                                                     PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2459:                                                     PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2460:   PetscAttrMPIPointerWithType(6,3);
2461: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2462:                                                        MPI_Request**,MPI_Request**,
2463:                                                        PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2464:                                                        PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2465:   PetscAttrMPIPointerWithType(6,3);

2467: PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2468: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2469: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);

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

2473: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);

2475: PETSC_EXTERN const char *const PetscSubcommTypes[];

2477: struct _n_PetscSubcomm {
2478:   MPI_Comm         parent;           /* parent communicator */
2479:   MPI_Comm         dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2480:   MPI_Comm         child;            /* the sub-communicator */
2481:   PetscMPIInt      n;                /* num of subcommunicators under the parent communicator */
2482:   PetscMPIInt      color;            /* color of processors belong to this communicator */
2483:   PetscMPIInt      *subsize;         /* size of subcommunicator[color] */
2484:   PetscSubcommType type;
2485:   char             *subcommprefix;
2486: };

2488: PETSC_STATIC_INLINE MPI_Comm PetscSubcommParent(PetscSubcomm scomm) {return scomm->parent;}
2489: PETSC_STATIC_INLINE MPI_Comm PetscSubcommChild(PetscSubcomm scomm) {return scomm->child;}
2490: PETSC_STATIC_INLINE MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm) {return scomm->dupparent;}
2491: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2492: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2493: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2494: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2495: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2496: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2497: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2498: PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm,const char[]);
2499: PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm,MPI_Comm*);
2500: PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm,MPI_Comm*);
2501: PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm,MPI_Comm*);

2503: PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt,PetscHeap*);
2504: PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap,PetscInt,PetscInt);
2505: PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap,PetscInt*,PetscInt*);
2506: PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap,PetscInt*,PetscInt*);
2507: PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap,PetscInt,PetscInt);
2508: PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2509: PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap*);
2510: PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap,PetscViewer);

2512: PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2513: PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm,PetscShmComm*);
2514: PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2515: PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2516: PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm,MPI_Comm*);

2518: /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2519: PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm,PetscInt,PetscOmpCtrl*);
2520: PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl,MPI_Comm*,MPI_Comm*,PetscBool*);
2521: PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl*);
2522: PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2523: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2524: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);

2526: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2527: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2528: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2529: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2530: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2531: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2532: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2533: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);


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

2541: extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2542: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted*);
2543: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted*);
2544: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted,const char*,const char*,PetscBool*);

2546: #include <stdarg.h>
2547: PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list);
2548: PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE*,const char[],va_list);

2550: PETSC_EXTERN PetscSegBuffer PetscCitationsList;

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

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

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

2561:    Level: intermediate

2563:    Not available from Fortran

2565:      Options Database:
2566: .     -citations [filename]   - print out the bibtex entries for the given computation
2567: @*/
2568: PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2569: {
2570:   size_t         len;
2571:   char           *vstring;

2575:   if (set && *set) return(0);
2576:   PetscStrlen(cit,&len);
2577:   PetscSegBufferGet(PetscCitationsList,len,&vstring);
2578:   PetscArraycpy(vstring,cit,len);
2579:   if (set) *set = PETSC_TRUE;
2580:   return(0);
2581: }

2583: PETSC_EXTERN PetscErrorCode PetscURLShorten(const char[],char[],size_t);
2584: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm,char[],char[],size_t);
2585: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm,const char[],char[],size_t);
2586: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm,const char[],const char []);

2588: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm,char[],char[],size_t);
2589: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm,const char[],char[],char[],size_t);

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

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

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


2600: #if defined(PETSC_USE_DEBUG)
2601: /*
2602:    Verify that all processes in the communicator have called this from the same line of code
2603:  */
2604: PETSC_EXTERN PetscErrorCode PetscAllreduceBarrierCheck(MPI_Comm,PetscMPIInt,int,const char*,const char *);

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

2611:    Collective

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

2617:    Input Parameters:
2618: +  indata - pointer to the input data to be reduced
2619: .  count - the number of MPI data items in indata and outdata
2620: .  datatype - the MPI datatype, for example MPI_INT
2621: .  op - the MPI operation, for example MPI_SUM
2622: -  comm - the MPI communicator on which the operation occurs

2624:    Output Parameter:
2625: .  outdata - the reduced values

2627:    Notes:
2628:    In optimized mode this directly calls MPI_Allreduce()

2630:    Level: developer

2632: .seealso: MPI_Allreduce()
2633: M*/
2634: #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))
2635: #else
2636: #define MPIU_Allreduce(a,b,c,d,e,fcomm) MPI_Allreduce(a,b,c,d,e,fcomm)
2637: #endif

2639: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
2640: PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint,PetscMPIInt,MPI_Info,MPI_Comm,void*,MPI_Win*);
2641: PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win,PetscMPIInt,MPI_Aint*,PetscMPIInt*,void*);
2642: #endif

2644: /*
2645:     Returned from PETSc functions that are called from MPI, such as related to attributes
2646: */
2647: PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CLASS;
2648: PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CODE;

2650: #endif