Actual source code: petscsys.h

petsc-master 2020-01-28
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: /*MC
435:    PetscMalloc - Allocates memory, One should use PetscNew(), PetscMalloc1() or PetscCalloc1() usually instead of this

437:    Synopsis:
438:     #include <petscsys.h>
439:    PetscErrorCode PetscMalloc(size_t m,void **result)

441:    Not Collective

443:    Input Parameter:
444: .  m - number of bytes to allocate

446:    Output Parameter:
447: .  result - memory allocated

449:    Level: beginner

451:    Notes:
452:    Memory is always allocated at least double aligned

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

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

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

461: /*MC
462:    PetscRealloc - Rellocates memory

464:    Synopsis:
465:     #include <petscsys.h>
466:    PetscErrorCode PetscRealloc(size_t m,void **result)

468:    Not Collective

470:    Input Parameters:
471: +  m - number of bytes to allocate
472: -  result - previous memory

474:    Output Parameter:
475: .  result - new memory allocated

477:    Level: developer

479:    Notes:
480:    Memory is always allocated at least double aligned

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

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

487: /*MC
488:    PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment

490:    Synopsis:
491:     #include <petscsys.h>
492:    void *PetscAddrAlign(void *addr)

494:    Not Collective

496:    Input Parameters:
497: .  addr - address to align (any pointer type)

499:    Level: developer

501: .seealso: PetscMallocAlign()

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

506: /*MC
507:    PetscMalloc1 - Allocates an array of memory aligned to PETSC_MEMALIGN

509:    Synopsis:
510:     #include <petscsys.h>
511:    PetscErrorCode PetscMalloc1(size_t m1,type **r1)

513:    Not Collective

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

518:    Output Parameter:
519: .  r1 - memory allocated

521:    Note:
522:    This uses the sizeof() of the memory type requested to determine the total memory to be allocated, therefore you should not
523:          multiply the number of elements requested by the sizeof() the type. For example use
524: $  PetscInt *id;
525: $  PetscMalloc1(10,&id);
526:         not
527: $  PetscInt *id;
528: $  PetscMalloc1(10*sizeof(PetscInt),&id);

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

532:    Level: beginner

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

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

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

542:    Synopsis:
543:     #include <petscsys.h>
544:    PetscErrorCode PetscCalloc1(size_t m1,type **r1)

546:    Not Collective

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

551:    Output Parameter:
552: .  r1 - memory allocated

554:    Notes:
555:    See PetsMalloc1() for more details on usage.

557:    Level: beginner

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

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

564: /*MC
565:    PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN

567:    Synopsis:
568:     #include <petscsys.h>
569:    PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)

571:    Not Collective

573:    Input Parameter:
574: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
575: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

577:    Output Parameter:
578: +  r1 - memory allocated in first chunk
579: -  r2 - memory allocated in second chunk

581:    Level: developer

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

585: M*/
586: #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))

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

591:    Synopsis:
592:     #include <petscsys.h>
593:    PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)

595:    Not Collective

597:    Input Parameter:
598: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
599: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

601:    Output Parameter:
602: +  r1 - memory allocated in first chunk
603: -  r2 - memory allocated in second chunk

605:    Level: developer

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

609: M*/
610: #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))

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

615:    Synopsis:
616:     #include <petscsys.h>
617:    PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

619:    Not Collective

621:    Input Parameter:
622: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
623: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
624: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

626:    Output Parameter:
627: +  r1 - memory allocated in first chunk
628: .  r2 - memory allocated in second chunk
629: -  r3 - memory allocated in third chunk

631:    Level: developer

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

635: M*/
636: #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))

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

641:    Synopsis:
642:     #include <petscsys.h>
643:    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

645:    Not Collective

647:    Input Parameter:
648: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
649: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
650: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

652:    Output Parameter:
653: +  r1 - memory allocated in first chunk
654: .  r2 - memory allocated in second chunk
655: -  r3 - memory allocated in third chunk

657:    Level: developer

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

661: M*/
662: #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))

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

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

671:    Not Collective

673:    Input Parameter:
674: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
675: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
676: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
677: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

679:    Output Parameter:
680: +  r1 - memory allocated in first chunk
681: .  r2 - memory allocated in second chunk
682: .  r3 - memory allocated in third chunk
683: -  r4 - memory allocated in fourth chunk

685:    Level: developer

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

689: M*/
690: #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))

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

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

699:    Not Collective

701:    Input Parameters:
702: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
703: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
704: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
705: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

707:    Output Parameters:
708: +  r1 - memory allocated in first chunk
709: .  r2 - memory allocated in second chunk
710: .  r3 - memory allocated in third chunk
711: -  r4 - memory allocated in fourth chunk

713:    Level: developer

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

717: M*/
718: #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))

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

723:    Synopsis:
724:     #include <petscsys.h>
725:    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)

727:    Not Collective

