Actual source code: petscsys.h

  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 that
 11:    PETSc's makefiles add to the compiler rules.
 12:    For --prefix installs the directory ${PETSC_ARCH} does not exist and petscconf.h is in the same
 13:    directory as the other PETSc include files.
 14: */
 15: #include <petscconf.h>
 16: #include <petscconf_poison.h>
 17: #include <petscfix.h>

 19: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP) || defined(PETSC_HAVE_KOKKOS)
 20:    #define PETSC_HAVE_DEVICE
 21: #endif

 23: #if defined(PETSC_DESIRE_FEATURE_TEST_MACROS)
 24: /*
 25:    Feature test macros must be included before headers defined by IEEE Std 1003.1-2001
 26:    We only turn these in PETSc source files that require them by setting PETSC_DESIRE_FEATURE_TEST_MACROS
 27: */
 28: #  if defined(PETSC__POSIX_C_SOURCE_200112L) && !defined(_POSIX_C_SOURCE)
 29: #    define _POSIX_C_SOURCE 200112L
 30: #  endif
 31: #  if defined(PETSC__BSD_SOURCE) && !defined(_BSD_SOURCE)
 32: #    define _BSD_SOURCE
 33: #  endif
 34: #  if defined(PETSC__DEFAULT_SOURCE) && !defined(_DEFAULT_SOURCE)
 35: #    define _DEFAULT_SOURCE
 36: #  endif
 37: #  if defined(PETSC__GNU_SOURCE) && !defined(_GNU_SOURCE)
 38: #    define _GNU_SOURCE
 39: #  endif
 40: #endif

 42: #include <petscsystypes.h>

 44: /* ========================================================================== */
 45: /*
 46:    This facilitates using the C version of PETSc from C++ and the C++ version from C.
 47: */
 48: #if defined(__cplusplus)
 49: #  define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_CXX
 50: #else
 51: #  define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_C
 52: #endif

 54: /* ========================================================================== */
 55: /*
 56:    Since PETSc manages its own extern "C" handling users should never include PETSc include
 57:    files within extern "C". This will generate a compiler error if a user does put the include
 58:    file within an extern "C".
 59: */
 60: #if defined(__cplusplus)
 61: void assert_never_put_petsc_headers_inside_an_extern_c(int); void assert_never_put_petsc_headers_inside_an_extern_c(double);
 62: #endif

 64: #if defined(__cplusplus)
 65: #  define PETSC_RESTRICT PETSC_CXX_RESTRICT
 66: #else
 67: #  define PETSC_RESTRICT PETSC_C_RESTRICT
 68: #endif

 70: #if defined(__cplusplus)
 71: #  define PETSC_INLINE PETSC_CXX_INLINE
 72: #else
 73: #  define PETSC_INLINE PETSC_C_INLINE
 74: #endif

 76: #define PETSC_STATIC_INLINE static PETSC_INLINE

 78: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES) /* For Win32 shared libraries */
 79: #  define  __declspec(dllexport)
 80: #  define PETSC_DLLIMPORT __declspec(dllimport)
 81: #  define PETSC_VISIBILITY_INTERNAL
 82: #elif defined(PETSC_USE_VISIBILITY_CXX) && defined(__cplusplus)
 83: #  define  __attribute__((visibility ("default")))
 84: #  define PETSC_DLLIMPORT __attribute__((visibility ("default")))
 85: #  define PETSC_VISIBILITY_INTERNAL __attribute__((visibility ("hidden")))
 86: #elif defined(PETSC_USE_VISIBILITY_C) && !defined(__cplusplus)
 87: #  define  __attribute__((visibility ("default")))
 88: #  define PETSC_DLLIMPORT __attribute__((visibility ("default")))
 89: #  define PETSC_VISIBILITY_INTERNAL __attribute__((visibility ("hidden")))
 90: #else
 91: #  define 
 92: #  define PETSC_DLLIMPORT
 93: #  define PETSC_VISIBILITY_INTERNAL
 94: #endif

 96: #if defined(petsc_EXPORTS)      /* CMake defines this when building the shared library */
 97: #  define PETSC_VISIBILITY_PUBLIC 
 98: #else  /* Win32 users need this to import symbols from petsc.dll */
 99: #  define PETSC_VISIBILITY_PUBLIC PETSC_DLLIMPORT
100: #endif

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

117: #if defined(__cplusplus) && defined(PETSC_HAVE_CXX_DIALECT_CXX11)
118: #  define PETSC_NULLPTR             nullptr
119: #  define PETSC_CONSTEXPR           constexpr
120: #  define PETSC_NOEXCEPT            noexcept
121: #  define PETSC_NOEXCEPT_ARG(cond_) noexcept(cond_)
122: #  define PETSC_CXX_DEFAULT(func_)  func_ = default
123: #else
124: #  define PETSC_NULLPTR             NULL
125: #  define PETSC_CONSTEXPR
126: #  define PETSC_NOEXCEPT
127: #  define PETSC_NOEXCEPT_ARG(cond_)
128: #  define PETSC_CXX_DEFAULT(func_)
129: #endif

131: /* C++17 features */
132: /* We met cases that the host CXX compiler (say mpicxx) supports C++17, but nvcc does not agree, even with -ccbin mpicxx! */
133: #if defined(__cplusplus) && defined(PETSC_HAVE_CXX_DIALECT_CXX17) && (!defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_CUDA_DIALECT_CXX17))
134: #  define PETSC_NODISCARD    [[nodiscard]]
135: #else
136: #  define PETSC_NODISCARD
137: #endif

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

142: /* ========================================================================== */

144: /*
145:     Defines the interface to MPI allowing the use of all MPI functions.

147:     PETSc does not use the C++ binding of MPI at ALL. The following flag
148:     makes sure the C++ bindings are not included. The C++ bindings REQUIRE
149:     putting mpi.h before ANY C++ include files, we cannot control this
150:     with all PETSc users. Users who want to use the MPI C++ bindings can include
151:     mpicxx.h directly in their code
152: */
153: #if !defined(MPICH_SKIP_MPICXX)
154: #  define MPICH_SKIP_MPICXX 1
155: #endif
156: #if !defined(OMPI_SKIP_MPICXX)
157: #  define OMPI_SKIP_MPICXX 1
158: #endif
159: #if defined(PETSC_HAVE_MPIUNI)
160: #  include <petsc/mpiuni/mpi.h>
161: #else
162: #  include <mpi.h>
163: #endif

165: /*MC
166:   PetscDefined - determine whether a boolean macro is defined

168:   Notes:
169:   The prefix "PETSC_" is added to the argument.

171:   Typical usage is within normal code,

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

175:   but can also be used in the preprocessor,

177: $   #if PetscDefined(USE_DEBUG)
178: $     ...
179: $   #else

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

184:   To avoid prepending "PETSC_", say to add custom checks in user code, one can use e.g.

186: $  #define FooDefined(d) PetscDefined_(FOO_ ## d)

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

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

197:   Level: developer
198: M*/
199: #if !defined(PETSC_SKIP_VARIADIC_MACROS)
200: #  define PetscDefined_arg_1    shift,
201: #  define PetscDefined_arg_     shift,
202: #  define PetscDefined__take_second_expanded(ignored, val, ...) val
203: #  define PetscDefined__take_second_expand(args) PetscDefined__take_second_expanded args
204: #  define PetscDefined__take_second(...) PetscDefined__take_second_expand((__VA_ARGS__))
205: #  define PetscDefined___(arg1_or_junk) PetscDefined__take_second(arg1_or_junk 1, 0, at_)
206: #  define PetscDefined__(value) PetscDefined___(PetscDefined_arg_ ## value)
207: #  define PetscDefined_(d)      PetscDefined__(d)
208: #  define PetscDefined(d)       PetscDefined_(PETSC_ ## d)
209: #endif

