Actual source code: petscsys.h

petsc-master 2020-07-10
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: /*MC
138:   PetscDefined - determine whether a boolean macro is defined

140:   Notes:
141:   Typical usage is within normal code,

143: $   if (PetscDefined(USE_DEBUG)) { ... }

145:   but can also be used in the preprocessor,

147: $   #if PetscDefined(USE_DEBUG)
148: $     ...
149: $   #else

151:   Either way evaluates true if PETSC_USE_DEBUG is defined (merely defined or defined to 1) or undefined.  This macro
152:   should not be used if its argument may be defined to a non-empty value other than 1.

154:   Developer Notes:
155:   Getting something that works in C and CPP for an arg that may or may not be defined is tricky.  Here, if we have
156:   "#define PETSC_HAVE_BOOGER 1" we match on the placeholder define, insert the "0," for arg1 and generate the triplet
157:   (0, 1, 0).  Then the last step cherry picks the 2nd arg (a one).  When PETSC_HAVE_BOOGER is not defined, we generate
158:   a (... 1, 0) pair, and when the last step cherry picks the 2nd arg, we get a zero.

160:   Our extra expansion via PetscDefined__take_second_expand() is needed with MSVC, which has a nonconforming
161:   implementation of variadic macros.

163:   Level: developer
164: M*/
165: #if !defined(PETSC_SKIP_VARIADIC_MACROS)
166: #  define PetscDefined_arg_1    shift,
167: #  define PetscDefined_arg_     shift,
168: #  define PetscDefined__take_second_expanded(ignored, val, ...) val
169: #  define PetscDefined__take_second_expand(args) PetscDefined__take_second_expanded args
170: #  define PetscDefined__take_second(...) PetscDefined__take_second_expand((__VA_ARGS__))
171: #  define PetscDefined___(arg1_or_junk) PetscDefined__take_second(arg1_or_junk 1, 0, at_)
172: #  define PetscDefined__(value) PetscDefined___(PetscDefined_arg_ ## value)
173: #  define PetscDefined_(d)      PetscDefined__(d)
174: #  define PetscDefined(d)       PetscDefined_(PETSC_ ## d)
175: #endif

177: /*
178:    Perform various sanity checks that the correct mpi.h is being included at compile time.
179:    This usually happens because
180:       * either an unexpected mpi.h is in the default compiler path (i.e. in /usr/include) or
181:       * an extra include path -I/something (which contains the unexpected mpi.h) is being passed to the compiler
182: */
183: #if defined(PETSC_HAVE_MPIUNI)
184: #  if !defined(MPIUNI_H)
185: #    error "PETSc was configured with --with-mpi=0 but now appears to be compiling using a different mpi.h"
186: #  endif
187: #elif defined(PETSC_HAVE_I_MPI_NUMVERSION)
188: #  if !defined(I_MPI_NUMVERSION)
189: #    error "PETSc was configured with I_MPI but now appears to be compiling using a non-I_MPI mpi.h"
190: #  elif I_MPI_NUMVERSION != PETSC_HAVE_I_MPI_NUMVERSION
191: #    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"
192: #  endif
193: #elif defined(PETSC_HAVE_MVAPICH2_NUMVERSION)
194: #  if !defined(MVAPICH2_NUMVERSION)
195: #    error "PETSc was configured with MVAPICH2 but now appears to be compiling using a non-MVAPICH2 mpi.h"
196: #  elif MVAPICH2_NUMVERSION != PETSC_HAVE_MVAPICH2_NUMVERSION
197: #    error "PETSc was configured with one MVAPICH2 mpi.h version but now appears to be compiling using a different MVAPICH2 mpi.h version"
198: #  endif
199: #elif defined(PETSC_HAVE_MPICH_NUMVERSION)
200: #  if !defined(MPICH_NUMVERSION) || defined(MVAPICH2_NUMVERSION) || defined(I_MPI_NUMVERSION)
201: #    error "PETSc was configured with MPICH but now appears to be compiling using a non-MPICH mpi.h"
202: #  elif (MPICH_NUMVERSION/100000 != PETSC_HAVE_MPICH_NUMVERSION/100000) || (MPICH_NUMVERSION%100000/1000 < PETSC_HAVE_MPICH_NUMVERSION%100000/1000)
203: #    error "PETSc was configured with one MPICH mpi.h version but now appears to be compiling using a different MPICH mpi.h version"
204: #  endif
205: #elif defined(PETSC_HAVE_OMPI_MAJOR_VERSION)
206: #  if !defined(OMPI_MAJOR_VERSION)
207: #    error "PETSc was configured with OpenMPI but now appears to be compiling using a non-OpenMPI mpi.h"
208: #  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)
209: #    error "PETSc was configured with one OpenMPI mpi.h version but now appears to be compiling using a different OpenMPI mpi.h version"
210: #  endif
211: #elif defined(PETSC_HAVE_MSMPI_VERSION)
212: #  if !defined(MSMPI_VER)
213: #    error "PETSc was configured with MSMPI but now appears to be compiling using a non-MSMPI mpi.h"
214: #  elif (MSMPI_VER != PETSC_HAVE_MSMPI_VERSION)
215: #    error "PETSc was configured with one MSMPI mpi.h version but now appears to be compiling using a different MSMPI mpi.h version"
216: #  endif
217: #elif defined(OMPI_MAJOR_VERSION) || defined(MPICH_NUMVERSION) || defined(MSMPI_VER)
218: #  error "PETSc was configured with undetermined MPI - but now appears to be compiling using any of OpenMPI, MS-MPI or a MPICH variant"
219: #endif

221: /*
222:     Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler
223:     see the top of mpicxx.h in the MPICH2 distribution.
224: */
225: #include <stdio.h>

227: /* MSMPI on 32bit windows requires this yukky hack - that breaks MPI standard compliance */
228: #if !defined(MPIAPI)
229: #define MPIAPI
230: #endif

232: /*
233:    Support for Clang (>=3.2) matching type tag arguments with void* buffer types.
234:    This allows the compiler to detect cases where the MPI datatype argument passed to a MPI routine
235:    does not match the actual type of the argument being passed in
236: */
237: #if defined(__has_attribute) && defined(works_with_const_which_is_not_true)
238: #  if __has_attribute(argument_with_type_tag) && __has_attribute(pointer_with_type_tag) && __has_attribute(type_tag_for_datatype)
239: #    define PetscAttrMPIPointerWithType(bufno,typeno) __attribute__((pointer_with_type_tag(MPI,bufno,typeno)))
240: #    define PetscAttrMPITypeTag(type)                 __attribute__((type_tag_for_datatype(MPI,type)))
241: #    define PetscAttrMPITypeTagLayoutCompatible(type) __attribute__((type_tag_for_datatype(MPI,type,layout_compatible)))
242: #  endif
243: #endif
244: #if !defined(PetscAttrMPIPointerWithType)
245: #  define PetscAttrMPIPointerWithType(bufno,typeno)
246: #  define PetscAttrMPITypeTag(type)
247: #  define PetscAttrMPITypeTagLayoutCompatible(type)
248: #endif

250: PETSC_EXTERN MPI_Datatype MPIU_ENUM PetscAttrMPITypeTag(PetscEnum);
251: PETSC_EXTERN MPI_Datatype MPIU_BOOL PetscAttrMPITypeTag(PetscBool);