729:    Input Parameters:
730: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
731: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
732: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
733: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
734: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

736:    Output Parameters:
737: +  r1 - memory allocated in first chunk
738: .  r2 - memory allocated in second chunk
739: .  r3 - memory allocated in third chunk
740: .  r4 - memory allocated in fourth chunk
741: -  r5 - memory allocated in fifth chunk

743:    Level: developer

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

747: M*/
748: #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))

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

753:    Synopsis:
754:     #include <petscsys.h>
755:    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)

757:    Not Collective

759:    Input Parameters:
760: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
761: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
762: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
763: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
764: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

766:    Output Parameters:
767: +  r1 - memory allocated in first chunk
768: .  r2 - memory allocated in second chunk
769: .  r3 - memory allocated in third chunk
770: .  r4 - memory allocated in fourth chunk
771: -  r5 - memory allocated in fifth chunk

773:    Level: developer

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

777: M*/
778: #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))

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

783:    Synopsis:
784:     #include <petscsys.h>
785:    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)

787:    Not Collective

789:    Input Parameters:
790: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
791: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
792: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
793: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
794: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
795: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

797:    Output Parameteasr:
798: +  r1 - memory allocated in first chunk
799: .  r2 - memory allocated in second chunk
800: .  r3 - memory allocated in third chunk
801: .  r4 - memory allocated in fourth chunk
802: .  r5 - memory allocated in fifth chunk
803: -  r6 - memory allocated in sixth chunk

805:    Level: developer

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

809: M*/
810: #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))

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

815:    Synopsis:
816:     #include <petscsys.h>
817:    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)

819:    Not Collective

821:    Input Parameters:
822: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
823: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
824: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
825: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
826: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
827: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

829:    Output Parameters:
830: +  r1 - memory allocated in first chunk
831: .  r2 - memory allocated in second chunk
832: .  r3 - memory allocated in third chunk
833: .  r4 - memory allocated in fourth chunk
834: .  r5 - memory allocated in fifth chunk
835: -  r6 - memory allocated in sixth chunk

837:    Level: developer

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

841: M*/
842: #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))

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

847:    Synopsis:
848:     #include <petscsys.h>
849:    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)

851:    Not Collective

853:    Input Parameters:
854: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
855: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
856: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
857: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
858: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
859: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
860: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

862:    Output Parameters:
863: +  r1 - memory allocated in first chunk
864: .  r2 - memory allocated in second chunk
865: .  r3 - memory allocated in third chunk
866: .  r4 - memory allocated in fourth chunk
867: .  r5 - memory allocated in fifth chunk
868: .  r6 - memory allocated in sixth chunk
869: -  r7 - memory allocated in seventh chunk

871:    Level: developer

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

875: M*/
876: #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))

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

881:    Synopsis:
882:     #include <petscsys.h>
883:    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)

885:    Not Collective

887:    Input Parameters:
888: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
889: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
890: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
891: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
892: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
893: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
894: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

896:    Output Parameters:
897: +  r1 - memory allocated in first chunk
898: .  r2 - memory allocated in second chunk
899: .  r3 - memory allocated in third chunk
900: .  r4 - memory allocated in fourth chunk
901: .  r5 - memory allocated in fifth chunk
902: .  r6 - memory allocated in sixth chunk
903: -  r7 - memory allocated in seventh chunk

905:    Level: developer

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

909: M*/
910: #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))

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

915:    Synopsis:
916:     #include <petscsys.h>
917:    PetscErrorCode PetscNew(type **result)

919:    Not Collective

921:    Output Parameter:
922: .  result - memory allocated, sized to match pointer type

924:    Level: beginner

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

928: M*/
929: #define PetscNew(b)      PetscCalloc1(1,(b))

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

935:    Synopsis:
936:     #include <petscsys.h>
937:    PetscErrorCode PetscNewLog(PetscObject obj,type **result)

939:    Not Collective

941:    Input Parameter:
942: .  obj - object memory is logged to

944:    Output Parameter:
945: .  result - memory allocated, sized to match pointer type

947:    Level: developer

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

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