211: /*
212:    Perform various sanity checks that the correct mpi.h is being included at compile time.
213:    This usually happens because
214:       * either an unexpected mpi.h is in the default compiler path (i.e. in /usr/include) or
215:       * an extra include path -I/something (which contains the unexpected mpi.h) is being passed to the compiler
216: */
217: #if defined(PETSC_HAVE_MPIUNI)
218: #  if !defined(MPIUNI_H)
219: #    error "PETSc was configured with --with-mpi=0 but now appears to be compiling using a different mpi.h"
220: #  endif
221: #elif defined(PETSC_HAVE_I_MPI_NUMVERSION)
222: #  if !defined(I_MPI_NUMVERSION)
223: #    error "PETSc was configured with I_MPI but now appears to be compiling using a non-I_MPI mpi.h"
224: #  elif I_MPI_NUMVERSION != PETSC_HAVE_I_MPI_NUMVERSION
225: #    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"
226: #  endif
227: #elif defined(PETSC_HAVE_MVAPICH2_NUMVERSION)
228: #  if !defined(MVAPICH2_NUMVERSION)
229: #    error "PETSc was configured with MVAPICH2 but now appears to be compiling using a non-MVAPICH2 mpi.h"
230: #  elif MVAPICH2_NUMVERSION != PETSC_HAVE_MVAPICH2_NUMVERSION
231: #    error "PETSc was configured with one MVAPICH2 mpi.h version but now appears to be compiling using a different MVAPICH2 mpi.h version"
232: #  endif
233: #elif defined(PETSC_HAVE_MPICH_NUMVERSION)
234: #  if !defined(MPICH_NUMVERSION) || defined(MVAPICH2_NUMVERSION) || defined(I_MPI_NUMVERSION)
235: #    error "PETSc was configured with MPICH but now appears to be compiling using a non-MPICH mpi.h"
236: #  elif (MPICH_NUMVERSION/100000 != PETSC_HAVE_MPICH_NUMVERSION/100000) || (MPICH_NUMVERSION%100000/1000 < PETSC_HAVE_MPICH_NUMVERSION%100000/1000)
237: #    error "PETSc was configured with one MPICH mpi.h version but now appears to be compiling using a different MPICH mpi.h version"
238: #  endif
239: #elif defined(PETSC_HAVE_OMPI_MAJOR_VERSION)
240: #  if !defined(OMPI_MAJOR_VERSION)
241: #    error "PETSc was configured with OpenMPI but now appears to be compiling using a non-OpenMPI mpi.h"
242: #  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)
243: #    error "PETSc was configured with one OpenMPI mpi.h version but now appears to be compiling using a different OpenMPI mpi.h version"
244: #  endif
245: #elif defined(PETSC_HAVE_MSMPI_VERSION)
246: #  if !defined(MSMPI_VER)
247: #    error "PETSc was configured with MSMPI but now appears to be compiling using a non-MSMPI mpi.h"
248: #  elif (MSMPI_VER != PETSC_HAVE_MSMPI_VERSION)
249: #    error "PETSc was configured with one MSMPI mpi.h version but now appears to be compiling using a different MSMPI mpi.h version"
250: #  endif
251: #elif defined(OMPI_MAJOR_VERSION) || defined(MPICH_NUMVERSION) || defined(MSMPI_VER)
252: #  error "PETSc was configured with undetermined MPI - but now appears to be compiling using any of OpenMPI, MS-MPI or a MPICH variant"
253: #endif

255: /*
256:     Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler
257:     see the top of mpicxx.h in the MPICH2 distribution.
258: */
259: #include <stdio.h>

261: /* MSMPI on 32bit windows requires this yukky hack - that breaks MPI standard compliance */
262: #if !defined(MPIAPI)
263: #define MPIAPI
264: #endif

266: /*
267:    Support for Clang (>=3.2) matching type tag arguments with void* buffer types.
268:    This allows the compiler to detect cases where the MPI datatype argument passed to a MPI routine
269:    does not match the actual type of the argument being passed in
270: */
271: #if defined(__has_attribute) && defined(works_with_const_which_is_not_true)
272: #  if __has_attribute(argument_with_type_tag) && __has_attribute(pointer_with_type_tag) && __has_attribute(type_tag_for_datatype)
273: #    define PetscAttrMPIPointerWithType(bufno,typeno) __attribute__((pointer_with_type_tag(MPI,bufno,typeno)))
274: #    define PetscAttrMPITypeTag(type)                 __attribute__((type_tag_for_datatype(MPI,type)))
275: #    define PetscAttrMPITypeTagLayoutCompatible(type) __attribute__((type_tag_for_datatype(MPI,type,layout_compatible)))
276: #  endif
277: #endif
278: #if !defined(PetscAttrMPIPointerWithType)
279: #  define PetscAttrMPIPointerWithType(bufno,typeno)
280: #  define PetscAttrMPITypeTag(type)
281: #  define PetscAttrMPITypeTagLayoutCompatible(type)
282: #endif

284: PETSC_EXTERN MPI_Datatype MPIU_ENUM PetscAttrMPITypeTag(PetscEnum);
285: PETSC_EXTERN MPI_Datatype MPIU_BOOL PetscAttrMPITypeTag(PetscBool);

287: /*MC
288:    MPIU_INT - MPI datatype corresponding to PetscInt

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

293:    Level: beginner

295: .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_SCALAR, MPIU_COMPLEX
296: M*/

298: #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 */
299: #  define MPIU_INT64 MPI_INT64_T
300: #  define PetscInt64_FMT PRId64
301: #elif (PETSC_SIZEOF_LONG_LONG == 8)
302: #  define MPIU_INT64 MPI_LONG_LONG_INT
303: #  define PetscInt64_FMT "lld"
304: #elif defined(PETSC_HAVE___INT64)
305: #  define MPIU_INT64 MPI_INT64_T
306: #  define PetscInt64_FMT "ld"
307: #else
308: #  error "cannot determine PetscInt64 type"
309: #endif

311: PETSC_EXTERN MPI_Datatype MPIU_FORTRANADDR;

313: #if defined(PETSC_USE_64BIT_INDICES)
314: #  define MPIU_INT MPIU_INT64
315: #  define PetscInt_FMT PetscInt64_FMT
316: #else
317: #  define MPIU_INT MPI_INT
318: #  define PetscInt_FMT "d"
319: #endif

321: /*
322:     For the rare cases when one needs to send a size_t object with MPI
323: */
324: PETSC_EXTERN MPI_Datatype MPIU_SIZE_T;

326: /*
327:       You can use PETSC_STDOUT as a replacement of stdout. You can also change
328:     the value of PETSC_STDOUT to redirect all standard output elsewhere
329: */
330: PETSC_EXTERN FILE* PETSC_STDOUT;

332: /*
333:       You can use PETSC_STDERR as a replacement of stderr. You can also change
334:     the value of PETSC_STDERR to redirect all standard error elsewhere
335: */
336: PETSC_EXTERN FILE* PETSC_STDERR;

338: /*MC
339:     PetscUnlikely - hints the compiler that the given condition is usually FALSE

341:     Synopsis:
342: #include <petscsys.h>
343:     PetscBool  PetscUnlikely(PetscBool  cond)

345:     Not Collective

347:     Input Parameters:
348: .   cond - condition or expression

350:     Notes:
351:     This returns the same truth value, it is only a hint to compilers that the resulting
352:     branch is unlikely.

354:     Level: advanced

356: .seealso: PetscUnlikelyDebug(), PetscLikely(), CHKERRQ
357: M*/

359: /*MC
360:     PetscLikely - hints the compiler that the given condition is usually TRUE

362:     Synopsis:
363: #include <petscsys.h>
364:     PetscBool  PetscLikely(PetscBool  cond)

366:     Not Collective

368:     Input Parameters:
369: .   cond - condition or expression

371:     Notes:
372:     This returns the same truth value, it is only a hint to compilers that the resulting
373:     branch is likely.

375:     Level: advanced

377: .seealso: PetscUnlikely()
378: M*/
379: #if defined(PETSC_HAVE_BUILTIN_EXPECT)
380: #  define PetscUnlikely(cond)   __builtin_expect(!!(cond),0)
381: #  define PetscLikely(cond)     __builtin_expect(!!(cond),1)
382: #else
383: #  define PetscUnlikely(cond)   (cond)
384: #  define PetscLikely(cond)     (cond)
385: #endif

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

390:     Synopsis:
391: #include <petscsys.h>
392:     PetscBool  PetscUnlikelyDebug(PetscBool  cond)

394:     Not Collective

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

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

403:     Level: advanced

405: .seealso: PetscUnlikely(), CHKERRQ, SETERRQ
406: M*/
407: #define PetscUnlikelyDebug(cond) (PetscDefined(USE_DEBUG) && PetscUnlikely(cond))

409: /* PetscPragmaSIMD - from CeedPragmaSIMD */

411: #if defined(__NEC__)
412: #  define PetscPragmaSIMD _Pragma("_NEC ivdep")
413: #elif defined(__INTEL_COMPILER) && !defined(_WIN32)
414: #  define PetscPragmaSIMD _Pragma("vector")
415: #elif defined(__GNUC__) && __GNUC__ >= 5 && !defined(__PGI)
416: #  define PetscPragmaSIMD _Pragma("GCC ivdep")
417: #elif defined(_OPENMP) && _OPENMP >= 201307
418: #  if defined(_MSC_VER)
419: #    define PetscPragmaSIMD __pragma(omp simd)
420: #  else
421: #    define PetscPragmaSIMD _Pragma("omp simd")
422: #  endif
423: #elif defined(PETSC_HAVE_CRAY_VECTOR)
424: #  define PetscPragmaSIMD _Pragma("_CRI ivdep")
425: #else
426: #  define PetscPragmaSIMD
427: #endif

429: /*
430:     Declare extern C stuff after including external header files
431: */