253: /*MC
254:    MPIU_INT - MPI datatype corresponding to PetscInt

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

259:    Level: beginner

261: .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_SCALAR, MPIU_COMPLEX
262: M*/

264: #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 */
265: #  define MPIU_INT64 MPI_INT64_T
266: #  define PetscInt64_FMT PRId64
267: #elif (PETSC_SIZEOF_LONG_LONG == 8)
268: #  define MPIU_INT64 MPI_LONG_LONG_INT
269: #  define PetscInt64_FMT "lld"
270: #elif defined(PETSC_HAVE___INT64)
271: #  define MPIU_INT64 MPI_INT64_T
272: #  define PetscInt64_FMT "ld"
273: #else
274: #  error "cannot determine PetscInt64 type"
275: #endif

277: PETSC_EXTERN MPI_Datatype MPIU_FORTRANADDR;

279: #if defined(PETSC_USE_64BIT_INDICES)
280: #  define MPIU_INT MPIU_INT64
281: #  define PetscInt_FMT PetscInt64_FMT
282: #else
283: #  define MPIU_INT MPI_INT
284: #  define PetscInt_FMT "d"
285: #endif

287: /*
288:     For the rare cases when one needs to send a size_t object with MPI
289: */
290: PETSC_EXTERN MPI_Datatype MPIU_SIZE_T;

292: /*
293:       You can use PETSC_STDOUT as a replacement of stdout. You can also change
294:     the value of PETSC_STDOUT to redirect all standard output elsewhere
295: */
296: PETSC_EXTERN FILE* PETSC_STDOUT;

298: /*
299:       You can use PETSC_STDERR as a replacement of stderr. You can also change
300:     the value of PETSC_STDERR to redirect all standard error elsewhere
301: */
302: PETSC_EXTERN FILE* PETSC_STDERR;

304: /*MC
305:     PetscUnlikely - hints the compiler that the given condition is usually FALSE

307:     Synopsis:
308:  #include <petscsys.h>
309:     PetscBool  PetscUnlikely(PetscBool  cond)

311:     Not Collective

313:     Input Parameters:
314: .   cond - condition or expression

316:     Notes:
317:     This returns the same truth value, it is only a hint to compilers that the resulting
318:     branch is unlikely.

320:     Level: advanced

322: .seealso: PetscUnlikelyDebug(), PetscLikely(), CHKERRQ
323: M*/

325: /*MC
326:     PetscLikely - hints the compiler that the given condition is usually TRUE

328:     Synopsis:
329:  #include <petscsys.h>
330:     PetscBool  PetscLikely(PetscBool  cond)

332:     Not Collective

334:     Input Parameters:
335: .   cond - condition or expression

337:     Notes:
338:     This returns the same truth value, it is only a hint to compilers that the resulting
339:     branch is likely.

341:     Level: advanced

343: .seealso: PetscUnlikely()
344: M*/
345: #if defined(PETSC_HAVE_BUILTIN_EXPECT)
346: #  define PetscUnlikely(cond)   __builtin_expect(!!(cond),0)
347: #  define PetscLikely(cond)     __builtin_expect(!!(cond),1)
348: #else
349: #  define PetscUnlikely(cond)   (cond)
350: #  define PetscLikely(cond)     (cond)
351: #endif

353: /*MC
354:     PetscUnlikelyDebug - hints the compiler that the given condition is usually FALSE, eliding the check in optimized mode

356:     Synopsis:
357:  #include <petscsys.h>
358:     PetscBool  PetscUnlikelyDebug(PetscBool  cond)

360:     Not Collective

362:     Input Parameters:
363: .   cond - condition or expression

365:     Notes:
366:     This returns the same truth value, it is only a hint to compilers that the resulting
367:     branch is unlikely.  When compiled in optimized mode, it always returns false.

369:     Level: advanced

371: .seealso: PetscUnlikely(), CHKERRQ, SETERRQ
372: M*/
373: #define PetscUnlikelyDebug(cond) (PetscDefined(USE_DEBUG) && PetscUnlikely(cond))

375: /*
376:     Declare extern C stuff after including external header files
377: */

379: PETSC_EXTERN const char *const PetscBools[];

381: /*
382:     Defines elementary mathematics functions and constants.
383: */
384:  #include <petscmath.h>

386: PETSC_EXTERN const char *const PetscCopyModes[];

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

391:    Level: beginner

393:    Note:
394:    Accepted by many PETSc functions to not set a parameter and instead use
395:           some default

397:    Fortran Notes:
398:    This macro does not exist in Fortran; you must use PETSC_NULL_INTEGER,
399:           PETSC_NULL_DOUBLE_PRECISION etc

401: .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_DETERMINE

403: M*/
404: #define PETSC_IGNORE NULL

406: /* This is deprecated */
407: #define PETSC_NULL NULL

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

413:    Level: beginner

415: .seealso: PETSC_DEFAULT, PETSC_IGNORE, PETSC_DETERMINE

417: M*/
418: #define PETSC_DECIDE -1

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

424:    Level: beginner

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

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

432: M*/
433: #define PETSC_DETERMINE PETSC_DECIDE

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

439:    Level: beginner

441:    Fortran Notes:
442:    You need to use PETSC_DEFAULT_INTEGER or PETSC_DEFAULT_REAL.

444: .seealso: PETSC_DECIDE, PETSC_IGNORE, PETSC_DETERMINE

446: M*/
447: #define PETSC_DEFAULT -2

449: /*MC
450:     PETSC_COMM_WORLD - the equivalent of the MPI_COMM_WORLD communicator which represents
451:            all the processes that PETSc knows about.

453:    Level: beginner

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

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

464: .seealso: PETSC_COMM_SELF

466: M*/
467: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;

469: /*MC
470:     PETSC_COMM_SELF - This is always MPI_COMM_SELF

472:    Level: beginner

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

477: .seealso: PETSC_COMM_WORLD

479: M*/
480: #define PETSC_COMM_SELF MPI_COMM_SELF

482: /*MC
483:     PETSC_MPI_THREAD_REQUIRED - the required threading support used if PETSc initializes
484:            MPI with MPI_Init_thread.

486:    Level: beginner

488:    Notes:
489:    By default PETSC_MPI_THREAD_REQUIRED equals MPI_THREAD_FUNNELED.

491: .seealso: PetscInitialize()

493: M*/
494: PETSC_EXTERN PetscMPIInt PETSC_MPI_THREAD_REQUIRED;

496: PETSC_EXTERN PetscBool PetscBeganMPI;
497: PETSC_EXTERN PetscBool PetscErrorHandlingInitialized;
498: PETSC_EXTERN PetscBool PetscInitializeCalled;
499: PETSC_EXTERN PetscBool PetscFinalizeCalled;
500: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;
501: PETSC_EXTERN PetscBool PetscCUDASynchronize;

503: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm));
504: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*);
505: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm*);

507: #if defined(PETSC_HAVE_CUDA)
508: PETSC_EXTERN PetscErrorCode PetscCUDAInitialize(MPI_Comm);
509: #endif

511: #if defined(PETSC_HAVE_ELEMENTAL)
512: PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
513: PETSC_EXTERN PetscErrorCode PetscElementalInitialized(PetscBool*);
514: PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
515: #endif

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