954: /*MC
955:    PetscFree - Frees memory

957:    Synopsis:
958:     #include <petscsys.h>
959:    PetscErrorCode PetscFree(void *memory)

961:    Not Collective

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

966:    Level: beginner

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

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

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

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

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

981:    Synopsis:
982:     #include <petscsys.h>
983:    PetscErrorCode PetscFree2(void *memory1,void *memory2)

985:    Not Collective

987:    Input Parameters:
988: +   memory1 - memory to free
989: -   memory2 - 2nd memory to free

991:    Level: developer

993:    Notes:
994:     Memory must have been obtained with PetscMalloc2()

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

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

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

1004:    Synopsis:
1005:     #include <petscsys.h>
1006:    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)

1008:    Not Collective

1010:    Input Parameters:
1011: +   memory1 - memory to free
1012: .   memory2 - 2nd memory to free
1013: -   memory3 - 3rd memory to free

1015:    Level: developer

1017:    Notes:
1018:     Memory must have been obtained with PetscMalloc3()

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

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

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

1028:    Synopsis:
1029:     #include <petscsys.h>
1030:    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)

1032:    Not Collective

1034:    Input Parameters:
1035: +   m1 - memory to free
1036: .   m2 - 2nd memory to free
1037: .   m3 - 3rd memory to free
1038: -   m4 - 4th memory to free

1040:    Level: developer

1042:    Notes:
1043:     Memory must have been obtained with PetscMalloc4()

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

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

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

1053:    Synopsis:
1054:     #include <petscsys.h>
1055:    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)

1057:    Not Collective

1059:    Input Parameters:
1060: +   m1 - memory to free
1061: .   m2 - 2nd memory to free
1062: .   m3 - 3rd memory to free
1063: .   m4 - 4th memory to free
1064: -   m5 - 5th memory to free

1066:    Level: developer

1068:    Notes:
1069:     Memory must have been obtained with PetscMalloc5()

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

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

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

1079:    Synopsis:
1080:     #include <petscsys.h>
1081:    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)

1083:    Not Collective

1085:    Input Parameters:
1086: +   m1 - memory to free
1087: .   m2 - 2nd memory to free
1088: .   m3 - 3rd memory to free
1089: .   m4 - 4th memory to free
1090: .   m5 - 5th memory to free
1091: -   m6 - 6th memory to free


1094:    Level: developer

1096:    Notes:
1097:     Memory must have been obtained with PetscMalloc6()

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

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

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

1107:    Synopsis:
1108:     #include <petscsys.h>
1109:    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)

1111:    Not Collective

1113:    Input Parameters:
1114: +   m1 - memory to free
1115: .   m2 - 2nd memory to free
1116: .   m3 - 3rd memory to free
1117: .   m4 - 4th memory to free
1118: .   m5 - 5th memory to free
1119: .   m6 - 6th memory to free
1120: -   m7 - 7th memory to free


1123:    Level: developer

1125:    Notes:
1126:     Memory must have been obtained with PetscMalloc7()

1128: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1129:           PetscMalloc7()

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

1134: PETSC_EXTERN PetscErrorCode PetscMallocA(int,PetscBool,int,const char *,const char *,size_t,void *,...);
1135: PETSC_EXTERN PetscErrorCode PetscFreeA(int,int,const char *,const char *,void *,...);
1136: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,PetscBool,int,const char[],const char[],void**);
1137: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1138: PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t,int,const char[],const char[],void**);
1139: PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1140: 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 **));
1141: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);

1143: /*
1144:   Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1145: */
1146: PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1147: PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);

1149: #define MPIU_PETSCLOGDOUBLE  MPI_DOUBLE
1150: #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION

1152: /*
1153:    Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1154: */
1155: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1156: PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1157: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1158: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1159: PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1160: PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int,PetscLogDouble*);
1161: PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool,PetscBool);
1162: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*,PetscBool*,PetscBool*);
1163: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1164: PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1165: PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool*);

1167: PETSC_EXTERN const char *const PetscDataTypes[];
1168: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1169: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1170: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1171: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);

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

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

1209: PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1210: PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1211: PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);
1212: PETSC_EXTERN PetscErrorCode PetscStrInList(const char[],const char[],char,PetscBool*);

1214: PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1215: PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);

1217: /*
1218:    These are MPI operations for MPI_Allreduce() etc
1219: */
1220: PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1221: #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1222: PETSC_EXTERN MPI_Op MPIU_SUM;
1223: #else
1224: #define MPIU_SUM MPI_SUM
1225: #endif
1226: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1227: PETSC_EXTERN MPI_Op MPIU_MAX;
1228: PETSC_EXTERN MPI_Op MPIU_MIN;
1229: #else
1230: #define MPIU_MAX MPI_MAX
1231: #define MPIU_MIN MPI_MIN
1232: #endif
1233: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);

1235: PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1236: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);

1238: PETSC_EXTERN const char *const PetscFileModes[];

1240: /*
1241:     Defines PETSc error handling.
1242: */
1243:  #include <petscerror.h>

1245: #define PETSC_SMALLEST_CLASSID  1211211
1246: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1247: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1248: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);
1249: PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject,PetscObjectId*);
1250: PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject,PetscObjectId,PetscBool*);


1253: /*
1254:    Routines that get memory usage information from the OS
1255: */
1256: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1257: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1258: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1259: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);

1261: PETSC_EXTERN PetscErrorCode PetscInfoAllow(PetscBool ,const char []);
1262: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);

1264: /*
1265:    Initialization of PETSc
1266: */
1267: PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1268: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1269: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1270: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1271: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1272: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1273: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1274: PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1275: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1276: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);