433: PETSC_EXTERN const char *const PetscBools[];

435: PETSC_EXTERN PetscBool PETSC_RUNNING_ON_VALGRIND;
436: /*
437:     Defines elementary mathematics functions and constants.
438: */
439: #include <petscmath.h>

441: PETSC_EXTERN const char *const PetscCopyModes[];

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

446:    Level: beginner

448:    Note:
449:    Accepted by many PETSc functions to not set a parameter and instead use
450:           some default

452:    Fortran Notes:
453:    This macro does not exist in Fortran; you must use PETSC_NULL_INTEGER,
454:           PETSC_NULL_DOUBLE_PRECISION etc

456: .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_DETERMINE

458: M*/
459: #define PETSC_IGNORE NULL

461: /* This is deprecated */
462: #define PETSC_NULL NULL

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

468:    Level: beginner

470: .seealso: PETSC_DEFAULT, PETSC_IGNORE, PETSC_DETERMINE

472: M*/
473: #define PETSC_DECIDE -1

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

479:    Level: beginner

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

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

487: M*/
488: #define PETSC_DETERMINE PETSC_DECIDE

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

494:    Level: beginner

496:    Fortran Notes:
497:    You need to use PETSC_DEFAULT_INTEGER or PETSC_DEFAULT_REAL.

499: .seealso: PETSC_DECIDE, PETSC_IGNORE, PETSC_DETERMINE

501: M*/
502: #define PETSC_DEFAULT -2

504: /*MC
505:     PETSC_COMM_WORLD - the equivalent of the MPI_COMM_WORLD communicator which represents
506:            all the processes that PETSc knows about.

508:    Level: beginner

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

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

519: .seealso: PETSC_COMM_SELF

521: M*/
522: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;

524: /*MC
525:     PETSC_COMM_SELF - This is always MPI_COMM_SELF

527:    Level: beginner

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

532: .seealso: PETSC_COMM_WORLD

534: M*/
535: #define PETSC_COMM_SELF MPI_COMM_SELF

537: /*MC
538:     PETSC_MPI_THREAD_REQUIRED - the required threading support used if PETSc initializes
539:            MPI with MPI_Init_thread.

541:    Level: beginner

543:    Notes:
544:    By default PETSC_MPI_THREAD_REQUIRED equals MPI_THREAD_FUNNELED.

546: .seealso: PetscInitialize()

548: M*/
549: PETSC_EXTERN PetscMPIInt PETSC_MPI_THREAD_REQUIRED;

551: PETSC_EXTERN PetscBool PetscBeganMPI;
552: PETSC_EXTERN PetscBool PetscErrorHandlingInitialized;
553: PETSC_EXTERN PetscBool PetscInitializeCalled;
554: PETSC_EXTERN PetscBool PetscFinalizeCalled;
555: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;

557: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm));
558: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*);
559: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm*);

561: #if defined(PETSC_HAVE_CUDA)
562: PETSC_EXTERN PetscBool      PetscCUDASynchronize;
563: PETSC_EXTERN PetscErrorCode PetscCUDAInitialize(MPI_Comm,PetscInt);
564: PETSC_EXTERN PetscErrorCode PetscCUDAInitializeCheck(void);
565: #endif

567: #if defined(PETSC_HAVE_HIP)
568: PETSC_EXTERN PetscBool      PetscHIPSynchronize;
569: PETSC_EXTERN PetscErrorCode PetscHIPInitialize(MPI_Comm,PetscInt);
570: PETSC_EXTERN PetscErrorCode PetscHIPInitializeCheck(void);
571: #endif

573: #if defined(PETSC_HAVE_KOKKOS)
574: PETSC_EXTERN PetscErrorCode PetscKokkosInitializeCheck(void);  /* Initialize Kokkos if not yet. */
575: #endif

577: #if defined(PETSC_HAVE_NVSHMEM)
578: PETSC_EXTERN PetscBool      PetscBeganNvshmem;
579: PETSC_EXTERN PetscBool      PetscNvshmemInitialized;
580: PETSC_EXTERN PetscErrorCode PetscNvshmemFinalize(void);
581: #endif

583: #if defined(PETSC_HAVE_ELEMENTAL)
584: PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
585: PETSC_EXTERN PetscErrorCode PetscElementalInitialized(PetscBool*);
586: PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
587: #endif

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

592:    Synopsis:
593: #include <petscsys.h>
594:    PetscErrorCode PetscMalloc(size_t m,void **result)

596:    Not Collective

598:    Input Parameter:
599: .  m - number of bytes to allocate

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

604:    Level: beginner

606:    Notes:
607:    Memory is always allocated at least double aligned

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

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

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

616: /*MC
617:    PetscRealloc - Rellocates memory

619:    Synopsis:
620: #include <petscsys.h>
621:    PetscErrorCode PetscRealloc(size_t m,void **result)

623:    Not Collective

625:    Input Parameters:
626: +  m - number of bytes to allocate
627: -  result - previous memory

629:    Output Parameter:
630: .  result - new memory allocated

632:    Level: developer

634:    Notes:
635:    Memory is always allocated at least double aligned

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

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

642: /*MC
643:    PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment

645:    Synopsis:
646: #include <petscsys.h>
647:    void *PetscAddrAlign(void *addr)

649:    Not Collective

651:    Input Parameters:
652: .  addr - address to align (any pointer type)

654:    Level: developer

656: .seealso: PetscMallocAlign()

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

661: /*MC
662:    PetscCalloc - Allocates a cleared (zeroed) memory region aligned to PETSC_MEMALIGN

664:    Synopsis:
665: #include <petscsys.h>
666:    PetscErrorCode PetscCalloc(size_t m,void **result)

668:    Not Collective

670:    Input Parameter:
671: .  m - number of bytes to allocate

673:    Output Parameter:
674: .  result - memory allocated

676:    Level: beginner

678:    Notes:
679:    Memory is always allocated at least double aligned. This macro is useful in allocating memory pointed by void pointers

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

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

685: M*/
686: #define PetscCalloc(m,result)  PetscMallocA(1,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m),(result))

688: /*MC
689:    PetscMalloc1 - Allocates an array of memory aligned to PETSC_MEMALIGN

691:    Synopsis:
692: #include <petscsys.h>
693:    PetscErrorCode PetscMalloc1(size_t m1,type **r1)

695:    Not Collective

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

700:    Output Parameter:
701: .  r1 - memory allocated

703:    Note:
704:    This uses the sizeof() of the memory type requested to determine the total memory to be allocated, therefore you should not
705:          multiply the number of elements requested by the sizeof() the type. For example use
706: $  PetscInt *id;
707: $  PetscMalloc1(10,&id);
708:         not
709: $  PetscInt *id;
710: $  PetscMalloc1(10*sizeof(PetscInt),&id);

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

714:    Level: beginner

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

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

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

724:    Synopsis:
725: #include <petscsys.h>
726:    PetscErrorCode PetscCalloc1(size_t m1,type **r1)

728:    Not Collective

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

733:    Output Parameter:
734: .  r1 - memory allocated

736:    Notes:
737:    See PetsMalloc1() for more details on usage.

739:    Level: beginner

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

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

746: /*MC
747:    PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN

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

753:    Not Collective

755:    Input Parameters:
756: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
757: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

759:    Output Parameters:
760: +  r1 - memory allocated in first chunk
761: -  r2 - memory allocated in second chunk

763:    Level: developer

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

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

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

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

777:    Not Collective

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

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

787:    Level: developer

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

791: M*/
792: #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))

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

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

801:    Not Collective

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

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

813:    Level: developer

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

817: M*/
818: #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))

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

823:    Synopsis:
824: #include <petscsys.h>
825:    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

827:    Not Collective

829:    Input Parameters:
830: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
831: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
832: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

834:    Output Parameters:
835: +  r1 - memory allocated in first chunk
836: .  r2 - memory allocated in second chunk
837: -  r3 - memory allocated in third chunk

839:    Level: developer

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

843: M*/
844: #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))

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

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

853:    Not Collective

855:    Input Parameters:
856: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
857: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
858: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
859: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

861:    Output Parameters:
862: +  r1 - memory allocated in first chunk
863: .  r2 - memory allocated in second chunk
864: .  r3 - memory allocated in third chunk
865: -  r4 - memory allocated in fourth chunk

867:    Level: developer

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

871: M*/
872: #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))

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

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

881:    Not Collective

883:    Input Parameters:
884: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
885: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
886: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
887: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

889:    Output Parameters:
890: +  r1 - memory allocated in first chunk
891: .  r2 - memory allocated in second chunk
892: .  r3 - memory allocated in third chunk
893: -  r4 - memory allocated in fourth chunk

895:    Level: developer

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

899: M*/
900: #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))

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

905:    Synopsis:
906: #include <petscsys.h>
907:    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)

909:    Not Collective