520:    Synopsis:
521:  #include <petscsys.h>
522:    PetscErrorCode PetscMalloc(size_t m,void **result)

524:    Not Collective

526:    Input Parameter:
527: .  m - number of bytes to allocate

529:    Output Parameter:
530: .  result - memory allocated

532:    Level: beginner

534:    Notes:
535:    Memory is always allocated at least double aligned

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

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

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

544: /*MC
545:    PetscRealloc - Rellocates memory

547:    Synopsis:
548:  #include <petscsys.h>
549:    PetscErrorCode PetscRealloc(size_t m,void **result)

551:    Not Collective

553:    Input Parameters:
554: +  m - number of bytes to allocate
555: -  result - previous memory

557:    Output Parameter:
558: .  result - new memory allocated

560:    Level: developer

562:    Notes:
563:    Memory is always allocated at least double aligned

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

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

570: /*MC
571:    PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment

573:    Synopsis:
574:  #include <petscsys.h>
575:    void *PetscAddrAlign(void *addr)

577:    Not Collective

579:    Input Parameters:
580: .  addr - address to align (any pointer type)

582:    Level: developer

584: .seealso: PetscMallocAlign()

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

589: /*MC
590:    PetscMalloc1 - Allocates an array of memory aligned to PETSC_MEMALIGN

592:    Synopsis:
593:  #include <petscsys.h>
594:    PetscErrorCode PetscMalloc1(size_t m1,type **r1)

596:    Not Collective

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

601:    Output Parameter:
602: .  r1 - memory allocated

604:    Note:
605:    This uses the sizeof() of the memory type requested to determine the total memory to be allocated, therefore you should not
606:          multiply the number of elements requested by the sizeof() the type. For example use
607: $  PetscInt *id;
608: $  PetscMalloc1(10,&id);
609:         not
610: $  PetscInt *id;
611: $  PetscMalloc1(10*sizeof(PetscInt),&id);

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

615:    Level: beginner

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

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

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

625:    Synopsis:
626:  #include <petscsys.h>
627:    PetscErrorCode PetscCalloc1(size_t m1,type **r1)

629:    Not Collective

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

634:    Output Parameter:
635: .  r1 - memory allocated

637:    Notes:
638:    See PetsMalloc1() for more details on usage.

640:    Level: beginner

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

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

647: /*MC
648:    PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN

650:    Synopsis:
651:  #include <petscsys.h>
652:    PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)

654:    Not Collective

656:    Input Parameter:
657: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
658: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

660:    Output Parameter:
661: +  r1 - memory allocated in first chunk
662: -  r2 - memory allocated in second chunk

664:    Level: developer

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

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

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

674:    Synopsis:
675:  #include <petscsys.h>
676:    PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)

678:    Not Collective

680:    Input Parameter:
681: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
682: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

684:    Output Parameter:
685: +  r1 - memory allocated in first chunk
686: -  r2 - memory allocated in second chunk

688:    Level: developer

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

692: M*/
693: #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))

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

698:    Synopsis:
699:  #include <petscsys.h>
700:    PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

702:    Not Collective

704:    Input Parameter:
705: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
706: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
707: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

709:    Output Parameter:
710: +  r1 - memory allocated in first chunk
711: .  r2 - memory allocated in second chunk
712: -  r3 - memory allocated in third chunk

714:    Level: developer

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

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

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

724:    Synopsis:
725:  #include <petscsys.h>
726:    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

728:    Not Collective

730:    Input Parameter:
731: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
732: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
733: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

735:    Output Parameter:
736: +  r1 - memory allocated in first chunk
737: .  r2 - memory allocated in second chunk
738: -  r3 - memory allocated in third chunk

740:    Level: developer

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

744: M*/
745: #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))

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

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

754:    Not Collective

756:    Input Parameter:
757: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
758: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
759: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
760: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

762:    Output Parameter:
763: +  r1 - memory allocated in first chunk
764: .  r2 - memory allocated in second chunk
765: .  r3 - memory allocated in third chunk
766: -  r4 - memory allocated in fourth chunk

768:    Level: developer

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

772: M*/
773: #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))

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

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

782:    Not Collective

784:    Input Parameters:
785: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
786: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
787: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
788: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

790:    Output Parameters:
791: +  r1 - memory allocated in first chunk
792: .  r2 - memory allocated in second chunk
793: .  r3 - memory allocated in third chunk
794: -  r4 - memory allocated in fourth chunk

796:    Level: developer

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

800: M*/
801: #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))

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

806:    Synopsis:
807:  #include <petscsys.h>
808:    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)

810:    Not Collective

812:    Input Parameters:
813: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
814: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
815: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
816: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
817: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

819:    Output Parameters:
820: +  r1 - memory allocated in first chunk
821: .  r2 - memory allocated in second chunk
822: .  r3 - memory allocated in third chunk
823: .  r4 - memory allocated in fourth chunk
824: -  r5 - memory allocated in fifth chunk

826:    Level: developer

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

830: M*/
831: #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))

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

836:    Synopsis:
837:  #include <petscsys.h>
838:    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)

840:    Not Collective

842:    Input Parameters:
843: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
844: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
845: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
846: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
847: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

849:    Output Parameters:
850: +  r1 - memory allocated in first chunk
851: .  r2 - memory allocated in second chunk
852: .  r3 - memory allocated in third chunk
853: .  r4 - memory allocated in fourth chunk
854: -  r5 - memory allocated in fifth chunk

856:    Level: developer

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

860: M*/
861: #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))

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

866:    Synopsis:
867:  #include <petscsys.h>
868:    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)

870:    Not Collective

872:    Input Parameters:
873: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
874: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
875: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
876: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
877: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
878: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

880:    Output Parameteasr:
881: +  r1 - memory allocated in first chunk
882: .  r2 - memory allocated in second chunk
883: .  r3 - memory allocated in third chunk
884: .  r4 - memory allocated in fourth chunk
885: .  r5 - memory allocated in fifth chunk
886: -  r6 - memory allocated in sixth chunk

888:    Level: developer

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

892: M*/
893: #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))

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

898:    Synopsis:
899:  #include <petscsys.h>
900:    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)

902:    Not Collective

904:    Input Parameters:
905: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
906: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
907: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
908: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
909: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
910: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

912:    Output Parameters:
913: +  r1 - memory allocated in first chunk
914: .  r2 - memory allocated in second chunk
915: .  r3 - memory allocated in third chunk
916: .  r4 - memory allocated in fourth chunk
917: .  r5 - memory allocated in fifth chunk
918: -  r6 - memory allocated in sixth chunk

920:    Level: developer

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

924: M*/
925: #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))

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

930:    Synopsis:
931:  #include <petscsys.h>
932:    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)

934:    Not Collective

936:    Input Parameters:
937: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
938: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
939: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
940: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
941: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
942: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
943: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

945:    Output Parameters:
946: +  r1 - memory allocated in first chunk
947: .  r2 - memory allocated in second chunk
948: .  r3 - memory allocated in third chunk
949: .  r4 - memory allocated in fourth chunk
950: .  r5 - memory allocated in fifth chunk
951: .  r6 - memory allocated in sixth chunk
952: -  r7 - memory allocated in seventh chunk

954:    Level: developer

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

958: M*/
959: #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))

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

964:    Synopsis:
965:  #include <petscsys.h>
966:    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)

968:    Not Collective