1278: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1279: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);

1281: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1282: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1283: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1284: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);

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


1289: /*
1290:      These are so that in extern C code we can caste function pointers to non-extern C
1291:    function pointers. Since the regular C++ code expects its function pointers to be C++
1292: */
1293: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1294: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1295: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);

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

1326:  #include <petscviewertypes.h>
1327:  #include <petscoptions.h>

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

1331: PETSC_EXTERN PetscErrorCode PetscMemoryShowUsage(PetscViewer,const char[]);
1332: PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer,const char[]);
1333: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1334: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1335: #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1336: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1337: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1338: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1339: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1340: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1341: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1342: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1343: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1344: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,PetscObject,const char[]);
1345: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1346: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1347: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject,const char[],PetscBool *);
1348: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1349: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1350: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1351: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);

1353: #if defined(PETSC_HAVE_SAWS)
1354: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1355: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1356: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1357: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1358: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1359: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1360: PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1361: PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1362: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1363: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);

1365: #else
1366: #define PetscSAWsBlock()                        0
1367: #define PetscObjectSAWsViewOff(obj)             0
1368: #define PetscObjectSAWsSetBlock(obj,flg)        0
1369: #define PetscObjectSAWsBlock(obj)               0
1370: #define PetscObjectSAWsGrantAccess(obj)         0
1371: #define PetscObjectSAWsTakeAccess(obj)          0
1372: #define PetscStackViewSAWs()                    0
1373: #define PetscStackSAWsViewOff()                 0
1374: #define PetscStackSAWsTakeAccess()
1375: #define PetscStackSAWsGrantAccess()

1377: #endif

1379: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1380: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1381: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);

1383: #if defined(PETSC_USE_DEBUG)
1384: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1385: #endif
1386: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);

1388: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1389: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1390: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1391: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1392: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1393: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);

1395: /*
1396:     Dynamic library lists. Lists of names of routines in objects or in dynamic
1397:   link libraries that will be loaded as needed.
1398: */

1400: #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1401: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1402: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1403: #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1404: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1405: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[],const char[]);
1406: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1407: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1408: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);

1410: PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1411: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1412: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1413: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1414: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1415: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1416: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1417: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);

1419: /*
1420:      Useful utility routines
1421: */
1422: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1423: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1424: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1425: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1426: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1427: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);
1428: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm,PetscInt[],PetscInt[]);
1429: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm,PetscReal[],PetscReal[]);

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

1434:     Notes:
1435:     This is useful in cases like
1436: $     int        *a;
1437: $     PetscBool  flag = PetscNot(a)
1438:      where !a would not return a PetscBool because we cannot provide a cast from int to PetscBool in C.

1440:     Level: beginner

1442:     .seealso : PetscBool, PETSC_TRUE, PETSC_FALSE
1443: M*/
1444: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)

1446: /*MC
1447:     PetscHelpPrintf - Prints help messages.

1449:    Synopsis:
1450:     #include <petscsys.h>
1451:      PetscErrorCode (*PetscHelpPrintf)(const char format[],...);

1453:     Not Collective

1455:     Input Parameters:
1456: .   format - the usual printf() format string

1458:    Level: developer

1460:     Fortran Note:
1461:     This routine is not supported in Fortran.


1464: .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1465: M*/
1466: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);

1468: /*
1469:      Defines PETSc profiling.
1470: */
1471:  #include <petsclog.h>

1473: /*
1474:       Simple PETSc parallel IO for ASCII printing
1475: */
1476: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1477: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1478: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1479: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1480: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1481: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1482: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);
1483: PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[],size_t,const char*,PetscInt,const PetscReal[]);

1485: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1486: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1487: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);

1489: PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char*,size_t*);
1490: PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char*,char *);

1492: #if defined(PETSC_HAVE_POPEN)
1493: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1494: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*);
1495: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1496: #endif

1498: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1499: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1500: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*);
1501: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1502: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1503: PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1504: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);

1506: PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1507: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void **);
1508: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void *);
1509: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1510: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer *);
1511: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void*));
1512: PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void*);

1514: /*
1515:    For use in debuggers
1516: */
1517: PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1518: PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1519: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1520: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1521: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);

1523: #include <stddef.h>
1524: #include <string.h>             /* for memcpy, memset */
1525: #include <stdlib.h>

1527: #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1528: #include <xmmintrin.h>
1529: #endif

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

1536:    Not Collective

1538:    Input Parameters:
1539: +  b - pointer to initial memory space
1540: .  a - pointer to copy space
1541: -  n - length (in bytes) of space to copy

1543:    Level: intermediate

1545:    Note:
1546:    PetscArraymove() is preferred
1547:    This routine is analogous to memmove().

1549:    Developers Note: This is inlined for performance

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

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

1584:    Not Collective