911:    Input Parameters:
912: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
913: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
914: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
915: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
916: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

918:    Output Parameters:
919: +  r1 - memory allocated in first chunk
920: .  r2 - memory allocated in second chunk
921: .  r3 - memory allocated in third chunk
922: .  r4 - memory allocated in fourth chunk
923: -  r5 - memory allocated in fifth chunk

925:    Level: developer

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

929: M*/
930: #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))

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

935:    Synopsis:
936: #include <petscsys.h>
937:    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)

939:    Not Collective

941:    Input Parameters:
942: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
943: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
944: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
945: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
946: -  m5 - number of elements to allocate in 5th chunk  (may be zero)

948:    Output Parameters:
949: +  r1 - memory allocated in first chunk
950: .  r2 - memory allocated in second chunk
951: .  r3 - memory allocated in third chunk
952: .  r4 - memory allocated in fourth chunk
953: -  r5 - memory allocated in fifth chunk

955:    Level: developer

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

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

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

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

969:    Not Collective

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

979:    Output Parameteasr:
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

987:    Level: developer

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

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

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

997:    Synopsis:
998: #include <petscsys.h>
999:    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)

1001:    Not Collective

1003:    Input Parameters:
1004: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
1005: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
1006: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
1007: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
1008: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
1009: -  m6 - number of elements to allocate in 6th chunk  (may be zero)

1011:    Output Parameters:
1012: +  r1 - memory allocated in first chunk
1013: .  r2 - memory allocated in second chunk
1014: .  r3 - memory allocated in third chunk
1015: .  r4 - memory allocated in fourth chunk
1016: .  r5 - memory allocated in fifth chunk
1017: -  r6 - memory allocated in sixth chunk

1019:    Level: developer

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

1023: M*/
1024: #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))

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

1029:    Synopsis:
1030: #include <petscsys.h>
1031:    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)

1033:    Not Collective

1035:    Input Parameters:
1036: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
1037: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
1038: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
1039: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
1040: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
1041: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
1042: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

1044:    Output Parameters:
1045: +  r1 - memory allocated in first chunk
1046: .  r2 - memory allocated in second chunk
1047: .  r3 - memory allocated in third chunk
1048: .  r4 - memory allocated in fourth chunk
1049: .  r5 - memory allocated in fifth chunk
1050: .  r6 - memory allocated in sixth chunk
1051: -  r7 - memory allocated in seventh chunk

1053:    Level: developer

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

1057: M*/
1058: #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))

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

1063:    Synopsis:
1064: #include <petscsys.h>
1065:    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)

1067:    Not Collective

1069:    Input Parameters:
1070: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
1071: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
1072: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
1073: .  m4 - number of elements to allocate in 4th chunk  (may be zero)
1074: .  m5 - number of elements to allocate in 5th chunk  (may be zero)
1075: .  m6 - number of elements to allocate in 6th chunk  (may be zero)
1076: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

1078:    Output Parameters:
1079: +  r1 - memory allocated in first chunk
1080: .  r2 - memory allocated in second chunk
1081: .  r3 - memory allocated in third chunk
1082: .  r4 - memory allocated in fourth chunk
1083: .  r5 - memory allocated in fifth chunk
1084: .  r6 - memory allocated in sixth chunk
1085: -  r7 - memory allocated in seventh chunk

1087:    Level: developer

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

1091: M*/
1092: #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))

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

1097:    Synopsis:
1098: #include <petscsys.h>
1099:    PetscErrorCode PetscNew(type **result)

1101:    Not Collective

1103:    Output Parameter:
1104: .  result - memory allocated, sized to match pointer type

1106:    Level: beginner

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

1110: M*/
1111: #define PetscNew(b)      PetscCalloc1(1,(b))

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

1117:    Synopsis:
1118: #include <petscsys.h>
1119:    PetscErrorCode PetscNewLog(PetscObject obj,type **result)

1121:    Not Collective

1123:    Input Parameter:
1124: .  obj - object memory is logged to

1126:    Output Parameter:
1127: .  result - memory allocated, sized to match pointer type

1129:    Level: developer

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

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

1136: /*MC
1137:    PetscFree - Frees memory

1139:    Synopsis:
1140: #include <petscsys.h>
1141:    PetscErrorCode PetscFree(void *memory)

1143:    Not Collective

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

1148:    Level: beginner

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

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

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

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

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

1163:    Synopsis:
1164: #include <petscsys.h>
1165:    PetscErrorCode PetscFree2(void *memory1,void *memory2)

1167:    Not Collective

1169:    Input Parameters:
1170: +   memory1 - memory to free
1171: -   memory2 - 2nd memory to free

1173:    Level: developer

1175:    Notes:
1176:     Memory must have been obtained with PetscMalloc2()

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

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

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

1186:    Synopsis:
1187: #include <petscsys.h>
1188:    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)

1190:    Not Collective

1192:    Input Parameters:
1193: +   memory1 - memory to free
1194: .   memory2 - 2nd memory to free
1195: -   memory3 - 3rd memory to free

1197:    Level: developer

1199:    Notes:
1200:     Memory must have been obtained with PetscMalloc3()

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

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

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

1210:    Synopsis:
1211: #include <petscsys.h>
1212:    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)

1214:    Not Collective

1216:    Input Parameters:
1217: +   m1 - memory to free
1218: .   m2 - 2nd memory to free
1219: .   m3 - 3rd memory to free
1220: -   m4 - 4th memory to free

1222:    Level: developer

1224:    Notes:
1225:     Memory must have been obtained with PetscMalloc4()

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

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

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

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

1239:    Not Collective

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

1248:    Level: developer

1250:    Notes:
1251:     Memory must have been obtained with PetscMalloc5()

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

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

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

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

1265:    Not Collective

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

1275:    Level: developer

1277:    Notes:
1278:     Memory must have been obtained with PetscMalloc6()

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

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

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

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

1292:    Not Collective

1294:    Input Parameters:
1295: +   m1 - memory to free
1296: .   m2 - 2nd memory to free
1297: .   m3 - 3rd memory to free
1298: .   m4 - 4th memory to free
1299: .   m5 - 5th memory to free
1300: .   m6 - 6th memory to free
1301: -   m7 - 7th memory to free

1303:    Level: developer

1305:    Notes:
1306:     Memory must have been obtained with PetscMalloc7()

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

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

1314: PETSC_EXTERN PetscErrorCode PetscMallocA(int,PetscBool,int,const char *,const char *,size_t,void *,...);
1315: PETSC_EXTERN PetscErrorCode PetscFreeA(int,int,const char *,const char *,void *,...);
1316: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,PetscBool,int,const char[],const char[],void**);
1317: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1318: PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t,int,const char[],const char[],void**);
1319: PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1320: 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 **));
1321: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);

1323: /*
1324:   Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1325: */
1326: PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1327: PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);
1328: #if defined(PETSC_HAVE_CUDA)
1329: PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void);
1330: PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void);
1331: #endif
1332: #if defined(PETSC_HAVE_HIP)
1333: PETSC_EXTERN PetscErrorCode PetscMallocSetHIPHost(void);
1334: PETSC_EXTERN PetscErrorCode PetscMallocResetHIPHost(void);
1335: #endif

1337: #define MPIU_PETSCLOGDOUBLE  MPI_DOUBLE
1338: #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION

1340: /*
1341:    Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1342: */
1343: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1344: PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1345: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1346: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1347: PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1348: PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int,PetscLogDouble*);
1349: PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool,PetscBool);
1350: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*,PetscBool*,PetscBool*);
1351: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1352: PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1353: PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool*);
1354: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool);
1355: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool*);

1357: PETSC_EXTERN const char *const PetscDataTypes[];
1358: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1359: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1360: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1361: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);

1363: /*
1364:     Basic memory and string operations. These are usually simple wrappers
1365:    around the basic Unix system calls, but a few of them have additional
1366:    functionality and/or error checking.
1367: */
1368: PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1369: PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1370: PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***);
1371: PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1372: PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool  *);
1373: PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool  *);
1374: PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1375: PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1376: PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1377: PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1378: PETSC_EXTERN PetscErrorCode PetscStrlcat(char[],const char[],size_t);
1379: PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1380: PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1381: PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1382: PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1383: PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1384: PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1385: PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1386: PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1387: PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*);
1388: PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1389: PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1390: PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1391: PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1392: PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt,const char *const*,char***);
1393: PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt,char***);
1394: PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);

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

1398: PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1399: PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1400: PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);

1402: PETSC_EXTERN PetscErrorCode PetscStrInList(const char[],const char[],char,PetscBool*);
1403: PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1404: PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);