970:    Input Parameters:
971: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
972: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
973: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
974: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
975: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
976: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
977: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

979:    Output Parameters:
980: +  r1 - memory allocated in first chunk
981: .  r2 - memory allocated in second chunk
982: .  r3 - memory allocated in third chunk
983: .  r4 - memory allocated in fourth chunk
984: .  r5 - memory allocated in fifth chunk
985: .  r6 - memory allocated in sixth chunk
986: -  r7 - memory allocated in seventh chunk

988:    Level: developer

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

992: M*/
993: #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))

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

998:    Synopsis:
999:  #include <petscsys.h>
1000:    PetscErrorCode PetscNew(type **result)

1002:    Not Collective

1004:    Output Parameter:
1005: .  result - memory allocated, sized to match pointer type

1007:    Level: beginner

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

1011: M*/
1012: #define PetscNew(b)      PetscCalloc1(1,(b))

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

1018:    Synopsis:
1019:  #include <petscsys.h>
1020:    PetscErrorCode PetscNewLog(PetscObject obj,type **result)

1022:    Not Collective

1024:    Input Parameter:
1025: .  obj - object memory is logged to

1027:    Output Parameter:
1028: .  result - memory allocated, sized to match pointer type

1030:    Level: developer

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

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

1037: /*MC
1038:    PetscFree - Frees memory

1040:    Synopsis:
1041:  #include <petscsys.h>
1042:    PetscErrorCode PetscFree(void *memory)

1044:    Not Collective

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

1049:    Level: beginner

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

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

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

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

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

1064:    Synopsis:
1065:  #include <petscsys.h>
1066:    PetscErrorCode PetscFree2(void *memory1,void *memory2)

1068:    Not Collective

1070:    Input Parameters:
1071: +   memory1 - memory to free
1072: -   memory2 - 2nd memory to free

1074:    Level: developer

1076:    Notes:
1077:     Memory must have been obtained with PetscMalloc2()

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

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

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

1087:    Synopsis:
1088:  #include <petscsys.h>
1089:    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)

1091:    Not Collective

1093:    Input Parameters:
1094: +   memory1 - memory to free
1095: .   memory2 - 2nd memory to free
1096: -   memory3 - 3rd memory to free

1098:    Level: developer

1100:    Notes:
1101:     Memory must have been obtained with PetscMalloc3()

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

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

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

1111:    Synopsis:
1112:  #include <petscsys.h>
1113:    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)

1115:    Not Collective

1117:    Input Parameters:
1118: +   m1 - memory to free
1119: .   m2 - 2nd memory to free
1120: .   m3 - 3rd memory to free
1121: -   m4 - 4th memory to free

1123:    Level: developer

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

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

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

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

1136:    Synopsis:
1137:  #include <petscsys.h>
1138:    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)

1140:    Not Collective

1142:    Input Parameters:
1143: +   m1 - memory to free
1144: .   m2 - 2nd memory to free
1145: .   m3 - 3rd memory to free
1146: .   m4 - 4th memory to free
1147: -   m5 - 5th memory to free

1149:    Level: developer

1151:    Notes:
1152:     Memory must have been obtained with PetscMalloc5()

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

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

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

1162:    Synopsis:
1163:  #include <petscsys.h>
1164:    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)

1166:    Not Collective

1168:    Input Parameters:
1169: +   m1 - memory to free
1170: .   m2 - 2nd memory to free
1171: .   m3 - 3rd memory to free
1172: .   m4 - 4th memory to free
1173: .   m5 - 5th memory to free
1174: -   m6 - 6th memory to free


1177:    Level: developer

1179:    Notes:
1180:     Memory must have been obtained with PetscMalloc6()

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

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

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

1190:    Synopsis:
1191:  #include <petscsys.h>
1192:    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)

1194:    Not Collective

1196:    Input Parameters:
1197: +   m1 - memory to free
1198: .   m2 - 2nd memory to free
1199: .   m3 - 3rd memory to free
1200: .   m4 - 4th memory to free
1201: .   m5 - 5th memory to free
1202: .   m6 - 6th memory to free
1203: -   m7 - 7th memory to free


1206:    Level: developer

1208:    Notes:
1209:     Memory must have been obtained with PetscMalloc7()

1211: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1212:           PetscMalloc7()

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

1217: PETSC_EXTERN PetscErrorCode PetscMallocA(int,PetscBool,int,const char *,const char *,size_t,void *,...);
1218: PETSC_EXTERN PetscErrorCode PetscFreeA(int,int,const char *,const char *,void *,...);
1219: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,PetscBool,int,const char[],const char[],void**);
1220: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1221: PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t,int,const char[],const char[],void**);
1222: PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1223: 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 **));
1224: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);

1226: /*
1227:   Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1228: */
1229: PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1230: PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);
1231: #if defined(PETSC_HAVE_CUDA)
1232: PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void);
1233: PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void);
1234: #endif

1236: #define MPIU_PETSCLOGDOUBLE  MPI_DOUBLE
1237: #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION

1239: /*
1240:    Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1241: */
1242: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1243: PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1244: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1245: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1246: PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1247: PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int,PetscLogDouble*);
1248: PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool,PetscBool);
1249: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*,PetscBool*,PetscBool*);
1250: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1251: PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1252: PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool*);

1254: PETSC_EXTERN const char *const PetscDataTypes[];
1255: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1256: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1257: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1258: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);

1260: /*
1261:     Basic memory and string operations. These are usually simple wrappers
1262:    around the basic Unix system calls, but a few of them have additional
1263:    functionality and/or error checking.
1264: */
1265: PETSC_EXTERN PetscErrorCode PetscBitMemcpy(void*,PetscInt,const void*,PetscInt,PetscInt,PetscDataType);
1266: PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1267: PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1268: PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***);
1269: PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1270: PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool  *);
1271: PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool  *);
1272: PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1273: PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1274: PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1275: PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1276: PETSC_EXTERN PetscErrorCode PetscStrlcat(char[],const char[],size_t);
1277: PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1278: PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1279: PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1280: PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1281: PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1282: PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1283: PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1284: PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1285: PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*);
1286: PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1287: PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1288: PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1289: PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1290: PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt,const char *const*,char***);
1291: PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt,char***);
1292: PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);

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

1296: PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1297: PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1298: PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);

1300: PETSC_EXTERN PetscErrorCode PetscStrInList(const char[],const char[],char,PetscBool*);
1301: PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1302: PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);

1304: /*
1305:    These are MPI operations for MPI_Allreduce() etc
1306: */
1307: PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1308: #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1309: PETSC_EXTERN MPI_Op MPIU_SUM;
1310: #else
1311: #define MPIU_SUM MPI_SUM
1312: #endif
1313: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1314: PETSC_EXTERN MPI_Op MPIU_MAX;
1315: PETSC_EXTERN MPI_Op MPIU_MIN;
1316: #else
1317: #define MPIU_MAX MPI_MAX
1318: #define MPIU_MIN MPI_MIN
1319: #endif
1320: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);

1322: PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1323: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);

1325: PETSC_EXTERN const char *const PetscFileModes[];

1327: /*
1328:     Defines PETSc error handling.
1329: */
1330:  #include <petscerror.h>

1332: #define PETSC_SMALLEST_CLASSID  1211211
1333: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1334: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1335: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);
1336: PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject,PetscObjectId*);
1337: PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject,PetscObjectId,PetscBool*);