1586:    Input Parameters:
1587: +  b - pointer to initial memory space
1588: -  n - length (in bytes) of space to copy

1590:    Output Parameter:
1591: .  a - pointer to copy space

1593:    Level: intermediate

1595:    Compile Option:
1596:     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1597:                                   for memory copies on double precision values.
1598:     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1599:                                   for memory copies on double precision values.
1600:     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1601:                                   for memory copies on double precision values.

1603:    Note:
1604:    Prefer PetscArraycpy()
1605:    This routine is analogous to memcpy().
1606:    Not available from Fortran

1608:    Developer Note: this is inlined for fastest performance

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

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

1655: /*@C
1656:    PetscMemzero - Zeros the specified memory.

1658:    Not Collective

1660:    Input Parameters:
1661: +  a - pointer to beginning memory location
1662: -  n - length (in bytes) of memory to initialize

1664:    Level: intermediate

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

1670:    Not available from Fortran
1671:    Prefer PetscArrayzero()

1673:    Developer Note: this is inlined for fastest performance

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

1707: /*MC
1708:    PetscArraycmp - Compares two arrays in memory.

1710:    Synopsis:
1711:     #include <petscsys.h>
1712:     PetscErrorCode PetscArraycmp(const anytype *str1,const anytype *str2,size_t cnt,PetscBool *e)

1714:    Not Collective

1716:    Input Parameters:
1717: +  str1 - First array
1718: .  str2 - Second array
1719: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1721:    Output Parameters:
1722: .   e - PETSC_TRUE if equal else PETSC_FALSE.

1724:    Level: intermediate

1726:    Note:
1727:    This routine is a preferred replacement to PetscMemcmp()
1728:    The arrays must be of the same type

1730: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy(),
1731:           PetscArraymove()
1732: M*/
1733: #define  PetscArraycmp(str1,str2,cnt,e) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcmp(str1,str2,(size_t)(cnt)*sizeof(*(str1)),e))

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

1739:    Synopsis:
1740:     #include <petscsys.h>
1741:     PetscErrorCode PetscArraymove(anytype *str1,const anytype *str2,size_t cnt)

1743:    Not Collective

1745:    Input Parameters:
1746: +  str1 - First array
1747: .  str2 - Second array
1748: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1750:    Level: intermediate

1752:    Note:
1753:    This routine is a preferred replacement to PetscMemmove()
1754:    The arrays must be of the same type

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

1760: /*MC
1761:    PetscArraycpy - Copies from one array in memory to another

1763:    Synopsis:
1764:     #include <petscsys.h>
1765:     PetscErrorCode PetscArraycpy(anytype *str1,const anytype *str2,size_t cnt)

1767:    Not Collective

1769:    Input Parameters:
1770: +  str1 - First array (destination)
1771: .  str2 - Second array (source)
1772: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1774:    Level: intermediate

1776:    Note:
1777:    This routine is a preferred replacement to PetscMemcpy()
1778:    The arrays must be of the same type

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

1784: /*MC
1785:    PetscArrayzero - Zeros an array in memory.

1787:    Synopsis:
1788:     #include <petscsys.h>
1789:     PetscErrorCode PetscArrayzero(anytype *str1,size_t cnt)

1791:    Not Collective

1793:    Input Parameters:
1794: +  str1 - array
1795: -  cnt  - Count of the array, not in bytes, but number of entries in the array

1797:    Level: intermediate

1799:    Note:
1800:    This routine is a preferred replacement to PetscMemzero()

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

1806: /*MC
1807:    PetscPrefetchBlock - Prefetches a block of memory

1809:    Synopsis:
1810:     #include <petscsys.h>
1811:     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)

1813:    Not Collective

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

1821:    Level: developer

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

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

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

1836: M*/
1837: #define PetscPrefetchBlock(a,n,rw,t) do {                               \
1838:     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
1839:     for ( ; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
1840:   } while (0)

1842: /*
1843:       Determine if some of the kernel computation routines use
1844:    Fortran (rather than C) for the numerical calculations. On some machines
1845:    and compilers (like complex numbers) the Fortran version of the routines
1846:    is faster than the C/C++ versions. The flag --with-fortran-kernels
1847:    should be used with ./configure to turn these on.
1848: */
1849: #if defined(PETSC_USE_FORTRAN_KERNELS)

1851: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1852: #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1853: #endif

1855: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
1856: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
1857: #endif

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

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

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

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

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

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

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

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

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

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

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

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

1907: #endif

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

1924: /* --------------------------------------------------------------------*/

1926: /*MC
1927:     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
1928:         communication

1930:    Level: beginner

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

1934: .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
1935: M*/

1937: #if defined(PETSC_HAVE_MPIIO)
1938: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1939: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
1940: #endif

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