1406: /*
1407:    These are MPI operations for MPI_Allreduce() etc
1408: */
1409: PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1410: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1411: PETSC_EXTERN MPI_Op MPIU_SUM;
1412: PETSC_EXTERN MPI_Op MPIU_MAX;
1413: PETSC_EXTERN MPI_Op MPIU_MIN;
1414: #else
1415: #define MPIU_SUM MPI_SUM
1416: #define MPIU_MAX MPI_MAX
1417: #define MPIU_MIN MPI_MIN
1418: #endif
1419: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);

1421: PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1422: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);

1424: PETSC_EXTERN const char *const PetscFileModes[];

1426: /*
1427:     Defines PETSc error handling.
1428: */
1429: #include <petscerror.h>

1431: #define PETSC_SMALLEST_CLASSID  1211211
1432: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1433: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1434: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);
1435: PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject,PetscObjectId*);
1436: PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject,PetscObjectId,PetscBool*);

1438: /*
1439:    Routines that get memory usage information from the OS
1440: */
1441: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1442: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1443: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1444: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);

1446: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);

1448: /*
1449:    Initialization of PETSc
1450: */
1451: PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1452: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1453: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1454: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1455: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1456: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1457: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1458: PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1459: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1460: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);

1462: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1463: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);

1465: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1466: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1467: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1468: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);

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

1472: /*
1473:      These are so that in extern C code we can caste function pointers to non-extern C
1474:    function pointers. Since the regular C++ code expects its function pointers to be C++
1475: */
1476: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1477: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1478: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);

1480: /*
1481:     Functions that can act on any PETSc object.
1482: */
1483: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1484: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1485: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1486: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1487: PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1488: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1489: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1490: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1491: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1492: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1493: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1494: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1495: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1496: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1497: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt*);
1498: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1499: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1500: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject*);
1501: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void));
1502: #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d))
1503: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1504: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1505: PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1506: PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject,PetscObject);
1507: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt*);

1509: #include <petscviewertypes.h>
1510: #include <petscoptions.h>

1512: PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer,PetscBool,PetscLogDouble);
1513: PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool*);

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

1517: PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer,const char[]);
1518: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1519: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1520: #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1521: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1522: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1523: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1524: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1525: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1526: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1527: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1528: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1529: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,PetscObject,const char[]);
1530: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1531: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1532: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject,const char[],PetscBool *);
1533: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1534: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1535: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1536: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);

1538: #if defined(PETSC_HAVE_SAWS)
1539: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1540: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1541: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1542: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1543: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1544: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1545: PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1546: PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1547: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1548: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);

1550: #else
1551: #define PetscSAWsBlock()                        0
1552: #define PetscObjectSAWsViewOff(obj)             0
1553: #define PetscObjectSAWsSetBlock(obj,flg)        0
1554: #define PetscObjectSAWsBlock(obj)               0
1555: #define PetscObjectSAWsGrantAccess(obj)         0
1556: #define PetscObjectSAWsTakeAccess(obj)          0
1557: #define PetscStackViewSAWs()                    0
1558: #define PetscStackSAWsViewOff()                 0
1559: #define PetscStackSAWsTakeAccess()
1560: #define PetscStackSAWsGrantAccess()

1562: #endif

1564: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1565: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1566: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);
1567: PETSC_EXTERN PetscErrorCode PetscDLAddr(void (*)(void), char **);
1568: #ifdef PETSC_HAVE_CXX
1569: PETSC_EXTERN PetscErrorCode PetscDemangleSymbol(const char *, char **);
1570: #endif

1572: #if defined(PETSC_USE_DEBUG)
1573: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1574: #endif
1575: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);

1577: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1578: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1579: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1580: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1581: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1582: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);

1584: /*
1585:     Dynamic library lists. Lists of names of routines in objects or in dynamic
1586:   link libraries that will be loaded as needed.
1587: */

1589: #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1590: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1591: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1592: #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1593: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1594: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[],const char[]);
1595: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1596: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1597: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);

1599: PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1600: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1601: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1602: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1603: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1604: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1605: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1606: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);

1608: /*
1609:      Useful utility routines
1610: */
1611: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1612: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1613: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm,PetscInt*,PetscInt*);
1614: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1615: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1616: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1617: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);
1618: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm,PetscInt[2],PetscInt[2]);
1619: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm,PetscReal[2],PetscReal[2]);

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

1624:     Notes:
1625:     This is useful in cases like
1626: $     int        *a;
1627: $     PetscBool  flag = PetscNot(a)
1628:      where !a would not return a PetscBool because we cannot provide a cast from int to PetscBool in C.

1630:     Level: beginner

1632:     .seealso : PetscBool, PETSC_TRUE, PETSC_FALSE
1633: M*/
1634: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)

1636: /*MC
1637:     PetscHelpPrintf - Prints help messages.

1639:    Synopsis:
1640: #include <petscsys.h>
1641:      PetscErrorCode (*PetscHelpPrintf)(MPI_Comm comm, const char format[],args);

1643:     Collective on comm

1645:     Input Parameters:
1646: +  comm - the MPI communicator over which the help message is printed
1647: .  format - the usual printf() format string
1648: -  args - arguments to be printed

1650:    Level: developer

1652:    Fortran Note:
1653:      This routine is not supported in Fortran.

1655:    Note:
1656:      You can change how help messages are printed by replacing the function pointer with a function that does not simply write to stdout.

1658:       To use, write your own function, for example,
1659: $PetscErrorCode mypetschelpprintf(MPI_Comm comm,const char format[],....)
1660: ${
1661: $ return(0);
1662: $}
1663: then do the assigment
1664: $    PetscHelpPrintf = mypetschelpprintf;
1665:    You can do the assignment before PetscInitialize().

1667:   The default routine used is called PetscHelpPrintfDefault().

1669: .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1670: M*/
1671: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);

1673: /*
1674:      Defines PETSc profiling.
1675: */
1676: #include <petsclog.h>

1678: /*
1679:       Simple PETSc parallel IO for ASCII printing
1680: */
1681: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1682: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1683: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1684: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1685: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1686: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1687: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);
1688: PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[],size_t,const char*,PetscInt,const PetscReal[]);

1690: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1691: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1692: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);

1694: PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char*,size_t*);
1695: PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char*,char *);

1697: #if defined(PETSC_HAVE_POPEN)
1698: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1699: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*);
1700: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1701: #endif

1703: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1704: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1705: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*);
1706: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1707: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1708: PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1709: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);

1711: PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1712: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void**);
1713: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void*);
1714: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1715: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer*);
1716: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer,PetscErrorCode (*)(void*));
1717: PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void*);

1719: /*
1720:    For use in debuggers
1721: */
1722: PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1723: PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1724: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1725: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1726: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);

1728: #include <stddef.h>
1729: #include <string.h>             /* for memcpy, memset */
1730: #include <stdlib.h>

1732: #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1733: #include <xmmintrin.h>
1734: #endif

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

1741:    Not Collective

1743:    Input Parameters:
1744: +  b - pointer to initial memory space
1745: .  a - pointer to copy space
1746: -  n - length (in bytes) of space to copy

1748:    Level: intermediate

1750:    Note:
1751:    PetscArraymove() is preferred
1752:    This routine is analogous to memmove().

1754:    Developers Note: This is inlined for performance

1756: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscStrallocpy(),
1757:           PetscArraymove()
1758: @*/
1759: PETSC_STATIC_INLINE PetscErrorCode PetscMemmove(void *a,const void *b,size_t n)
1760: {
1762:   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to null pointer");
1763:   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1764: #if !defined(PETSC_HAVE_MEMMOVE)
1765:   if (a < b) {
1766:     if (a <= b - n) memcpy(a,b,n);
1767:     else {
1768:       memcpy(a,b,(int)(b - a));
1769:       PetscMemmove(b,b + (int)(b - a),n - (int)(b - a));
1770:     }
1771:   } else {
1772:     if (b <= a - n) memcpy(a,b,n);
1773:     else {
1774:       memcpy(b + n,b + (n - (int)(a - b)),(int)(a - b));
1775:       PetscMemmove(a,b,n - (int)(a - b));
1776:     }
1777:   }
1778: #else
1779:   memmove((char*)(a),(char*)(b),n);
1780: #endif
1781:   return(0);
1782: }

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

1789:    Not Collective

1791:    Input Parameters:
1792: +  b - pointer to initial memory space
1793: -  n - length (in bytes) of space to copy

1795:    Output Parameter:
1796: .  a - pointer to copy space

1798:    Level: intermediate

1800:    Compile Option:
1801:     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1802:                                   for memory copies on double precision values.
1803:     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1804:                                   for memory copies on double precision values.
1805:     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1806:                                   for memory copies on double precision values.

1808:    Note:
1809:    Prefer PetscArraycpy()
1810:    This routine is analogous to memcpy().
1811:    Not available from Fortran

1813:    Developer Note: this is inlined for fastest performance

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