1340: /*
1341:    Routines that get memory usage information from the OS
1342: */
1343: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1344: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1345: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1346: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);

1348: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);

1350: /*
1351:    Initialization of PETSc
1352: */
1353: PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1354: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1355: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1356: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1357: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1358: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1359: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1360: PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1361: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1362: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);

1364: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1365: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);

1367: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1368: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1369: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1370: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);

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


1375: /*
1376:      These are so that in extern C code we can caste function pointers to non-extern C
1377:    function pointers. Since the regular C++ code expects its function pointers to be C++
1378: */
1379: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1380: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1381: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);

1383: /*
1384:     Functions that can act on any PETSc object.
1385: */
1386: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1387: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1388: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1389: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1390: PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1391: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1392: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1393: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1394: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1395: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1396: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1397: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1398: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1399: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1400: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt*);
1401: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1402: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1403: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject*);
1404: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void));
1405: #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d))
1406: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1407: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1408: PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1409: PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject,PetscObject);
1410: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt*);

1412:  #include <petscviewertypes.h>
1413:  #include <petscoptions.h>

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

1417: PETSC_EXTERN PetscErrorCode PetscMemoryShowUsage(PetscViewer,const char[]);
1418: PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer,const char[]);
1419: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1420: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1421: #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1422: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1423: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1424: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1425: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1426: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1427: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1428: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1429: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1430: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,PetscObject,const char[]);
1431: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1432: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1433: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject,const char[],PetscBool *);
1434: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1435: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1436: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1437: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);

1439: #if defined(PETSC_HAVE_SAWS)
1440: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1441: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1442: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1443: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1444: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1445: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1446: PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1447: PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1448: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1449: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);

1451: #else
1452: #define PetscSAWsBlock()                        0
1453: #define PetscObjectSAWsViewOff(obj)             0
1454: #define PetscObjectSAWsSetBlock(obj,flg)        0
1455: #define PetscObjectSAWsBlock(obj)               0
1456: #define PetscObjectSAWsGrantAccess(obj)         0
1457: #define PetscObjectSAWsTakeAccess(obj)          0
1458: #define PetscStackViewSAWs()                    0
1459: #define PetscStackSAWsViewOff()                 0
1460: #define PetscStackSAWsTakeAccess()
1461: #define PetscStackSAWsGrantAccess()

1463: #endif

1465: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1466: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1467: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);

1469: #if defined(PETSC_USE_DEBUG)
1470: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1471: #endif
1472: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);

1474: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1475: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1476: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1477: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1478: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1479: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);

1481: /*
1482:     Dynamic library lists. Lists of names of routines in objects or in dynamic
1483:   link libraries that will be loaded as needed.
1484: */

1486: #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1487: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1488: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1489: #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1490: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1491: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[],const char[]);
1492: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1493: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1494: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);

1496: PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1497: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1498: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1499: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1500: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1501: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1502: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1503: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);

1505: /*
1506:      Useful utility routines
1507: */
1508: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1509: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1510: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm,PetscInt*,PetscInt*);
1511: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1512: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1513: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1514: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);
1515: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm,PetscInt[],PetscInt[]);
1516: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm,PetscReal[],PetscReal[]);

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

1521:     Notes:
1522:     This is useful in cases like
1523: $     int        *a;
1524: $     PetscBool  flag = PetscNot(a)
1525:      where !a would not return a PetscBool because we cannot provide a cast from int to PetscBool in C.

1527:     Level: beginner

1529:     .seealso : PetscBool, PETSC_TRUE, PETSC_FALSE
1530: M*/
1531: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)

1533: /*MC
1534:     PetscHelpPrintf - Prints help messages.

1536:    Synopsis:
1537:  #include <petscsys.h>
1538:      PetscErrorCode (*PetscHelpPrintf)(const char format[],...);

1540:     Not Collective

1542:     Input Parameters:
1543: .   format - the usual printf() format string

1545:    Level: developer

1547:     Fortran Note:
1548:     This routine is not supported in Fortran.


1551: .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1552: M*/
1553: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);

1555: /*
1556:      Defines PETSc profiling.
1557: */
1558:  #include <petsclog.h>

1560: /*
1561:       Simple PETSc parallel IO for ASCII printing
1562: */
1563: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1564: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1565: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1566: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1567: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1568: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1569: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);
1570: PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[],size_t,const char*,PetscInt,const PetscReal[]);

1572: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1573: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1574: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);

1576: PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char*,size_t*);
1577: PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char*,char *);

1579: #if defined(PETSC_HAVE_POPEN)
1580: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1581: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*);
1582: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1583: #endif

1585: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1586: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1587: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*);
1588: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1589: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1590: PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1591: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);

1593: PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1594: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void**);
1595: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void*);
1596: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1597: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer*);
1598: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer,PetscErrorCode (*)(void*));
1599: PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void*);

1601: /*
1602:    For use in debuggers
1603: */
1604: PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1605: PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1606: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1607: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1608: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);

1610: #include <stddef.h>
1611: #include <string.h>             /* for memcpy, memset */
1612: #include <stdlib.h>

1614: #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1615: #include <xmmintrin.h>
1616: #endif

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

1623:    Not Collective

1625:    Input Parameters:
1626: +  b - pointer to initial memory space
1627: .  a - pointer to copy space
1628: -  n - length (in bytes) of space to copy

1630:    Level: intermediate

1632:    Note:
1633:    PetscArraymove() is preferred
1634:    This routine is analogous to memmove().

1636:    Developers Note: This is inlined for performance

1638: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscStrallocpy(),
1639:           PetscArraymove()
1640: @*/
1641: PETSC_STATIC_INLINE PetscErrorCode PetscMemmove(void *a,const void *b,size_t n)
1642: {
1644:   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to null pointer");
1645:   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1646: #if !defined(PETSC_HAVE_MEMMOVE)
1647:   if (a < b) {
1648:     if (a <= b - n) memcpy(a,b,n);
1649:     else {
1650:       memcpy(a,b,(int)(b - a));
1651:       PetscMemmove(b,b + (int)(b - a),n - (int)(b - a));
1652:     }
1653:   } else {
1654:     if (b <= a - n) memcpy(a,b,n);
1655:     else {
1656:       memcpy(b + n,b + (n - (int)(a - b)),(int)(a - b));
1657:       PetscMemmove(a,b,n - (int)(a - b));
1658:     }
1659:   }
1660: #else
1661:   memmove((char*)(a),(char*)(b),n);
1662: #endif
1663:   return(0);
1664: }

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

1671:    Not Collective

1673:    Input Parameters:
1674: +  b - pointer to initial memory space
1675: -  n - length (in bytes) of space to copy

1677:    Output Parameter:
1678: .  a - pointer to copy space

1680:    Level: intermediate

1682:    Compile Option:
1683:     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1684:                                   for memory copies on double precision values.
1685:     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1686:                                   for memory copies on double precision values.
1687:     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1688:                                   for memory copies on double precision values.

1690:    Note:
1691:    Prefer PetscArraycpy()
1692:    This routine is analogous to memcpy().
1693:    Not available from Fortran

1695:    Developer Note: this is inlined for fastest performance

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