1944: /* Limit MPI to 32-bits */
1945: #define PETSC_MPI_INT_MAX  2147483647
1946: #define PETSC_MPI_INT_MIN -2147483647
1947: /* Limit BLAS to 32-bits */
1948: #define PETSC_BLAS_INT_MAX  2147483647
1949: #define PETSC_BLAS_INT_MIN -2147483647

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

1955:    Not Collective

1957:    Input Parameter:
1958: .     a - the PetscInt64 value

1960:    Output Parameter:
1961: .     b - the resulting PetscInt value

1963:    Level: advanced

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

1967:    Not available from Fortran

1969: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscBLASIntCast()
1970: @*/
1971: PETSC_STATIC_INLINE PetscErrorCode PetscIntCast(PetscInt64 a,PetscInt *b)
1972: {
1974: #if !defined(PETSC_USE_64BIT_INDICES)
1975:   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");
1976: #endif
1977:   *b =  (PetscInt)(a);
1978:   return(0);
1979: }

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

1985:    Not Collective

1987:    Input Parameter:
1988: .     a - the PetscInt value

1990:    Output Parameter:
1991: .     b - the resulting PetscBLASInt value

1993:    Level: advanced

1995:    Not available from Fortran

1997: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscIntCast()
1998: @*/
1999: PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
2000: {
2002: #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
2003:   *b = 0;
2004:   if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK");
2005: #endif
2006:   *b =  (PetscBLASInt)(a);
2007:   return(0);
2008: }

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

2014:    Not Collective

2016:    Input Parameter:
2017: .     a - the PetscInt value

2019:    Output Parameter:
2020: .     b - the resulting PetscMPIInt value

2022:    Level: advanced

2024:    Not available from Fortran

2026: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntCast()
2027: @*/
2028: PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
2029: {
2031: #if defined(PETSC_USE_64BIT_INDICES)
2032:   *b = 0;
2033:   if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for MPI");
2034: #endif
2035:   *b =  (PetscMPIInt)(a);
2036:   return(0);
2037: }

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

2041: /*@C

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

2045:    Not Collective

2047:    Input Parameter:
2048: +     a - the PetscReal value
2049: -     b - the second value

2051:    Returns:
2052:       the result as a PetscInt value

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

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

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

2063:    Not available from Fortran

2065:    Level: advanced

2067: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2068: @*/
2069: PETSC_STATIC_INLINE PetscInt PetscRealIntMultTruncate(PetscReal a,PetscInt b)
2070: {
2071:   PetscInt64 r;

2073:   r  =  (PetscInt64) (a*(PetscReal)b);
2074:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2075:   return (PetscInt) r;
2076: }

2078: /*@C

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

2082:    Not Collective

2084:    Input Parameter:
2085: +     a - the PetscInt value
2086: -     b - the second value

2088:    Returns:
2089:       the result as a PetscInt value

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

2095:    Not available from Fortran

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

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

2101:    Level: advanced

2103: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2104: @*/
2105: PETSC_STATIC_INLINE PetscInt PetscIntMultTruncate(PetscInt a,PetscInt b)
2106: {
2107:   PetscInt64 r;

2109:   r  =  PetscInt64Mult(a,b);
2110:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2111:   return (PetscInt) r;
2112: }

2114: /*@C

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

2118:    Not Collective

2120:    Input Parameter:
2121: +     a - the PetscInt value
2122: -     b - the second value

2124:    Returns:
2125:      the result as a PetscInt value

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

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

2133:    Not available from Fortran

2135:    Level: advanced

2137: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2138: @*/
2139: PETSC_STATIC_INLINE PetscInt PetscIntSumTruncate(PetscInt a,PetscInt b)
2140: {
2141:   PetscInt64 r;

2143:   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2144:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2145:   return (PetscInt) r;
2146: }

2148: /*@C

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

2152:    Not Collective

2154:    Input Parameter:
2155: +     a - the PetscInt value
2156: -     b - the second value

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

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

2164:    Not available from Fortran

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

2168:    Level: advanced

2170: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64(), PetscIntSumError()
2171: @*/
2172: PETSC_STATIC_INLINE PetscErrorCode PetscIntMultError(PetscInt a,PetscInt b,PetscInt *result)
2173: {
2174:   PetscInt64 r;

2177:   r  =  PetscInt64Mult(a,b);
2178: #if !defined(PETSC_USE_64BIT_INDICES)
2179:   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);
2180: #endif
2181:   if (result) *result = (PetscInt) r;
2182:   return(0);
2183: }

2185: /*@C

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

2189:    Not Collective

2191:    Input Parameter:
2192: +     a - the PetscInt value
2193: -     b - the second value

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

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

2201:    Not available from Fortran

2203:    Level: advanced

2205: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2206: @*/
2207: PETSC_STATIC_INLINE PetscErrorCode PetscIntSumError(PetscInt a,PetscInt b,PetscInt *result)
2208: {
2209:   PetscInt64 r;

2212:   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2213: #if !defined(PETSC_USE_64BIT_INDICES)
2214:   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);
2215: #endif
2216:   if (result) *result = (PetscInt) r;
2217:   return(0);
2218: }