1817: @*/
1818: PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n)
1819: {
1820: #if defined(PETSC_USE_DEBUG)
1821:   size_t al = (size_t) a,bl = (size_t) b;
1822:   size_t nl = (size_t) n;
1824:   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1825:   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1826: #else
1828: #endif
1829:   if (a != b && n > 0) {
1830: #if defined(PETSC_USE_DEBUG)
1831:     if ((al > bl && (al - bl) < nl) || (bl - al) < nl)  SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1832:               or make sure your copy regions and lengths are correct. \n\
1833:               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1834: #endif
1835: #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1836:    if (!(a % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1837:       size_t len = n/sizeof(PetscScalar);
1838: #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1839:       PetscBLASInt   one = 1,blen;
1841:       PetscBLASIntCast(len,&blen);
1842:       PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one));
1843: #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1844:       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1845: #else
1846:       size_t      i;
1847:       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1848:       for (i=0; i<len; i++) y[i] = x[i];
1849: #endif
1850:     } else {
1851:       memcpy((char*)(a),(char*)(b),n);
1852:     }
1853: #else
1854:     memcpy((char*)(a),(char*)(b),n);
1855: #endif
1856:   }
1857:   return(0);
1858: }

1860: /*@C
1861:    PetscMemzero - Zeros the specified memory.

1863:    Not Collective

1865:    Input Parameters:
1866: +  a - pointer to beginning memory location
1867: -  n - length (in bytes) of memory to initialize

1869:    Level: intermediate

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

1875:    Not available from Fortran
1876:    Prefer PetscArrayzero()

1878:    Developer Note: this is inlined for fastest performance

1880: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy()
1881: @*/
1882: PETSC_STATIC_INLINE PetscErrorCode PetscMemzero(void *a,size_t n)
1883: {
1884:   if (n > 0) {
1885: #if defined(PETSC_USE_DEBUG)
1886:     if (!a) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer with %zu bytes",n);
1887: #endif
1888: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1889:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1890:       size_t      i,len = n/sizeof(PetscScalar);
1891:       PetscScalar *x = (PetscScalar*)a;
1892:       for (i=0; i<len; i++) x[i] = 0.0;
1893:     } else {
1894: #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1895:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1896:       PetscInt len = n/sizeof(PetscScalar);
1897:       fortranzero_(&len,(PetscScalar*)a);
1898:     } else {
1899: #endif
1900: #if defined(PETSC_PREFER_BZERO)
1901:       bzero((char *)a,n);
1902: #else
1903:       memset((char*)a,0,n);
1904: #endif
1905: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1906:     }
1907: #endif
1908:   }
1909:   return 0;
1910: }

1912: /*MC
1913:    PetscArraycmp - Compares two arrays in memory.

1915:    Synopsis:
1916: #include <petscsys.h>
1917:     PetscErrorCode PetscArraycmp(const anytype *str1,const anytype *str2,size_t cnt,PetscBool *e)

1919:    Not Collective

1921:    Input Parameters:
1922: +  str1 - First array
1923: .  str2 - Second array
1924: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1926:    Output Parameters:
1927: .   e - PETSC_TRUE if equal else PETSC_FALSE.

1929:    Level: intermediate

1931:    Note:
1932:    This routine is a preferred replacement to PetscMemcmp()
1933:    The arrays must be of the same type

1935: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy(),
1936:           PetscArraymove()
1937: M*/
1938: #define  PetscArraycmp(str1,str2,cnt,e) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcmp(str1,str2,(size_t)(cnt)*sizeof(*(str1)),e))

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

1944:    Synopsis:
1945: #include <petscsys.h>
1946:     PetscErrorCode PetscArraymove(anytype *str1,const anytype *str2,size_t cnt)

1948:    Not Collective

1950:    Input Parameters:
1951: +  str1 - First array
1952: .  str2 - Second array
1953: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1955:    Level: intermediate

1957:    Note:
1958:    This routine is a preferred replacement to PetscMemmove()
1959:    The arrays must be of the same type

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

1965: /*MC
1966:    PetscArraycpy - Copies from one array in memory to another

1968:    Synopsis:
1969: #include <petscsys.h>
1970:     PetscErrorCode PetscArraycpy(anytype *str1,const anytype *str2,size_t cnt)

1972:    Not Collective

1974:    Input Parameters:
1975: +  str1 - First array (destination)
1976: .  str2 - Second array (source)
1977: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1979:    Level: intermediate

1981:    Note:
1982:    This routine is a preferred replacement to PetscMemcpy()
1983:    The arrays must be of the same type

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

1989: /*MC
1990:    PetscArrayzero - Zeros an array in memory.

1992:    Synopsis:
1993: #include <petscsys.h>
1994:     PetscErrorCode PetscArrayzero(anytype *str1,size_t cnt)

1996:    Not Collective

1998:    Input Parameters:
1999: +  str1 - array
2000: -  cnt  - Count of the array, not in bytes, but number of entries in the array

2002:    Level: intermediate

2004:    Note:
2005:    This routine is a preferred replacement to PetscMemzero()

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

2011: /*MC
2012:    PetscPrefetchBlock - Prefetches a block of memory

2014:    Synopsis:
2015: #include <petscsys.h>
2016:     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)

2018:    Not Collective

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

2026:    Level: developer

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

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

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

2041: M*/
2042: #define PetscPrefetchBlock(a,n,rw,t) do {                               \
2043:     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
2044:     for (; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
2045:   } while (0)

2047: /*
2048:       Determine if some of the kernel computation routines use
2049:    Fortran (rather than C) for the numerical calculations. On some machines
2050:    and compilers (like complex numbers) the Fortran version of the routines
2051:    is faster than the C/C++ versions. The flag --with-fortran-kernels
2052:    should be used with ./configure to turn these on.
2053: */
2054: #if defined(PETSC_USE_FORTRAN_KERNELS)

2056: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
2057: #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
2058: #endif

2060: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
2061: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
2062: #endif

2064: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
2065: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
2066: #endif

2068: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
2069: #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
2070: #endif

2072: #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
2073: #define PETSC_USE_FORTRAN_KERNEL_NORM
2074: #endif

2076: #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
2077: #define PETSC_USE_FORTRAN_KERNEL_MAXPY
2078: #endif

2080: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
2081: #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
2082: #endif

2084: #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
2085: #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
2086: #endif

2088: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
2089: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
2090: #endif

2092: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
2093: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
2094: #endif

2096: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
2097: #define PETSC_USE_FORTRAN_KERNEL_MDOT
2098: #endif

2100: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
2101: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
2102: #endif

2104: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
2105: #define PETSC_USE_FORTRAN_KERNEL_AYPX
2106: #endif

2108: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
2109: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
2110: #endif

2112: #endif

2114: /*
2115:     Macros for indicating code that should be compiled with a C interface,
2116:    rather than a C++ interface. Any routines that are dynamically loaded
2117:    (such as the PCCreate_XXX() routines) must be wrapped so that the name
2118:    mangler does not change the functions symbol name. This just hides the
2119:    ugly extern "C" {} wrappers.
2120: */
2121: #if defined(__cplusplus)
2122: #  define EXTERN_C_BEGIN extern "C" {
2123: #  define EXTERN_C_END }
2124: #else
2125: #  define EXTERN_C_BEGIN
2126: #  define EXTERN_C_END
2127: #endif

2129: /* --------------------------------------------------------------------*/

2131: /*MC
2132:     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
2133:         communication

2135:    Level: beginner

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

2139: .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
2140: M*/

2142: #if defined(PETSC_HAVE_MPIIO)
2143: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2144: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2145: PETSC_EXTERN PetscErrorCode MPIU_File_write_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2146: PETSC_EXTERN PetscErrorCode MPIU_File_read_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2147: PETSC_EXTERN PetscErrorCode MPIU_File_write_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2148: PETSC_EXTERN PetscErrorCode MPIU_File_read_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2149: #endif

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

2153: /* Limit MPI to 32-bits */
2154: #define PETSC_MPI_INT_MAX  2147483647
2155: #define PETSC_MPI_INT_MIN -2147483647
2156: /* Limit BLAS to 32-bits */
2157: #define PETSC_BLAS_INT_MAX  2147483647
2158: #define PETSC_BLAS_INT_MIN -2147483647
2159: #define PETSC_CUBLAS_INT_MAX  2147483647

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

2165:    Not Collective

2167:    Input Parameter:
2168: .     a - the PetscInt64 value

2170:    Output Parameter:
2171: .     b - the resulting PetscInt value

2173:    Level: advanced

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

2177:    Not available from Fortran

2179: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscBLASIntCast(), PetscIntMultError(), PetscIntSumError()
2180: @*/
2181: PETSC_STATIC_INLINE PetscErrorCode PetscIntCast(PetscInt64 a,PetscInt *b)
2182: {
2184: #if !defined(PETSC_USE_64BIT_INDICES)
2185:   if (a > PETSC_MAX_INT) { *b = 0; SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"%D is too big for PetscInt, you may need to ./configure using --with-64-bit-indices",a); }
2186: #endif
2187:   *b = (PetscInt)(a);
2188:   return(0);
2189: }

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