1699: @*/
1700: PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n)
1701: {
1702: #if defined(PETSC_USE_DEBUG)
1703:   size_t al = (size_t) a,bl = (size_t) b;
1704:   size_t nl = (size_t) n;
1706:   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1707:   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1708: #else
1710: #endif
1711:   if (a != b && n > 0) {
1712: #if defined(PETSC_USE_DEBUG)
1713:     if ((al > bl && (al - bl) < nl) || (bl - al) < nl)  SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1714:               or make sure your copy regions and lengths are correct. \n\
1715:               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1716: #endif
1717: #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1718:    if (!(a % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1719:       size_t len = n/sizeof(PetscScalar);
1720: #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1721:       PetscBLASInt   one = 1,blen;
1723:       PetscBLASIntCast(len,&blen);
1724:       PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one));
1725: #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1726:       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1727: #else
1728:       size_t      i;
1729:       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1730:       for (i=0; i<len; i++) y[i] = x[i];
1731: #endif
1732:     } else {
1733:       memcpy((char*)(a),(char*)(b),n);
1734:     }
1735: #else
1736:     memcpy((char*)(a),(char*)(b),n);
1737: #endif
1738:   }
1739:   return(0);
1740: }

1742: /*@C
1743:    PetscMemzero - Zeros the specified memory.

1745:    Not Collective

1747:    Input Parameters:
1748: +  a - pointer to beginning memory location
1749: -  n - length (in bytes) of memory to initialize

1751:    Level: intermediate

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

1757:    Not available from Fortran
1758:    Prefer PetscArrayzero()

1760:    Developer Note: this is inlined for fastest performance

1762: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy()
1763: @*/
1764: PETSC_STATIC_INLINE PetscErrorCode  PetscMemzero(void *a,size_t n)
1765: {
1766:   if (n > 0) {
1767: #if defined(PETSC_USE_DEBUG)
1768:     if (!a) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer with %zu bytes",n);
1769: #endif
1770: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1771:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1772:       size_t      i,len = n/sizeof(PetscScalar);
1773:       PetscScalar *x = (PetscScalar*)a;
1774:       for (i=0; i<len; i++) x[i] = 0.0;
1775:     } else {
1776: #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1777:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1778:       PetscInt len = n/sizeof(PetscScalar);
1779:       fortranzero_(&len,(PetscScalar*)a);
1780:     } else {
1781: #endif
1782: #if defined(PETSC_PREFER_BZERO)
1783:       bzero((char *)a,n);
1784: #else
1785:       memset((char*)a,0,n);
1786: #endif
1787: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1788:     }
1789: #endif
1790:   }
1791:   return 0;
1792: }

1794: /*MC
1795:    PetscArraycmp - Compares two arrays in memory.

1797:    Synopsis:
1798:  #include <petscsys.h>
1799:     PetscErrorCode PetscArraycmp(const anytype *str1,const anytype *str2,size_t cnt,PetscBool *e)

1801:    Not Collective

1803:    Input Parameters:
1804: +  str1 - First array
1805: .  str2 - Second array
1806: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1808:    Output Parameters:
1809: .   e - PETSC_TRUE if equal else PETSC_FALSE.

1811:    Level: intermediate

1813:    Note:
1814:    This routine is a preferred replacement to PetscMemcmp()
1815:    The arrays must be of the same type

1817: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy(),
1818:           PetscArraymove()
1819: M*/
1820: #define  PetscArraycmp(str1,str2,cnt,e) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcmp(str1,str2,(size_t)(cnt)*sizeof(*(str1)),e))

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

1826:    Synopsis:
1827:  #include <petscsys.h>
1828:     PetscErrorCode PetscArraymove(anytype *str1,const anytype *str2,size_t cnt)

1830:    Not Collective

1832:    Input Parameters:
1833: +  str1 - First array
1834: .  str2 - Second array
1835: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1837:    Level: intermediate

1839:    Note:
1840:    This routine is a preferred replacement to PetscMemmove()
1841:    The arrays must be of the same type

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

1847: /*MC
1848:    PetscArraycpy - Copies from one array in memory to another

1850:    Synopsis:
1851:  #include <petscsys.h>
1852:     PetscErrorCode PetscArraycpy(anytype *str1,const anytype *str2,size_t cnt)

1854:    Not Collective

1856:    Input Parameters:
1857: +  str1 - First array (destination)
1858: .  str2 - Second array (source)
1859: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1861:    Level: intermediate

1863:    Note:
1864:    This routine is a preferred replacement to PetscMemcpy()
1865:    The arrays must be of the same type

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

1871: /*MC
1872:    PetscArrayzero - Zeros an array in memory.

1874:    Synopsis:
1875:  #include <petscsys.h>
1876:     PetscErrorCode PetscArrayzero(anytype *str1,size_t cnt)

1878:    Not Collective

1880:    Input Parameters:
1881: +  str1 - array
1882: -  cnt  - Count of the array, not in bytes, but number of entries in the array

1884:    Level: intermediate

1886:    Note:
1887:    This routine is a preferred replacement to PetscMemzero()

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

1893: /*MC
1894:    PetscPrefetchBlock - Prefetches a block of memory

1896:    Synopsis:
1897:  #include <petscsys.h>
1898:     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)

1900:    Not Collective

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

1908:    Level: developer

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

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

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

1923: M*/
1924: #define PetscPrefetchBlock(a,n,rw,t) do {                               \
1925:     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
1926:     for ( ; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
1927:   } while (0)

1929: /*
1930:       Determine if some of the kernel computation routines use
1931:    Fortran (rather than C) for the numerical calculations. On some machines
1932:    and compilers (like complex numbers) the Fortran version of the routines
1933:    is faster than the C/C++ versions. The flag --with-fortran-kernels
1934:    should be used with ./configure to turn these on.
1935: */
1936: #if defined(PETSC_USE_FORTRAN_KERNELS)

1938: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1939: #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1940: #endif

1942: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
1943: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
1944: #endif

1946: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1947: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
1948: #endif

1950: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1951: #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
1952: #endif

1954: #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
1955: #define PETSC_USE_FORTRAN_KERNEL_NORM
1956: #endif

1958: #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
1959: #define PETSC_USE_FORTRAN_KERNEL_MAXPY
1960: #endif

1962: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1963: #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
1964: #endif

1966: #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1967: #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
1968: #endif

1970: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
1971: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
1972: #endif

1974: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1975: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
1976: #endif

1978: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
1979: #define PETSC_USE_FORTRAN_KERNEL_MDOT
1980: #endif

1982: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
1983: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
1984: #endif

1986: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
1987: #define PETSC_USE_FORTRAN_KERNEL_AYPX
1988: #endif

1990: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
1991: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
1992: #endif

1994: #endif

1996: /*
1997:     Macros for indicating code that should be compiled with a C interface,
1998:    rather than a C++ interface. Any routines that are dynamically loaded
1999:    (such as the PCCreate_XXX() routines) must be wrapped so that the name
2000:    mangler does not change the functions symbol name. This just hides the
2001:    ugly extern "C" {} wrappers.
2002: */
2003: #if defined(__cplusplus)
2004: #  define EXTERN_C_BEGIN extern "C" {
2005: #  define EXTERN_C_END }
2006: #else
2007: #  define EXTERN_C_BEGIN
2008: #  define EXTERN_C_END
2009: #endif

2011: /* --------------------------------------------------------------------*/