2220: /*
2221:      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2222: */
2223: #if defined(hz)
2224: #  undef hz
2225: #endif

2227: #include <limits.h>

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

2231: #define PETSC_BITS_PER_BYTE CHAR_BIT

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

2235: #if defined(PETSC_HAVE_SYS_PARAM_H)
2236: #  include <sys/param.h>
2237: #endif
2238: #if defined(PETSC_HAVE_SYS_TYPES_H)
2239: #  include <sys/types.h>
2240: #endif
2241: #if defined(MAXPATHLEN)
2242: #  define PETSC_MAX_PATH_LEN MAXPATHLEN
2243: #elif defined(MAX_PATH)
2244: #  define PETSC_MAX_PATH_LEN MAX_PATH
2245: #elif defined(_MAX_PATH)
2246: #  define PETSC_MAX_PATH_LEN _MAX_PATH
2247: #else
2248: #  define PETSC_MAX_PATH_LEN 4096
2249: #endif

2251: /*MC

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


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

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

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

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

2266:     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),
2267:     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
2268:     version number (x.y.z).

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

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

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

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

2278:     Level: intermediate

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

2282: M*/

2284: /*MC

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

2289: $
2290: $#include "petsc/finclude/petscXXX.h"
2291: $         use petscXXX

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

2295: $    XXX variablename
2296: $    type(tXXX) variablename

2298:     For example,

2300: $#include "petsc/finclude/petscvec.h"
2301: $         use petscvec
2302: $
2303: $    Vec b
2304: $    type(tVec) x

2306:     Level: beginner

2308: M*/

2310: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2311: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2312: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2313: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2314: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2315: PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2316: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2317: PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt*,PetscInt*,PetscInt*,PetscInt*);

2319: PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt,const PetscInt[],PetscBool*);
2320: PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt,const PetscMPIInt[],PetscBool*);
2321: PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt,const PetscReal[],PetscBool*);
2322: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2323: PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt,PetscInt[]);
2324: PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt*,PetscInt[]);
2325: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2327: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2328: PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt*);
2329: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2330: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2331: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2332: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt[],PetscInt[],PetscInt[]);
2333: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2334: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2335: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2336: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt,PetscMPIInt[],PetscInt[]);
2337: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2338: PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt,PetscInt[],void*,size_t,void*);
2339: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2340: PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt,PetscReal[],PetscInt[]);
2341: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2342: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2343: PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal,PetscInt,const PetscReal[],PetscReal,PetscInt*);
2344: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2345: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2346: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2347: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt[],const PetscInt[],PetscInt,const PetscInt[],const PetscInt[],PetscInt*,PetscInt**,PetscInt**);
2348: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt*,PetscInt**);
2349: PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt,const PetscMPIInt[],PetscInt,const PetscMPIInt[],PetscInt*,PetscMPIInt**);

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

2353: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2354: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);

2356: /*J
2357:     PetscRandomType - String with the name of a PETSc randomizer

2359:    Level: beginner

2361:    Notes:
2362:    To use SPRNG or RANDOM123 you must have ./configure PETSc
2363:    with the option --download-sprng or --download-random123

2365: .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2366: J*/
2367: typedef const char* PetscRandomType;
2368: #define PETSCRAND       "rand"
2369: #define PETSCRAND48     "rand48"
2370: #define PETSCSPRNG      "sprng"
2371: #define PETSCRANDER48   "rander48"
2372: #define PETSCRANDOM123  "random123"

2374: /* Logging support */
2375: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;

2377: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);

2379: /* Dynamic creation and loading functions */
2380: PETSC_EXTERN PetscFunctionList PetscRandomList;

2382: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2383: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2384: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2385: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType*);
2386: PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom,PetscObject,const char[]);
2387: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);

2389: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2390: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2391: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2392: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2393: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2394: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2395: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2396: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2397: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);

2399: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2400: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2401: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2402: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2403: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2404: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2405: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);
2406: PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2407: PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2408: PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);

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

2412: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscInt*,PetscDataType);
2413: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscInt*,PetscDataType);
2414: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,void*,PetscInt,PetscDataType,PetscBool);
2415: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,void*,PetscInt,PetscDataType,PetscBool);
2416: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2417: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2418: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2419: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2420: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2421: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2422: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2423: #if defined(PETSC_USE_SOCKET_VIEWER)
2424: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);
2425: #endif

2427: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2428: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2429: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);

2431: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2432: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool );
2433: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2434: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2435: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2436: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);