2195:    Not Collective

2197:    Input Parameter:
2198: .     a - the PetscInt value

2200:    Output Parameter:
2201: .     b - the resulting PetscBLASInt value

2203:    Level: advanced

2205:    Notes:
2206:       Not available from Fortran
2207:       Errors if the integer is negative since PETSc calls to BLAS/LAPACK never need to cast negative integer inputs

2209: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscIntCast()
2210: @*/
2211: PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
2212: {
2214: #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
2215:   if (a > PETSC_BLAS_INT_MAX) { *b = 0; SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"%D is too big for BLAS/LAPACK, which is restricted to 32 bit integers. Either you have an invalidly large integer error in your code or you must ./configure PETSc with -with-64-bit-blas-indices for the case you are running",a); }
2216: #endif
2217:   if (a < 0) { *b = 0; SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Passing negative integer to BLAS/LAPACK routine"); }
2218:   *b = (PetscBLASInt)(a);
2219:   return(0);
2220: }

2222: /*@C
2223:     PetscCuBLASIntCast - like PetscBLASIntCast(), but for PetscCuBLASInt.

2225:    Not Collective

2227:    Input Parameter:
2228: .     a - the PetscInt value

2230:    Output Parameter:
2231: .     b - the resulting PetscCuBLASInt value

2233:    Level: advanced

2235:    Notes:
2236:       Errors if the integer is negative since PETSc calls to cuBLAS and friends never need to cast negative integer inputs

2238: .seealso: PetscCuBLASInt, PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscMPIIntCast(), PetscIntCast()
2239: @*/
2240: PETSC_STATIC_INLINE PetscErrorCode PetscCuBLASIntCast(PetscInt a,PetscCuBLASInt *b)
2241: {
2243: #if defined(PETSC_USE_64BIT_INDICES)
2244:   if (a > PETSC_CUBLAS_INT_MAX) { *b = 0; SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"%D is too big for cuBLAS, which is restricted to 32 bit integers.",a); }
2245: #endif
2246:   if (a < 0) { *b = 0; SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Passing negative integer to cuBLAS routine"); }
2247:   *b = (PetscCuBLASInt)(a);
2248:   return(0);
2249: }

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

2255:    Not Collective

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

2260:    Output Parameter:
2261: .     b - the resulting PetscMPIInt value

2263:    Level: advanced

2265:    Not available from Fortran

2267: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntCast()
2268: @*/
2269: PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
2270: {
2272: #if defined(PETSC_USE_64BIT_INDICES)
2273:   if (a > PETSC_MPI_INT_MAX) { *b = 0; SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"%D is too big for MPI buffer length. We currently only support 32 bit integers",a); }
2274: #endif
2275:   *b = (PetscMPIInt)(a);
2276:   return(0);
2277: }

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

2281: /*@C

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

2285:    Not Collective

2287:    Input Parameters:
2288: +     a - the PetscReal value
2289: -     b - the second value

2291:    Returns:
2292:       the result as a PetscInt value

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

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

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

2303:    Not available from Fortran

2305:    Level: advanced

2307: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError(), PetscIntSumError()
2308: @*/
2309: PETSC_STATIC_INLINE PetscInt PetscRealIntMultTruncate(PetscReal a,PetscInt b)
2310: {
2311:   PetscInt64 r;

2313:   r  =  (PetscInt64) (a*(PetscReal)b);
2314:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2315:   return (PetscInt) r;
2316: }

2318: /*@C

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

2322:    Not Collective

2324:    Input Parameters:
2325: +     a - the PetscInt value
2326: -     b - the second value

2328:    Returns:
2329:       the result as a PetscInt value

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

2335:    Not available from Fortran

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

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

2341:    Level: advanced

2343: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError(), PetscIntSumError()
2344: @*/
2345: PETSC_STATIC_INLINE PetscInt PetscIntMultTruncate(PetscInt a,PetscInt b)
2346: {
2347:   PetscInt64 r;

2349:   r  =  PetscInt64Mult(a,b);
2350:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2351:   return (PetscInt) r;
2352: }

2354: /*@C

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

2358:    Not Collective

2360:    Input Parameters:
2361: +     a - the PetscInt value
2362: -     b - the second value

2364:    Returns:
2365:      the result as a PetscInt value

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

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

2373:    Not available from Fortran

2375:    Level: advanced

2377: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError()
2378: @*/
2379: PETSC_STATIC_INLINE PetscInt PetscIntSumTruncate(PetscInt a,PetscInt b)
2380: {
2381:   PetscInt64 r;

2383:   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2384:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2385:   return (PetscInt) r;
2386: }

2388: /*@C

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

2392:    Not Collective

2394:    Input Parameters:
2395: +     a - the PetscInt value
2396: -     b - the second value

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

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

2404:    Not available from Fortran

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

2408:    Level: advanced

2410: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64(), PetscIntSumError()
2411: @*/
2412: PETSC_STATIC_INLINE PetscErrorCode PetscIntMultError(PetscInt a,PetscInt b,PetscInt *result)
2413: {
2414:   PetscInt64 r;

2417:   r  =  PetscInt64Mult(a,b);
2418: #if !defined(PETSC_USE_64BIT_INDICES)
2419:   if (r > PETSC_MAX_INT) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Product of two integers %d %d overflow, either you have an invalidly large integer error in your code or you must ./configure PETSc with --with-64-bit-indices for the case you are running",a,b);
2420: #endif
2421:   if (result) *result = (PetscInt) r;
2422:   return(0);
2423: }

2425: /*@C

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

2429:    Not Collective

2431:    Input Parameters:
2432: +     a - the PetscInt value
2433: -     b - the second value

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

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

2441:    Not available from Fortran

2443:    Level: advanced

2445: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError()
2446: @*/
2447: PETSC_STATIC_INLINE PetscErrorCode PetscIntSumError(PetscInt a,PetscInt b,PetscInt *result)
2448: {
2449:   PetscInt64 r;

2452:   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2453: #if !defined(PETSC_USE_64BIT_INDICES)
2454:   if (r > PETSC_MAX_INT) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Sum of two integers %d %d overflow, either you have an invalidly large integer error in your code or you must ./configure PETSc with --with-64-bit-indices for the case you are running",a,b);
2455: #endif
2456:   if (result) *result = (PetscInt) r;
2457:   return(0);
2458: }

2460: /*
2461:      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2462: */
2463: #if defined(hz)
2464: #  undef hz
2465: #endif

2467: #include <limits.h>

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

2471: #define PETSC_BITS_PER_BYTE CHAR_BIT

2473: #if defined(PETSC_HAVE_SYS_TYPES_H)
2474: #  include <sys/types.h>
2475: #endif

2477: /*MC

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

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

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

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

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

2491:     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),
2492:     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
2493:     version number (x.y.z).

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

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

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

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

2503:     Level: intermediate

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

2507: M*/

2509: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2510: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2511: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2512: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2513: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2514: PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2515: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2516: PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt*,PetscInt*,PetscInt*,PetscInt*);

2518: PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt,const PetscInt[],PetscBool*);
2519: PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt,const PetscMPIInt[],PetscBool*);
2520: PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt,const PetscReal[],PetscBool*);
2521: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2522: PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt,PetscInt[]);
2523: PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt*,PetscInt[]);
2524: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2526: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2527: PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt*);
2528: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2529: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2530: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2531: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt[],PetscInt[],PetscInt[]);
2532: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2533: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2534: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2535: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt,PetscMPIInt[],PetscInt[]);
2536: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2537: PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt,PetscInt[],void*,size_t,void*);
2538: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2539: PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt,PetscReal[],PetscInt[]);
2540: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2541: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2542: PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal,PetscInt,const PetscReal[],PetscReal,PetscInt*);
2543: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2544: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2545: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2546: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt[],const PetscInt[],PetscInt,const PetscInt[],const PetscInt[],PetscInt*,PetscInt**,PetscInt**);
2547: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt*,PetscInt**);
2548: PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt,const PetscMPIInt[],PetscInt,const PetscMPIInt[],PetscInt*,PetscMPIInt**);
2549: PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *);

2551: PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt,void*,size_t,int(*)(const void*,const void*,void*),void*);
2552: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt,PetscInt[]);
2553: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt,PetscMPIInt[]);
2554: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt,PetscReal[]);
2555: PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt,void*,size_t,void*,size_t,int(*)(const void*,const void*,void*),void*);
2556: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt,PetscInt[],PetscInt[]);
2557: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt,PetscMPIInt[],PetscMPIInt[]);
2558: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt,PetscReal[],PetscInt[]);

2560: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2561: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);