2013: /*MC
2014:     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
2015:         communication

2017:    Level: beginner

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

2021: .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
2022: M*/

2024: #if defined(PETSC_HAVE_MPIIO)
2025: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2026: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2027: PETSC_EXTERN PetscErrorCode MPIU_File_write_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2028: PETSC_EXTERN PetscErrorCode MPIU_File_read_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2029: PETSC_EXTERN PetscErrorCode MPIU_File_write_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2030: PETSC_EXTERN PetscErrorCode MPIU_File_read_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2031: #endif

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

2035: /* Limit MPI to 32-bits */
2036: #define PETSC_MPI_INT_MAX  2147483647
2037: #define PETSC_MPI_INT_MIN -2147483647
2038: /* Limit BLAS to 32-bits */
2039: #define PETSC_BLAS_INT_MAX  2147483647
2040: #define PETSC_BLAS_INT_MIN -2147483647

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

2046:    Not Collective

2048:    Input Parameter:
2049: .     a - the PetscInt64 value

2051:    Output Parameter:
2052: .     b - the resulting PetscInt value

2054:    Level: advanced

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

2058:    Not available from Fortran

2060: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscBLASIntCast()
2061: @*/
2062: PETSC_STATIC_INLINE PetscErrorCode PetscIntCast(PetscInt64 a,PetscInt *b)
2063: {
2065: #if !defined(PETSC_USE_64BIT_INDICES)
2066:   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");
2067: #endif
2068:   *b =  (PetscInt)(a);
2069:   return(0);
2070: }

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

2076:    Not Collective

2078:    Input Parameter:
2079: .     a - the PetscInt value

2081:    Output Parameter:
2082: .     b - the resulting PetscBLASInt value

2084:    Level: advanced

2086:    Not available from Fortran

2088: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscIntCast()
2089: @*/
2090: PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
2091: {
2093: #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
2094:   *b = 0;
2095:   if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK");
2096: #endif
2097:   *b =  (PetscBLASInt)(a);
2098:   return(0);
2099: }

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

2105:    Not Collective

2107:    Input Parameter:
2108: .     a - the PetscInt value

2110:    Output Parameter:
2111: .     b - the resulting PetscMPIInt value

2113:    Level: advanced

2115:    Not available from Fortran

2117: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntCast()
2118: @*/
2119: PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
2120: {
2122: #if defined(PETSC_USE_64BIT_INDICES)
2123:   *b = 0;
2124:   if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for MPI");
2125: #endif
2126:   *b =  (PetscMPIInt)(a);
2127:   return(0);
2128: }

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

2132: /*@C

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

2136:    Not Collective

2138:    Input Parameter:
2139: +     a - the PetscReal value
2140: -     b - the second value

2142:    Returns:
2143:       the result as a PetscInt value

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

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

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

2154:    Not available from Fortran

2156:    Level: advanced

2158: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2159: @*/
2160: PETSC_STATIC_INLINE PetscInt PetscRealIntMultTruncate(PetscReal a,PetscInt b)
2161: {
2162:   PetscInt64 r;

2164:   r  =  (PetscInt64) (a*(PetscReal)b);
2165:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2166:   return (PetscInt) r;
2167: }

2169: /*@C

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

2173:    Not Collective

2175:    Input Parameter:
2176: +     a - the PetscInt value
2177: -     b - the second value

2179:    Returns:
2180:       the result as a PetscInt value

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

2186:    Not available from Fortran

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

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

2192:    Level: advanced

2194: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2195: @*/
2196: PETSC_STATIC_INLINE PetscInt PetscIntMultTruncate(PetscInt a,PetscInt b)
2197: {
2198:   PetscInt64 r;

2200:   r  =  PetscInt64Mult(a,b);
2201:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2202:   return (PetscInt) r;
2203: }

2205: /*@C

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

2209:    Not Collective

2211:    Input Parameter:
2212: +     a - the PetscInt value
2213: -     b - the second value

2215:    Returns:
2216:      the result as a PetscInt value

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

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

2224:    Not available from Fortran

2226:    Level: advanced

2228: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2229: @*/
2230: PETSC_STATIC_INLINE PetscInt PetscIntSumTruncate(PetscInt a,PetscInt b)
2231: {
2232:   PetscInt64 r;

2234:   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2235:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2236:   return (PetscInt) r;
2237: }

2239: /*@C

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

2243:    Not Collective

2245:    Input Parameter:
2246: +     a - the PetscInt value
2247: -     b - the second value

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

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

2255:    Not available from Fortran

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

2259:    Level: advanced

2261: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64(), PetscIntSumError()
2262: @*/
2263: PETSC_STATIC_INLINE PetscErrorCode PetscIntMultError(PetscInt a,PetscInt b,PetscInt *result)
2264: {
2265:   PetscInt64 r;

2268:   r  =  PetscInt64Mult(a,b);
2269: #if !defined(PETSC_USE_64BIT_INDICES)
2270:   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);
2271: #endif
2272:   if (result) *result = (PetscInt) r;
2273:   return(0);
2274: }

2276: /*@C

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

2280:    Not Collective

2282:    Input Parameter:
2283: +     a - the PetscInt value
2284: -     b - the second value

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

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

2292:    Not available from Fortran

2294:    Level: advanced

2296: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2297: @*/
2298: PETSC_STATIC_INLINE PetscErrorCode PetscIntSumError(PetscInt a,PetscInt b,PetscInt *result)
2299: {
2300:   PetscInt64 r;

2303:   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2304: #if !defined(PETSC_USE_64BIT_INDICES)
2305:   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);
2306: #endif
2307:   if (result) *result = (PetscInt) r;
2308:   return(0);
2309: }

2311: /*
2312:      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2313: */
2314: #if defined(hz)
2315: #  undef hz
2316: #endif

2318: #include <limits.h>

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

2322: #define PETSC_BITS_PER_BYTE CHAR_BIT

2324: #if defined(PETSC_HAVE_SYS_TYPES_H)
2325: #  include <sys/types.h>
2326: #endif

2328: /*MC

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


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

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

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

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

2343:     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),
2344:     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
2345:     version number (x.y.z).

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

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

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

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

2355:     Level: intermediate

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

2359: M*/

2361: /*MC

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

2366: $
2367: $#include "petsc/finclude/petscXXX.h"
2368: $         use petscXXX

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

2372: $    XXX variablename
2373: $    type(tXXX) variablename

2375:     For example,

2377: $#include "petsc/finclude/petscvec.h"
2378: $         use petscvec
2379: $
2380: $    Vec b
2381: $    type(tVec) x

2383:     Level: beginner

2385: M*/

2387: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2388: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2389: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2390: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2391: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2392: PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2393: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2394: PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt*,PetscInt*,PetscInt*,PetscInt*);

2396: PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt,const PetscInt[],PetscBool*);
2397: PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt,const PetscMPIInt[],PetscBool*);
2398: PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt,const PetscReal[],PetscBool*);
2399: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2400: PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt,PetscInt[]);
2401: PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt*,PetscInt[]);
2402: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2404: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2405: PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt*);
2406: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2407: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2408: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2409: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt[],PetscInt[],PetscInt[]);
2410: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2411: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2412: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2413: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt,PetscMPIInt[],PetscInt[]);
2414: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2415: PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt,PetscInt[],void*,size_t,void*);
2416: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2417: PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt,PetscReal[],PetscInt[]);
2418: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2419: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2420: PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal,PetscInt,const PetscReal[],PetscReal,PetscInt*);
2421: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2422: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2423: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2424: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt[],const PetscInt[],PetscInt,const PetscInt[],const PetscInt[],PetscInt*,PetscInt**,PetscInt**);
2425: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt*,PetscInt**);
2426: PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt,const PetscMPIInt[],PetscInt,const PetscMPIInt[],PetscInt*,PetscMPIInt**);

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