2438: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2439: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2440: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2441: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2442: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2443: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt*,const void*,PetscMPIInt*,PetscMPIInt**,void*)
2444:   PetscAttrMPIPointerWithType(6,3);
2445: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2446:                                                     PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2447:                                                     PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2448:   PetscAttrMPIPointerWithType(6,3);
2449: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2450:                                                        MPI_Request**,MPI_Request**,
2451:                                                        PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2452:                                                        PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2453:   PetscAttrMPIPointerWithType(6,3);

2455: PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2456: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2457: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);

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

2461: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);

2463: PETSC_EXTERN const char *const PetscSubcommTypes[];

2465: struct _n_PetscSubcomm {
2466:   MPI_Comm         parent;           /* parent communicator */
2467:   MPI_Comm         dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2468:   MPI_Comm         child;            /* the sub-communicator */
2469:   PetscMPIInt      n;                /* num of subcommunicators under the parent communicator */
2470:   PetscMPIInt      color;            /* color of processors belong to this communicator */
2471:   PetscMPIInt      *subsize;         /* size of subcommunicator[color] */
2472:   PetscSubcommType type;
2473:   char             *subcommprefix;
2474: };

2476: PETSC_STATIC_INLINE MPI_Comm PetscSubcommParent(PetscSubcomm scomm) {return scomm->parent;}
2477: PETSC_STATIC_INLINE MPI_Comm PetscSubcommChild(PetscSubcomm scomm) {return scomm->child;}
2478: PETSC_STATIC_INLINE MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm) {return scomm->dupparent;}
2479: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2480: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2481: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2482: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2483: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2484: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2485: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2486: PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm,const char[]);

2488: PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt,PetscHeap*);
2489: PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap,PetscInt,PetscInt);
2490: PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap,PetscInt*,PetscInt*);
2491: PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap,PetscInt*,PetscInt*);
2492: PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap,PetscInt,PetscInt);
2493: PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2494: PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap*);
2495: PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap,PetscViewer);

2497: PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2498: PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm,PetscShmComm*);
2499: PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2500: PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2501: PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm,MPI_Comm*);

2503: /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2504: PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm,PetscInt,PetscOmpCtrl*);
2505: PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl,MPI_Comm*,MPI_Comm*,PetscBool*);
2506: PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl*);
2507: PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2508: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2509: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);

2511: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2512: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2513: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2514: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2515: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2516: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2517: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2518: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);


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

2526: extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2527: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted*);
2528: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted*);
2529: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted,const char*,const char*,PetscBool*);

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

2535: PETSC_EXTERN PetscSegBuffer PetscCitationsList;

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

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

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

2546:    Level: intermediate

2548:    Not available from Fortran

2550:      Options Database:
2551: .     -citations [filename]   - print out the bibtex entries for the given computation
2552: @*/
2553: PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2554: {
2555:   size_t         len;
2556:   char           *vstring;

2560:   if (set && *set) return(0);
2561:   PetscStrlen(cit,&len);
2562:   PetscSegBufferGet(PetscCitationsList,len,&vstring);
2563:   PetscArraycpy(vstring,cit,len);
2564:   if (set) *set = PETSC_TRUE;
2565:   return(0);
2566: }

2568: PETSC_EXTERN PetscErrorCode PetscURLShorten(const char[],char[],size_t);
2569: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm,char[],char[],size_t);
2570: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm,const char[],char[],size_t);
2571: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm,const char[],const char []);

2573: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm,char[],char[],size_t);
2574: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm,const char[],char[],char[],size_t);

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

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

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


2585: #if defined(PETSC_USE_DEBUG)
2586: /*
2587:    Verify that all processes in the communicator have called this from the same line of code
2588:  */
2589: PETSC_EXTERN PetscErrorCode PetscAllreduceBarrierCheck(MPI_Comm,PetscMPIInt,int,const char*,const char *);

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

2596:    Collective

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

2602:    Input Parameters:
2603: +  indata - pointer to the input data to be reduced
2604: .  count - the number of MPI data items in indata and outdata
2605: .  datatype - the MPI datatype, for example MPI_INT
2606: .  op - the MPI operation, for example MPI_SUM
2607: -  comm - the MPI communicator on which the operation occurs

2609:    Output Parameter:
2610: .  outdata - the reduced values

2612:    Notes:
2613:    In optimized mode this directly calls MPI_Allreduce()

2615:    Level: developer

2617: .seealso: MPI_Allreduce()
2618: M*/
2619: #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))
2620: #else
2621: #define MPIU_Allreduce(a,b,c,d,e,fcomm) MPI_Allreduce(a,b,c,d,e,fcomm)
2622: #endif

2624: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
2625: PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint,PetscMPIInt,MPI_Info,MPI_Comm,void*,MPI_Win*);
2626: PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win,PetscMPIInt,MPI_Aint*,PetscMPIInt*,void*);
2627: #endif

2629: /*
2630:     Returned from PETSc functions that are called from MPI, such as related to attributes
2631: */
2632: PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CLASS;
2633: PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CODE;

2635: #endif