2563: /*J
2564:     PetscRandomType - String with the name of a PETSc randomizer

2566:    Level: beginner

2568:    Notes:
2569:    To use SPRNG or RANDOM123 you must have ./configure PETSc
2570:    with the option --download-sprng or --download-random123

2572: .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2573: J*/
2574: typedef const char* PetscRandomType;
2575: #define PETSCRAND       "rand"
2576: #define PETSCRAND48     "rand48"
2577: #define PETSCSPRNG      "sprng"
2578: #define PETSCRANDER48   "rander48"
2579: #define PETSCRANDOM123  "random123"
2580: #define PETSCCURAND     "curand"

2582: /* Logging support */
2583: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;

2585: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);

2587: /* Dynamic creation and loading functions */
2588: PETSC_EXTERN PetscFunctionList PetscRandomList;

2590: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2591: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom,PetscRandomType);
2592: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2593: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom,PetscRandomType*);
2594: PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom,PetscObject,const char[]);
2595: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);

2597: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2598: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2599: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2600: PETSC_EXTERN PetscErrorCode PetscRandomGetValues(PetscRandom,PetscInt,PetscScalar*);
2601: PETSC_EXTERN PetscErrorCode PetscRandomGetValuesReal(PetscRandom,PetscInt,PetscReal*);
2602: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2603: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2604: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2605: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2606: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2607: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);

2609: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2610: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2611: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2612: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2613: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2614: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2615: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);
2616: PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2617: PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2618: PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);

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

2622: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscInt*,PetscDataType);
2623: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscInt*,PetscDataType);
2624: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,const void*,PetscInt,PetscDataType);
2625: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,const void*,PetscInt,PetscDataType);
2626: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2627: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2628: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2629: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2630: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2631: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2632: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2633: #if defined(PETSC_USE_SOCKET_VIEWER)
2634: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);
2635: #endif

2637: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2638: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2639: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);

2641: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2642: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool);
2643: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2644: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2645: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2646: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2647: PETSC_EXTERN PetscErrorCode PetscWaitOnError(void);

2649: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2650: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2651: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2652: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2653: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2654: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt*,const void*,PetscMPIInt*,PetscMPIInt**,void*)
2655:   PetscAttrMPIPointerWithType(6,3);
2656: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2657:                                                     PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2658:                                                     PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2659:   PetscAttrMPIPointerWithType(6,3);
2660: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2661:                                                        MPI_Request**,MPI_Request**,
2662:                                                        PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2663:                                                        PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2664:   PetscAttrMPIPointerWithType(6,3);

2666: PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2667: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2668: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);

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

2672: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);

2674: PETSC_EXTERN const char *const PetscSubcommTypes[];

2676: struct _n_PetscSubcomm {
2677:   MPI_Comm         parent;           /* parent communicator */
2678:   MPI_Comm         dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2679:   MPI_Comm         child;            /* the sub-communicator */
2680:   PetscMPIInt      n;                /* num of subcommunicators under the parent communicator */
2681:   PetscMPIInt      color;            /* color of processors belong to this communicator */
2682:   PetscMPIInt      *subsize;         /* size of subcommunicator[color] */
2683:   PetscSubcommType type;
2684:   char             *subcommprefix;
2685: };

2687: PETSC_STATIC_INLINE MPI_Comm PetscSubcommParent(PetscSubcomm scomm) {return scomm->parent;}
2688: PETSC_STATIC_INLINE MPI_Comm PetscSubcommChild(PetscSubcomm scomm) {return scomm->child;}
2689: PETSC_STATIC_INLINE MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm) {return scomm->dupparent;}
2690: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2691: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2692: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2693: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2694: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2695: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2696: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2697: PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm,const char[]);
2698: PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm,MPI_Comm*);
2699: PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm,MPI_Comm*);
2700: PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm,MPI_Comm*);

2702: PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt,PetscHeap*);
2703: PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap,PetscInt,PetscInt);
2704: PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap,PetscInt*,PetscInt*);
2705: PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap,PetscInt*,PetscInt*);
2706: PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap,PetscInt,PetscInt);
2707: PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2708: PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap*);
2709: PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap,PetscViewer);

2711: PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2712: PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm,PetscShmComm*);
2713: PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2714: PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2715: PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm,MPI_Comm*);

2717: /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2718: PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm,PetscInt,PetscOmpCtrl*);
2719: PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl,MPI_Comm*,MPI_Comm*,PetscBool*);
2720: PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl*);
2721: PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2722: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2723: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);

2725: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2726: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2727: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2728: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2729: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2730: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2731: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2732: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);

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

2739: extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2740: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted*);
2741: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted*);
2742: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted,const char*,const char*,PetscBool*);

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

2748: PETSC_EXTERN PetscSegBuffer PetscCitationsList;

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

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

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

2759:    Level: intermediate

2761:    Not available from Fortran

2763:      Options Database:
2764: .     -citations [filename]   - print out the bibtex entries for the given computation
2765: @*/
2766: PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2767: {
2768:   size_t         len;
2769:   char           *vstring;

2773:   if (set && *set) return(0);
2774:   PetscStrlen(cit,&len);
2775:   PetscSegBufferGet(PetscCitationsList,len,&vstring);
2776:   PetscArraycpy(vstring,cit,len);
2777:   if (set) *set = PETSC_TRUE;
2778:   return(0);
2779: }

2781: PETSC_EXTERN PetscErrorCode PetscURLShorten(const char[],char[],size_t);
2782: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm,char[],char[],size_t);
2783: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm,const char[],char[],size_t);
2784: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm,const char[],const char []);

2786: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm,char[],char[],size_t);
2787: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm,const char[],char[],char[],size_t);

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

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

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

2797: #if defined(PETSC_USE_DEBUG)
2798: PETSC_STATIC_INLINE unsigned int PetscStrHash(const char *str)
2799: {
2800:   unsigned int c,hash = 5381;

2802:   while ((c = (unsigned int)*str++)) hash = ((hash << 5) + hash) + c; /* hash * 33 + c */
2803:   return hash;
2804: }

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

2811:    Collective

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

2817:    Input Parameters:
2818: +  a - pointer to the input data to be reduced
2819: .  c - the number of MPI data items in a and b
2820: .  d - the MPI datatype, for example MPI_INT
2821: .  e - the MPI operation, for example MPI_SUM
2822: -  fcomm - the MPI communicator on which the operation occurs

2824:    Output Parameter:
2825: .  b - the reduced values

2827:    Notes:
2828:      In optimized mode this directly calls MPI_Allreduce()

2830:      This is defined as a macro that can return error codes internally so it cannot be used in a subroutine that returns void.

2832:      The error code this returns should be checked with CHKERRMPI()

2834:    Level: developer

2836: .seealso: MPI_Allreduce()
2837: M*/
2838: #define MPIU_Allreduce(a,b,c,d,e,fcomm) MPI_SUCCESS; do {\
2839:   PetscErrorCode _4_ierr; \
2840:   PetscMPIInt a_b1[6],a_b2[6];\
2841:   a_b1[0] = -(PetscMPIInt)__LINE__;                          a_b1[1] = -a_b1[0];\
2842:   a_b1[2] = -(PetscMPIInt)PetscStrHash(PETSC_FUNCTION_NAME); a_b1[3] = -a_b1[2];\
2843:   a_b1[4] = -(PetscMPIInt)(c);                               a_b1[5] = -a_b1[4];\
2844:   _4_MPI_Allreduce(a_b1,a_b2,6,MPI_INT,MPI_MAX,fcomm);CHKERRMPI(_4_ierr);\
2845:   if (-a_b2[0] != a_b2[1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_Allreduce() called in different locations (code lines) on different processors");\
2846:   if (-a_b2[2] != a_b2[3]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_Allreduce() called in different locations (functions) on different processors");\
2847:   if (-a_b2[4] != a_b2[5]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_Allreduce() called with different counts %d on different processors",c);\
2848:   _4_MPI_Allreduce((a),(b),(c),d,e,(fcomm));CHKERRMPI(_4_ierr);\
2849:   } while (0)
2850: #else
2851: #define MPIU_Allreduce(a,b,c,d,e,fcomm) MPI_Allreduce((a),(b),(c),d,e,(fcomm))
2852: #endif

2854: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
2855: PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint,PetscMPIInt,MPI_Info,MPI_Comm,void*,MPI_Win*);
2856: PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win,PetscMPIInt,MPI_Aint*,PetscMPIInt*,void*);
2857: #endif

2859: /* this is a vile hack */
2860: #if defined(PETSC_HAVE_NECMPI)
2861: #define MPI_Type_free(a) (*(a) = MPI_DATATYPE_NULL,0);
2862: #endif

2864: /*
2865:     List of external packages and queries on it
2866: */
2867: PETSC_EXTERN PetscErrorCode  PetscHasExternalPackage(const char[],PetscBool*);

2869: #endif