2430: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2431: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);

2433: /*J
2434:     PetscRandomType - String with the name of a PETSc randomizer

2436:    Level: beginner

2438:    Notes:
2439:    To use SPRNG or RANDOM123 you must have ./configure PETSc
2440:    with the option --download-sprng or --download-random123

2442: .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2443: J*/
2444: typedef const char* PetscRandomType;
2445: #define PETSCRAND       "rand"
2446: #define PETSCRAND48     "rand48"
2447: #define PETSCSPRNG      "sprng"
2448: #define PETSCRANDER48   "rander48"
2449: #define PETSCRANDOM123  "random123"

2451: /* Logging support */
2452: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;

2454: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);

2456: /* Dynamic creation and loading functions */
2457: PETSC_EXTERN PetscFunctionList PetscRandomList;

2459: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2460: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2461: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2462: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType*);
2463: PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom,PetscObject,const char[]);
2464: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);

2466: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2467: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2468: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2469: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2470: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2471: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2472: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2473: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2474: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);

2476: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2477: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2478: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2479: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2480: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2481: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2482: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);
2483: PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2484: PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2485: PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);

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

2489: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscInt*,PetscDataType);
2490: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscInt*,PetscDataType);
2491: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,const void*,PetscInt,PetscDataType);
2492: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,const void*,PetscInt,PetscDataType);
2493: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2494: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2495: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2496: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2497: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2498: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2499: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2500: #if defined(PETSC_USE_SOCKET_VIEWER)
2501: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);
2502: #endif

2504: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2505: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2506: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);

2508: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2509: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool );
2510: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2511: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2512: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2513: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);

2515: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2516: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2517: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2518: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2519: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2520: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt*,const void*,PetscMPIInt*,PetscMPIInt**,void*)
2521:   PetscAttrMPIPointerWithType(6,3);
2522: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2523:                                                     PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2524:                                                     PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2525:   PetscAttrMPIPointerWithType(6,3);
2526: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2527:                                                        MPI_Request**,MPI_Request**,
2528:                                                        PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2529:                                                        PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2530:   PetscAttrMPIPointerWithType(6,3);

2532: PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2533: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2534: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);

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

2538: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);

2540: PETSC_EXTERN const char *const PetscSubcommTypes[];

2542: struct _n_PetscSubcomm {
2543:   MPI_Comm         parent;           /* parent communicator */
2544:   MPI_Comm         dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2545:   MPI_Comm         child;            /* the sub-communicator */
2546:   PetscMPIInt      n;                /* num of subcommunicators under the parent communicator */
2547:   PetscMPIInt      color;            /* color of processors belong to this communicator */
2548:   PetscMPIInt      *subsize;         /* size of subcommunicator[color] */
2549:   PetscSubcommType type;
2550:   char             *subcommprefix;
2551: };

2553: PETSC_STATIC_INLINE MPI_Comm PetscSubcommParent(PetscSubcomm scomm) {return scomm->parent;}
2554: PETSC_STATIC_INLINE MPI_Comm PetscSubcommChild(PetscSubcomm scomm) {return scomm->child;}
2555: PETSC_STATIC_INLINE MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm) {return scomm->dupparent;}
2556: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2557: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2558: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2559: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2560: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2561: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2562: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2563: PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm,const char[]);
2564: PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm,MPI_Comm*);
2565: PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm,MPI_Comm*);
2566: PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm,MPI_Comm*);

2568: PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt,PetscHeap*);
2569: PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap,PetscInt,PetscInt);
2570: PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap,PetscInt*,PetscInt*);
2571: PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap,PetscInt*,PetscInt*);
2572: PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap,PetscInt,PetscInt);
2573: PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2574: PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap*);
2575: PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap,PetscViewer);

2577: PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2578: PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm,PetscShmComm*);
2579: PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2580: PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2581: PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm,MPI_Comm*);

2583: /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2584: PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm,PetscInt,PetscOmpCtrl*);
2585: PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl,MPI_Comm*,MPI_Comm*,PetscBool*);
2586: PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl*);
2587: PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2588: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2589: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);

2591: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2592: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2593: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2594: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2595: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2596: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2597: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2598: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);


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

2606: extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2607: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted*);
2608: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted*);
2609: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted,const char*,const char*,PetscBool*);

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

2615: PETSC_EXTERN PetscSegBuffer PetscCitationsList;

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

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

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

2626:    Level: intermediate

2628:    Not available from Fortran

2630:      Options Database:
2631: .     -citations [filename]   - print out the bibtex entries for the given computation
2632: @*/
2633: PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2634: {
2635:   size_t         len;
2636:   char           *vstring;

2640:   if (set && *set) return(0);
2641:   PetscStrlen(cit,&len);
2642:   PetscSegBufferGet(PetscCitationsList,len,&vstring);
2643:   PetscArraycpy(vstring,cit,len);
2644:   if (set) *set = PETSC_TRUE;
2645:   return(0);
2646: }

2648: PETSC_EXTERN PetscErrorCode PetscURLShorten(const char[],char[],size_t);
2649: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm,char[],char[],size_t);
2650: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm,const char[],char[],size_t);
2651: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm,const char[],const char []);

2653: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm,char[],char[],size_t);
2654: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm,const char[],char[],char[],size_t);

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

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

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


2665: #if defined(PETSC_USE_DEBUG)
2666: /*
2667:    Verify that all processes in the communicator have called this from the same line of code
2668:  */
2669: PETSC_EXTERN PetscErrorCode PetscAllreduceBarrierCheck(MPI_Comm,PetscMPIInt,int,const char*,const char *);

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

2676:    Collective

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

2682:    Input Parameters:
2683: +  indata - pointer to the input data to be reduced
2684: .  count - the number of MPI data items in indata and outdata
2685: .  datatype - the MPI datatype, for example MPI_INT
2686: .  op - the MPI operation, for example MPI_SUM
2687: -  comm - the MPI communicator on which the operation occurs

2689:    Output Parameter:
2690: .  outdata - the reduced values

2692:    Notes:
2693:    In optimized mode this directly calls MPI_Allreduce()

2695:    Level: developer

2697: .seealso: MPI_Allreduce()
2698: M*/
2699: #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))
2700: #else
2701: #define MPIU_Allreduce(a,b,c,d,e,fcomm) MPI_Allreduce(a,b,c,d,e,fcomm)
2702: #endif

2704: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
2705: PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint,PetscMPIInt,MPI_Info,MPI_Comm,void*,MPI_Win*);
2706: PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win,PetscMPIInt,MPI_Aint*,PetscMPIInt*,void*);
2707: #endif

2709: /*
2710:     Returned from PETSc functions that are called from MPI, such as related to attributes
2711: */
2712: PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CLASS;
2713: PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CODE;

2715: #endif