Actual source code: petscsys.h

petsc-3.14.6 2021-03-30
Report Typos and Errors
  1: /*
  2:    This is the main PETSc include file (for C and C++).  It is included by all
  3:    other PETSc include files, so it almost never has to be specifically included.
  4: */
  5: #if !defined(PETSCSYS_H)
  6: #define PETSCSYS_H
  7: /* ========================================================================== */
  8: /*
  9:    petscconf.h is contained in ${PETSC_ARCH}/include/petscconf.h it is
 10:    found automatically by the compiler due to the -I${PETSC_DIR}/${PETSC_ARCH}/include.
 11:    For --prefix installs the ${PETSC_ARCH}/ does not exist and petscconf.h is in the same
 12:    directory as the other PETSc include files.
 13: */
 14: #include <petscconf.h>
 15: #include <petscfix.h>

 17: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP) || defined(PETSC_HAVE_KOKKOS)
 18:    #define PETSC_HAVE_DEVICE
 19: #endif

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

 40: #include <petscsystypes.h>

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

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

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

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

 74: #define PETSC_STATIC_INLINE static PETSC_INLINE

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

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

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

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

118: /* ========================================================================== */

120: /*
121:     Defines the interface to MPI allowing the use of all MPI functions.

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

141: /*MC
142:   PetscDefined - determine whether a boolean macro is defined

144:   Notes:
145:   The prefix "PETSC_" is added to the argument.

147:   Typical usage is within normal code,

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

151:   but can also be used in the preprocessor,

153: $   #if PetscDefined(USE_DEBUG)
154: $     ...
155: $   #else

157:   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
158:   should not be used if its argument may be defined to a non-empty value other than 1.

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

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

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

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

173:   Level: developer
174: M*/
175: #if !defined(PETSC_SKIP_VARIADIC_MACROS)
176: #  define PetscDefined_arg_1    shift,
177: #  define PetscDefined_arg_     shift,
178: #  define PetscDefined__take_second_expanded(ignored, val, ...) val
179: #  define PetscDefined__take_second_expand(args) PetscDefined__take_second_expanded args
180: #  define PetscDefined__take_second(...) PetscDefined__take_second_expand((__VA_ARGS__))
181: #  define PetscDefined___(arg1_or_junk) PetscDefined__take_second(arg1_or_junk 1, 0, at_)
182: #  define PetscDefined__(value) PetscDefined___(PetscDefined_arg_ ## value)
183: #  define PetscDefined_(d)      PetscDefined__(d)
184: #  define PetscDefined(d)       PetscDefined_(PETSC_ ## d)
185: #endif

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

231: /*
232:     Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler
233:     see the top of mpicxx.h in the MPICH2 distribution.
234: */
235: #include <stdio.h>

237: /* MSMPI on 32bit windows requires this yukky hack - that breaks MPI standard compliance */
238: #if !defined(MPIAPI)
239: #define MPIAPI
240: #endif

242: /*
243:    Support for Clang (>=3.2) matching type tag arguments with void* buffer types.
244:    This allows the compiler to detect cases where the MPI datatype argument passed to a MPI routine
245:    does not match the actual type of the argument being passed in
246: */
247: #if defined(__has_attribute) && defined(works_with_const_which_is_not_true)
248: #  if __has_attribute(argument_with_type_tag) && __has_attribute(pointer_with_type_tag) && __has_attribute(type_tag_for_datatype)
249: #    define PetscAttrMPIPointerWithType(bufno,typeno) __attribute__((pointer_with_type_tag(MPI,bufno,typeno)))
250: #    define PetscAttrMPITypeTag(type)                 __attribute__((type_tag_for_datatype(MPI,type)))
251: #    define PetscAttrMPITypeTagLayoutCompatible(type) __attribute__((type_tag_for_datatype(MPI,type,layout_compatible)))
252: #  endif
253: #endif
254: #if !defined(PetscAttrMPIPointerWithType)
255: #  define PetscAttrMPIPointerWithType(bufno,typeno)
256: #  define PetscAttrMPITypeTag(type)
257: #  define PetscAttrMPITypeTagLayoutCompatible(type)
258: #endif

260: PETSC_EXTERN MPI_Datatype MPIU_ENUM PetscAttrMPITypeTag(PetscEnum);
261: PETSC_EXTERN MPI_Datatype MPIU_BOOL PetscAttrMPITypeTag(PetscBool);

263: /*MC
264:    MPIU_INT - MPI datatype corresponding to PetscInt

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

269:    Level: beginner

271: .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_SCALAR, MPIU_COMPLEX
272: M*/

274: #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 */
275: #  define MPIU_INT64 MPI_INT64_T
276: #  define PetscInt64_FMT PRId64
277: #elif (PETSC_SIZEOF_LONG_LONG == 8)
278: #  define MPIU_INT64 MPI_LONG_LONG_INT
279: #  define PetscInt64_FMT "lld"
280: #elif defined(PETSC_HAVE___INT64)
281: #  define MPIU_INT64 MPI_INT64_T
282: #  define PetscInt64_FMT "ld"
283: #else
284: #  error "cannot determine PetscInt64 type"
285: #endif

287: PETSC_EXTERN MPI_Datatype MPIU_FORTRANADDR;

289: #if defined(PETSC_USE_64BIT_INDICES)
290: #  define MPIU_INT MPIU_INT64
291: #  define PetscInt_FMT PetscInt64_FMT
292: #else
293: #  define MPIU_INT MPI_INT
294: #  define PetscInt_FMT "d"
295: #endif

297: /*
298:     For the rare cases when one needs to send a size_t object with MPI
299: */
300: PETSC_EXTERN MPI_Datatype MPIU_SIZE_T;

302: /*
303:       You can use PETSC_STDOUT as a replacement of stdout. You can also change
304:     the value of PETSC_STDOUT to redirect all standard output elsewhere
305: */
306: PETSC_EXTERN FILE* PETSC_STDOUT;

308: /*
309:       You can use PETSC_STDERR as a replacement of stderr. You can also change
310:     the value of PETSC_STDERR to redirect all standard error elsewhere
311: */
312: PETSC_EXTERN FILE* PETSC_STDERR;

314: /*MC
315:     PetscUnlikely - hints the compiler that the given condition is usually FALSE

317:     Synopsis:
318: #include <petscsys.h>
319:     PetscBool  PetscUnlikely(PetscBool  cond)

321:     Not Collective

323:     Input Parameters:
324: .   cond - condition or expression

326:     Notes:
327:     This returns the same truth value, it is only a hint to compilers that the resulting
328:     branch is unlikely.

330:     Level: advanced

332: .seealso: PetscUnlikelyDebug(), PetscLikely(), CHKERRQ
333: M*/

335: /*MC
336:     PetscLikely - hints the compiler that the given condition is usually TRUE

338:     Synopsis:
339: #include <petscsys.h>
340:     PetscBool  PetscLikely(PetscBool  cond)

342:     Not Collective

344:     Input Parameters:
345: .   cond - condition or expression

347:     Notes:
348:     This returns the same truth value, it is only a hint to compilers that the resulting
349:     branch is likely.

351:     Level: advanced

353: .seealso: PetscUnlikely()
354: M*/
355: #if defined(PETSC_HAVE_BUILTIN_EXPECT)
356: #  define PetscUnlikely(cond)   __builtin_expect(!!(cond),0)
357: #  define PetscLikely(cond)     __builtin_expect(!!(cond),1)
358: #else
359: #  define PetscUnlikely(cond)   (cond)
360: #  define PetscLikely(cond)     (cond)
361: #endif

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

366:     Synopsis:
367: #include <petscsys.h>
368:     PetscBool  PetscUnlikelyDebug(PetscBool  cond)

370:     Not Collective

372:     Input Parameters:
373: .   cond - condition or expression

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

379:     Level: advanced

381: .seealso: PetscUnlikely(), CHKERRQ, SETERRQ
382: M*/
383: #define PetscUnlikelyDebug(cond) (PetscDefined(USE_DEBUG) && PetscUnlikely(cond))

385: /* PetscPragmaSIMD - from CeedPragmaSIMD */

387: #if defined(__INTEL_COMPILER) && !defined(_WIN32)
388: #  define PetscPragmaSIMD _Pragma("vector")
389: #elif defined(__GNUC__) && __GNUC__ >= 5 && !defined(__PGI)
390: #  define PetscPragmaSIMD _Pragma("GCC ivdep")
391: #elif defined(_OPENMP) && _OPENMP >= 201307 && !defined(_WIN32)
392: #  define PetscPragmaSIMD _Pragma("omp simd")
393: #elif defined(_OPENMP) && _OPENMP >= 201307 && defined(_WIN32)
394: #  define PetscPragmaSIMD __pragma(omp simd)
395: #elif defined(PETSC_HAVE_CRAY_VECTOR)
396: #  define PetscPragmaSIMD _Pragma("_CRI ivdep")
397: #else
398: #  define PetscPragmaSIMD
399: #endif

401: /*
402:     Declare extern C stuff after including external header files
403: */

405: PETSC_EXTERN const char *const PetscBools[];

407: /*
408:     Defines elementary mathematics functions and constants.
409: */
410: #include <petscmath.h>

412: PETSC_EXTERN const char *const PetscCopyModes[];

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

417:    Level: beginner

419:    Note:
420:    Accepted by many PETSc functions to not set a parameter and instead use
421:           some default

423:    Fortran Notes:
424:    This macro does not exist in Fortran; you must use PETSC_NULL_INTEGER,
425:           PETSC_NULL_DOUBLE_PRECISION etc

427: .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_DETERMINE

429: M*/
430: #define PETSC_IGNORE NULL

432: /* This is deprecated */
433: #define PETSC_NULL NULL

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

439:    Level: beginner

441: .seealso: PETSC_DEFAULT, PETSC_IGNORE, PETSC_DETERMINE

443: M*/
444: #define PETSC_DECIDE -1

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

450:    Level: beginner

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

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

458: M*/
459: #define PETSC_DETERMINE PETSC_DECIDE

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

465:    Level: beginner

467:    Fortran Notes:
468:    You need to use PETSC_DEFAULT_INTEGER or PETSC_DEFAULT_REAL.

470: .seealso: PETSC_DECIDE, PETSC_IGNORE, PETSC_DETERMINE

472: M*/
473: #define PETSC_DEFAULT -2

475: /*MC
476:     PETSC_COMM_WORLD - the equivalent of the MPI_COMM_WORLD communicator which represents
477:            all the processes that PETSc knows about.

479:    Level: beginner

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

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

490: .seealso: PETSC_COMM_SELF

492: M*/
493: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;

495: /*MC
496:     PETSC_COMM_SELF - This is always MPI_COMM_SELF

498:    Level: beginner

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

503: .seealso: PETSC_COMM_WORLD

505: M*/
506: #define PETSC_COMM_SELF MPI_COMM_SELF

508: /*MC
509:     PETSC_MPI_THREAD_REQUIRED - the required threading support used if PETSc initializes
510:            MPI with MPI_Init_thread.

512:    Level: beginner

514:    Notes:
515:    By default PETSC_MPI_THREAD_REQUIRED equals MPI_THREAD_FUNNELED.

517: .seealso: PetscInitialize()

519: M*/
520: PETSC_EXTERN PetscMPIInt PETSC_MPI_THREAD_REQUIRED;

522: PETSC_EXTERN PetscBool PetscBeganMPI;
523: PETSC_EXTERN PetscBool PetscErrorHandlingInitialized;
524: PETSC_EXTERN PetscBool PetscInitializeCalled;
525: PETSC_EXTERN PetscBool PetscFinalizeCalled;
526: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;

528: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm));
529: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*);
530: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm*);

532: #if defined(PETSC_HAVE_CUDA)
533: PETSC_EXTERN PetscBool      PetscCUDASynchronize;
534: PETSC_EXTERN PetscErrorCode PetscCUDAInitialize(MPI_Comm,PetscInt);
535: PETSC_EXTERN PetscErrorCode PetscCUDAInitializeCheck(void);
536: #endif

538: #if defined(PETSC_HAVE_HIP)
539: PETSC_EXTERN PetscBool      PetscHIPSynchronize;
540: PETSC_EXTERN PetscErrorCode PetscHIPInitialize(MPI_Comm,PetscInt);
541: PETSC_EXTERN PetscErrorCode PetscHIPInitializeCheck(void);
542: #endif

544: #if defined(PETSC_HAVE_ELEMENTAL)
545: PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
546: PETSC_EXTERN PetscErrorCode PetscElementalInitialized(PetscBool*);
547: PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
548: #endif

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

553:    Synopsis:
554: #include <petscsys.h>
555:    PetscErrorCode PetscMalloc(size_t m,void **result)

557:    Not Collective

559:    Input Parameter:
560: .  m - number of bytes to allocate

562:    Output Parameter:
563: .  result - memory allocated

565:    Level: beginner

567:    Notes:
568:    Memory is always allocated at least double aligned

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

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

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

577: /*MC
578:    PetscRealloc - Rellocates memory

580:    Synopsis:
581: #include <petscsys.h>
582:    PetscErrorCode PetscRealloc(size_t m,void **result)

584:    Not Collective

586:    Input Parameters:
587: +  m - number of bytes to allocate
588: -  result - previous memory

590:    Output Parameter:
591: .  result - new memory allocated

593:    Level: developer

595:    Notes:
596:    Memory is always allocated at least double aligned

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

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

603: /*MC
604:    PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment

606:    Synopsis:
607: #include <petscsys.h>
608:    void *PetscAddrAlign(void *addr)

610:    Not Collective

612:    Input Parameters:
613: .  addr - address to align (any pointer type)

615:    Level: developer

617: .seealso: PetscMallocAlign()

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

622: /*MC
623:    PetscMalloc1 - Allocates an array of memory aligned to PETSC_MEMALIGN

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

629:    Not Collective

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

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

637:    Note:
638:    This uses the sizeof() of the memory type requested to determine the total memory to be allocated, therefore you should not
639:          multiply the number of elements requested by the sizeof() the type. For example use
640: $  PetscInt *id;
641: $  PetscMalloc1(10,&id);
642:         not
643: $  PetscInt *id;
644: $  PetscMalloc1(10*sizeof(PetscInt),&id);

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

648:    Level: beginner

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

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

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

658:    Synopsis:
659: #include <petscsys.h>
660:    PetscErrorCode PetscCalloc1(size_t m1,type **r1)

662:    Not Collective

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

667:    Output Parameter:
668: .  r1 - memory allocated

670:    Notes:
671:    See PetsMalloc1() for more details on usage.

673:    Level: beginner

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

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

680: /*MC
681:    PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN

683:    Synopsis:
684: #include <petscsys.h>
685:    PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)

687:    Not Collective

689:    Input Parameter:
690: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
691: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

693:    Output Parameter:
694: +  r1 - memory allocated in first chunk
695: -  r2 - memory allocated in second chunk

697:    Level: developer

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

701: M*/
702: #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))

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

707:    Synopsis:
708: #include <petscsys.h>
709:    PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)

711:    Not Collective

713:    Input Parameter:
714: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
715: -  m2 - number of elements to allocate in 2nd chunk  (may be zero)

717:    Output Parameter:
718: +  r1 - memory allocated in first chunk
719: -  r2 - memory allocated in second chunk

721:    Level: developer

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

725: M*/
726: #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))

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

731:    Synopsis:
732: #include <petscsys.h>
733:    PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

735:    Not Collective

737:    Input Parameter:
738: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
739: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
740: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

742:    Output Parameter:
743: +  r1 - memory allocated in first chunk
744: .  r2 - memory allocated in second chunk
745: -  r3 - memory allocated in third chunk

747:    Level: developer

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

751: M*/
752: #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))

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

757:    Synopsis:
758: #include <petscsys.h>
759:    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)

761:    Not Collective

763:    Input Parameter:
764: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
765: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
766: -  m3 - number of elements to allocate in 3rd chunk  (may be zero)

768:    Output Parameter:
769: +  r1 - memory allocated in first chunk
770: .  r2 - memory allocated in second chunk
771: -  r3 - memory allocated in third chunk

773:    Level: developer

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

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

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

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

787:    Not Collective

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

795:    Output Parameter:
796: +  r1 - memory allocated in first chunk
797: .  r2 - memory allocated in second chunk
798: .  r3 - memory allocated in third chunk
799: -  r4 - memory allocated in fourth chunk

801:    Level: developer

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

805: M*/
806: #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))

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

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

815:    Not Collective

817:    Input Parameters:
818: +  m1 - number of elements to allocate in 1st chunk  (may be zero)
819: .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
820: .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
821: -  m4 - number of elements to allocate in 4th chunk  (may be zero)

823:    Output Parameters:
824: +  r1 - memory allocated in first chunk
825: .  r2 - memory allocated in second chunk
826: .  r3 - memory allocated in third chunk
827: -  r4 - memory allocated in fourth chunk

829:    Level: developer

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

833: M*/
834: #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))

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

839:    Synopsis:
840: #include <petscsys.h>
841:    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)

843:    Not Collective

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

852:    Output Parameters:
853: +  r1 - memory allocated in first chunk
854: .  r2 - memory allocated in second chunk
855: .  r3 - memory allocated in third chunk
856: .  r4 - memory allocated in fourth chunk
857: -  r5 - memory allocated in fifth chunk

859:    Level: developer

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

863: M*/
864: #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))

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

869:    Synopsis:
870: #include <petscsys.h>
871:    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)

873:    Not Collective

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

882:    Output Parameters:
883: +  r1 - memory allocated in first chunk
884: .  r2 - memory allocated in second chunk
885: .  r3 - memory allocated in third chunk
886: .  r4 - memory allocated in fourth chunk
887: -  r5 - memory allocated in fifth chunk

889:    Level: developer

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

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

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

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

903:    Not Collective

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

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

921:    Level: developer

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

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

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

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

935:    Not Collective

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

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

953:    Level: developer

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

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

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

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

967:    Not Collective

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

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

987:    Level: developer

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

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

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

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

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)
1010: -  m7 - number of elements to allocate in 7th chunk  (may be zero)

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

1021:    Level: developer

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

1025: M*/
1026: #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))

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

1031:    Synopsis:
1032: #include <petscsys.h>
1033:    PetscErrorCode PetscNew(type **result)

1035:    Not Collective

1037:    Output Parameter:
1038: .  result - memory allocated, sized to match pointer type

1040:    Level: beginner

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

1044: M*/
1045: #define PetscNew(b)      PetscCalloc1(1,(b))

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

1051:    Synopsis:
1052: #include <petscsys.h>
1053:    PetscErrorCode PetscNewLog(PetscObject obj,type **result)

1055:    Not Collective

1057:    Input Parameter:
1058: .  obj - object memory is logged to

1060:    Output Parameter:
1061: .  result - memory allocated, sized to match pointer type

1063:    Level: developer

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

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

1070: /*MC
1071:    PetscFree - Frees memory

1073:    Synopsis:
1074: #include <petscsys.h>
1075:    PetscErrorCode PetscFree(void *memory)

1077:    Not Collective

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

1082:    Level: beginner

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

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

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

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

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

1097:    Synopsis:
1098: #include <petscsys.h>
1099:    PetscErrorCode PetscFree2(void *memory1,void *memory2)

1101:    Not Collective

1103:    Input Parameters:
1104: +   memory1 - memory to free
1105: -   memory2 - 2nd memory to free

1107:    Level: developer

1109:    Notes:
1110:     Memory must have been obtained with PetscMalloc2()

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

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

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

1120:    Synopsis:
1121: #include <petscsys.h>
1122:    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)

1124:    Not Collective

1126:    Input Parameters:
1127: +   memory1 - memory to free
1128: .   memory2 - 2nd memory to free
1129: -   memory3 - 3rd memory to free

1131:    Level: developer

1133:    Notes:
1134:     Memory must have been obtained with PetscMalloc3()

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

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

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

1144:    Synopsis:
1145: #include <petscsys.h>
1146:    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)

1148:    Not Collective

1150:    Input Parameters:
1151: +   m1 - memory to free
1152: .   m2 - 2nd memory to free
1153: .   m3 - 3rd memory to free
1154: -   m4 - 4th memory to free

1156:    Level: developer

1158:    Notes:
1159:     Memory must have been obtained with PetscMalloc4()

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

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

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

1169:    Synopsis:
1170: #include <petscsys.h>
1171:    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)

1173:    Not Collective

1175:    Input Parameters:
1176: +   m1 - memory to free
1177: .   m2 - 2nd memory to free
1178: .   m3 - 3rd memory to free
1179: .   m4 - 4th memory to free
1180: -   m5 - 5th memory to free

1182:    Level: developer

1184:    Notes:
1185:     Memory must have been obtained with PetscMalloc5()

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

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

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

1195:    Synopsis:
1196: #include <petscsys.h>
1197:    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)

1199:    Not Collective

1201:    Input Parameters:
1202: +   m1 - memory to free
1203: .   m2 - 2nd memory to free
1204: .   m3 - 3rd memory to free
1205: .   m4 - 4th memory to free
1206: .   m5 - 5th memory to free
1207: -   m6 - 6th memory to free


1210:    Level: developer

1212:    Notes:
1213:     Memory must have been obtained with PetscMalloc6()

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

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

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

1223:    Synopsis:
1224: #include <petscsys.h>
1225:    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)

1227:    Not Collective

1229:    Input Parameters:
1230: +   m1 - memory to free
1231: .   m2 - 2nd memory to free
1232: .   m3 - 3rd memory to free
1233: .   m4 - 4th memory to free
1234: .   m5 - 5th memory to free
1235: .   m6 - 6th memory to free
1236: -   m7 - 7th memory to free


1239:    Level: developer

1241:    Notes:
1242:     Memory must have been obtained with PetscMalloc7()

1244: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1245:           PetscMalloc7()

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

1250: PETSC_EXTERN PetscErrorCode PetscMallocA(int,PetscBool,int,const char *,const char *,size_t,void *,...);
1251: PETSC_EXTERN PetscErrorCode PetscFreeA(int,int,const char *,const char *,void *,...);
1252: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,PetscBool,int,const char[],const char[],void**);
1253: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1254: PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t,int,const char[],const char[],void**);
1255: PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1256: 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 **));
1257: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);

1259: /*
1260:   Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1261: */
1262: PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1263: PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);
1264: #if defined(PETSC_HAVE_CUDA)
1265: PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void);
1266: PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void);
1267: #endif

1269: #define MPIU_PETSCLOGDOUBLE  MPI_DOUBLE
1270: #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION

1272: /*
1273:    Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1274: */
1275: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1276: PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1277: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1278: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1279: PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1280: PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int,PetscLogDouble*);
1281: PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool,PetscBool);
1282: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*,PetscBool*,PetscBool*);
1283: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1284: PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1285: PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool*);
1286: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool);
1287: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool*);

1289: PETSC_EXTERN const char *const PetscDataTypes[];
1290: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1291: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1292: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1293: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);

1295: /*
1296:     Basic memory and string operations. These are usually simple wrappers
1297:    around the basic Unix system calls, but a few of them have additional
1298:    functionality and/or error checking.
1299: */
1300: PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1301: PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1302: PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***);
1303: PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1304: PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool  *);
1305: PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool  *);
1306: PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1307: PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1308: PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1309: PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1310: PETSC_EXTERN PetscErrorCode PetscStrlcat(char[],const char[],size_t);
1311: PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1312: PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1313: PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1314: PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1315: PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1316: PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1317: PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1318: PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1319: PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*);
1320: PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1321: PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1322: PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1323: PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1324: PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt,const char *const*,char***);
1325: PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt,char***);
1326: PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);

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

1330: PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1331: PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1332: PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);

1334: PETSC_EXTERN PetscErrorCode PetscStrInList(const char[],const char[],char,PetscBool*);
1335: PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1336: PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);

1338: /*
1339:    These are MPI operations for MPI_Allreduce() etc
1340: */
1341: PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1342: #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1343: PETSC_EXTERN MPI_Op MPIU_SUM;
1344: #else
1345: #define MPIU_SUM MPI_SUM
1346: #endif
1347: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1348: PETSC_EXTERN MPI_Op MPIU_MAX;
1349: PETSC_EXTERN MPI_Op MPIU_MIN;
1350: #else
1351: #define MPIU_MAX MPI_MAX
1352: #define MPIU_MIN MPI_MIN
1353: #endif
1354: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);

1356: PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1357: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);

1359: PETSC_EXTERN const char *const PetscFileModes[];

1361: /*
1362:     Defines PETSc error handling.
1363: */
1364: #include <petscerror.h>

1366: #define PETSC_SMALLEST_CLASSID  1211211
1367: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1368: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1369: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);
1370: PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject,PetscObjectId*);
1371: PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject,PetscObjectId,PetscBool*);


1374: /*
1375:    Routines that get memory usage information from the OS
1376: */
1377: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1378: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1379: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1380: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);

1382: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);

1384: /*
1385:    Initialization of PETSc
1386: */
1387: PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1388: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1389: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1390: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1391: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1392: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1393: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1394: PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1395: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1396: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);

1398: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1399: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);

1401: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1402: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1403: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1404: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);

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


1409: /*
1410:      These are so that in extern C code we can caste function pointers to non-extern C
1411:    function pointers. Since the regular C++ code expects its function pointers to be C++
1412: */
1413: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1414: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1415: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);

1417: /*
1418:     Functions that can act on any PETSc object.
1419: */
1420: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1421: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1422: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1423: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1424: PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1425: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1426: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1427: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1428: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1429: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1430: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1431: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1432: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1433: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1434: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt*);
1435: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1436: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1437: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject*);
1438: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void));
1439: #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d))
1440: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1441: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1442: PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1443: PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject,PetscObject);
1444: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt*);

1446: #include <petscviewertypes.h>
1447: #include <petscoptions.h>

1449: PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer,PetscBool,PetscLogDouble);
1450: PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool*);

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

1454: PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer,const char[]);
1455: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1456: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1457: #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1458: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1459: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1460: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1461: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1462: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1463: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1464: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1465: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1466: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,PetscObject,const char[]);
1467: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1468: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1469: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject,const char[],PetscBool *);
1470: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1471: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1472: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1473: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);

1475: #if defined(PETSC_HAVE_SAWS)
1476: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1477: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1478: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1479: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1480: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1481: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1482: PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1483: PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1484: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1485: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);

1487: #else
1488: #define PetscSAWsBlock()                        0
1489: #define PetscObjectSAWsViewOff(obj)             0
1490: #define PetscObjectSAWsSetBlock(obj,flg)        0
1491: #define PetscObjectSAWsBlock(obj)               0
1492: #define PetscObjectSAWsGrantAccess(obj)         0
1493: #define PetscObjectSAWsTakeAccess(obj)          0
1494: #define PetscStackViewSAWs()                    0
1495: #define PetscStackSAWsViewOff()                 0
1496: #define PetscStackSAWsTakeAccess()
1497: #define PetscStackSAWsGrantAccess()

1499: #endif

1501: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1502: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1503: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);

1505: #if defined(PETSC_USE_DEBUG)
1506: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1507: #endif
1508: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);

1510: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1511: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1512: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1513: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1514: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1515: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);

1517: /*
1518:     Dynamic library lists. Lists of names of routines in objects or in dynamic
1519:   link libraries that will be loaded as needed.
1520: */

1522: #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1523: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1524: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1525: #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1526: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1527: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[],const char[]);
1528: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1529: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1530: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);

1532: PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1533: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1534: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1535: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1536: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1537: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1538: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1539: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);

1541: /*
1542:      Useful utility routines
1543: */
1544: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1545: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1546: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm,PetscInt*,PetscInt*);
1547: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1548: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1549: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1550: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);
1551: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm,PetscInt[],PetscInt[]);
1552: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm,PetscReal[],PetscReal[]);

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

1557:     Notes:
1558:     This is useful in cases like
1559: $     int        *a;
1560: $     PetscBool  flag = PetscNot(a)
1561:      where !a would not return a PetscBool because we cannot provide a cast from int to PetscBool in C.

1563:     Level: beginner

1565:     .seealso : PetscBool, PETSC_TRUE, PETSC_FALSE
1566: M*/
1567: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)

1569: /*MC
1570:     PetscHelpPrintf - Prints help messages.

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

1576:     Collective on comm

1578:     Input Parameters:
1579: +  comm - the MPI communicator over which the help message is printed
1580: .  format - the usual printf() format string
1581: -  args - arguments to be printed

1583:    Level: developer

1585:    Fortran Note:
1586:      This routine is not supported in Fortran.

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

1591:       To use, write your own function, for example,
1592: $PetscErrorCode mypetschelpprintf(MPI_Comm comm,const char format[],....)
1593: ${
1594: $ return(0);
1595: $}
1596: then do the assigment
1597: $    PetscHelpPrintf = mypetschelpprintf;
1598:    You can do the assignment before PetscInitialize().

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

1602: .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1603: M*/
1604: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);

1606: /*
1607:      Defines PETSc profiling.
1608: */
1609: #include <petsclog.h>

1611: /*
1612:       Simple PETSc parallel IO for ASCII printing
1613: */
1614: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1615: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1616: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1617: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1618: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1619: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1620: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);
1621: PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[],size_t,const char*,PetscInt,const PetscReal[]);

1623: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1624: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1625: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);

1627: PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char*,size_t*);
1628: PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char*,char *);

1630: #if defined(PETSC_HAVE_POPEN)
1631: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1632: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*);
1633: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1634: #endif

1636: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1637: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1638: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*);
1639: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1640: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1641: PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1642: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);

1644: PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1645: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void**);
1646: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void*);
1647: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1648: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer*);
1649: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer,PetscErrorCode (*)(void*));
1650: PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void*);

1652: /*
1653:    For use in debuggers
1654: */
1655: PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1656: PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1657: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1658: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1659: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);

1661: #include <stddef.h>
1662: #include <string.h>             /* for memcpy, memset */
1663: #include <stdlib.h>

1665: #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1666: #include <xmmintrin.h>
1667: #endif

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

1674:    Not Collective

1676:    Input Parameters:
1677: +  b - pointer to initial memory space
1678: .  a - pointer to copy space
1679: -  n - length (in bytes) of space to copy

1681:    Level: intermediate

1683:    Note:
1684:    PetscArraymove() is preferred
1685:    This routine is analogous to memmove().

1687:    Developers Note: This is inlined for performance

1689: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscStrallocpy(),
1690:           PetscArraymove()
1691: @*/
1692: PETSC_STATIC_INLINE PetscErrorCode PetscMemmove(void *a,const void *b,size_t n)
1693: {
1695:   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to null pointer");
1696:   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1697: #if !defined(PETSC_HAVE_MEMMOVE)
1698:   if (a < b) {
1699:     if (a <= b - n) memcpy(a,b,n);
1700:     else {
1701:       memcpy(a,b,(int)(b - a));
1702:       PetscMemmove(b,b + (int)(b - a),n - (int)(b - a));
1703:     }
1704:   } else {
1705:     if (b <= a - n) memcpy(a,b,n);
1706:     else {
1707:       memcpy(b + n,b + (n - (int)(a - b)),(int)(a - b));
1708:       PetscMemmove(a,b,n - (int)(a - b));
1709:     }
1710:   }
1711: #else
1712:   memmove((char*)(a),(char*)(b),n);
1713: #endif
1714:   return(0);
1715: }

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

1722:    Not Collective

1724:    Input Parameters:
1725: +  b - pointer to initial memory space
1726: -  n - length (in bytes) of space to copy

1728:    Output Parameter:
1729: .  a - pointer to copy space

1731:    Level: intermediate

1733:    Compile Option:
1734:     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1735:                                   for memory copies on double precision values.
1736:     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1737:                                   for memory copies on double precision values.
1738:     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1739:                                   for memory copies on double precision values.

1741:    Note:
1742:    Prefer PetscArraycpy()
1743:    This routine is analogous to memcpy().
1744:    Not available from Fortran

1746:    Developer Note: this is inlined for fastest performance

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

1750: @*/
1751: PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n)
1752: {
1753: #if defined(PETSC_USE_DEBUG)
1754:   size_t al = (size_t) a,bl = (size_t) b;
1755:   size_t nl = (size_t) n;
1757:   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1758:   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1759: #else
1761: #endif
1762:   if (a != b && n > 0) {
1763: #if defined(PETSC_USE_DEBUG)
1764:     if ((al > bl && (al - bl) < nl) || (bl - al) < nl)  SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1765:               or make sure your copy regions and lengths are correct. \n\
1766:               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1767: #endif
1768: #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1769:    if (!(a % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1770:       size_t len = n/sizeof(PetscScalar);
1771: #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1772:       PetscBLASInt   one = 1,blen;
1774:       PetscBLASIntCast(len,&blen);
1775:       PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one));
1776: #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1777:       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1778: #else
1779:       size_t      i;
1780:       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1781:       for (i=0; i<len; i++) y[i] = x[i];
1782: #endif
1783:     } else {
1784:       memcpy((char*)(a),(char*)(b),n);
1785:     }
1786: #else
1787:     memcpy((char*)(a),(char*)(b),n);
1788: #endif
1789:   }
1790:   return(0);
1791: }

1793: /*@C
1794:    PetscMemzero - Zeros the specified memory.

1796:    Not Collective

1798:    Input Parameters:
1799: +  a - pointer to beginning memory location
1800: -  n - length (in bytes) of memory to initialize

1802:    Level: intermediate

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

1808:    Not available from Fortran
1809:    Prefer PetscArrayzero()

1811:    Developer Note: this is inlined for fastest performance

1813: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy()
1814: @*/
1815: PETSC_STATIC_INLINE PetscErrorCode PetscMemzero(void *a,size_t n)
1816: {
1817:   if (n > 0) {
1818: #if defined(PETSC_USE_DEBUG)
1819:     if (!a) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer with %zu bytes",n);
1820: #endif
1821: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1822:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1823:       size_t      i,len = n/sizeof(PetscScalar);
1824:       PetscScalar *x = (PetscScalar*)a;
1825:       for (i=0; i<len; i++) x[i] = 0.0;
1826:     } else {
1827: #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1828:     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1829:       PetscInt len = n/sizeof(PetscScalar);
1830:       fortranzero_(&len,(PetscScalar*)a);
1831:     } else {
1832: #endif
1833: #if defined(PETSC_PREFER_BZERO)
1834:       bzero((char *)a,n);
1835: #else
1836:       memset((char*)a,0,n);
1837: #endif
1838: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1839:     }
1840: #endif
1841:   }
1842:   return 0;
1843: }

1845: /*MC
1846:    PetscArraycmp - Compares two arrays in memory.

1848:    Synopsis:
1849: #include <petscsys.h>
1850:     PetscErrorCode PetscArraycmp(const anytype *str1,const anytype *str2,size_t cnt,PetscBool *e)

1852:    Not Collective

1854:    Input Parameters:
1855: +  str1 - First array
1856: .  str2 - Second array
1857: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1859:    Output Parameters:
1860: .   e - PETSC_TRUE if equal else PETSC_FALSE.

1862:    Level: intermediate

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

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

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

1877:    Synopsis:
1878: #include <petscsys.h>
1879:     PetscErrorCode PetscArraymove(anytype *str1,const anytype *str2,size_t cnt)

1881:    Not Collective

1883:    Input Parameters:
1884: +  str1 - First array
1885: .  str2 - Second array
1886: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1888:    Level: intermediate

1890:    Note:
1891:    This routine is a preferred replacement to PetscMemmove()
1892:    The arrays must be of the same type

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

1898: /*MC
1899:    PetscArraycpy - Copies from one array in memory to another

1901:    Synopsis:
1902: #include <petscsys.h>
1903:     PetscErrorCode PetscArraycpy(anytype *str1,const anytype *str2,size_t cnt)

1905:    Not Collective

1907:    Input Parameters:
1908: +  str1 - First array (destination)
1909: .  str2 - Second array (source)
1910: -  cnt  - Count of the array, not in bytes, but number of entries in the arrays

1912:    Level: intermediate

1914:    Note:
1915:    This routine is a preferred replacement to PetscMemcpy()
1916:    The arrays must be of the same type

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

1922: /*MC
1923:    PetscArrayzero - Zeros an array in memory.

1925:    Synopsis:
1926: #include <petscsys.h>
1927:     PetscErrorCode PetscArrayzero(anytype *str1,size_t cnt)

1929:    Not Collective

1931:    Input Parameters:
1932: +  str1 - array
1933: -  cnt  - Count of the array, not in bytes, but number of entries in the array

1935:    Level: intermediate

1937:    Note:
1938:    This routine is a preferred replacement to PetscMemzero()

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

1944: /*MC
1945:    PetscPrefetchBlock - Prefetches a block of memory

1947:    Synopsis:
1948: #include <petscsys.h>
1949:     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)

1951:    Not Collective

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

1959:    Level: developer

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

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

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

1974: M*/
1975: #define PetscPrefetchBlock(a,n,rw,t) do {                               \
1976:     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
1977:     for (; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
1978:   } while (0)

1980: /*
1981:       Determine if some of the kernel computation routines use
1982:    Fortran (rather than C) for the numerical calculations. On some machines
1983:    and compilers (like complex numbers) the Fortran version of the routines
1984:    is faster than the C/C++ versions. The flag --with-fortran-kernels
1985:    should be used with ./configure to turn these on.
1986: */
1987: #if defined(PETSC_USE_FORTRAN_KERNELS)

1989: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1990: #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1991: #endif

1993: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
1994: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
1995: #endif

1997: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1998: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
1999: #endif

2001: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
2002: #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
2003: #endif

2005: #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
2006: #define PETSC_USE_FORTRAN_KERNEL_NORM
2007: #endif

2009: #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
2010: #define PETSC_USE_FORTRAN_KERNEL_MAXPY
2011: #endif

2013: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
2014: #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
2015: #endif

2017: #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
2018: #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
2019: #endif

2021: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
2022: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
2023: #endif

2025: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
2026: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
2027: #endif

2029: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
2030: #define PETSC_USE_FORTRAN_KERNEL_MDOT
2031: #endif

2033: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
2034: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
2035: #endif

2037: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
2038: #define PETSC_USE_FORTRAN_KERNEL_AYPX
2039: #endif

2041: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
2042: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
2043: #endif

2045: #endif

2047: /*
2048:     Macros for indicating code that should be compiled with a C interface,
2049:    rather than a C++ interface. Any routines that are dynamically loaded
2050:    (such as the PCCreate_XXX() routines) must be wrapped so that the name
2051:    mangler does not change the functions symbol name. This just hides the
2052:    ugly extern "C" {} wrappers.
2053: */
2054: #if defined(__cplusplus)
2055: #  define EXTERN_C_BEGIN extern "C" {
2056: #  define EXTERN_C_END }
2057: #else
2058: #  define EXTERN_C_BEGIN
2059: #  define EXTERN_C_END
2060: #endif

2062: /* --------------------------------------------------------------------*/

2064: /*MC
2065:     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
2066:         communication

2068:    Level: beginner

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

2072: .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
2073: M*/

2075: #if defined(PETSC_HAVE_MPIIO)
2076: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2077: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2078: PETSC_EXTERN PetscErrorCode MPIU_File_write_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2079: PETSC_EXTERN PetscErrorCode MPIU_File_read_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2080: PETSC_EXTERN PetscErrorCode MPIU_File_write_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2081: PETSC_EXTERN PetscErrorCode MPIU_File_read_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2082: #endif

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

2086: /* Limit MPI to 32-bits */
2087: #define PETSC_MPI_INT_MAX  2147483647
2088: #define PETSC_MPI_INT_MIN -2147483647
2089: /* Limit BLAS to 32-bits */
2090: #define PETSC_BLAS_INT_MAX  2147483647
2091: #define PETSC_BLAS_INT_MIN -2147483647

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

2097:    Not Collective

2099:    Input Parameter:
2100: .     a - the PetscInt64 value

2102:    Output Parameter:
2103: .     b - the resulting PetscInt value

2105:    Level: advanced

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

2109:    Not available from Fortran

2111: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscBLASIntCast(), PetscIntMultError(), PetscIntSumError()
2112: @*/
2113: PETSC_STATIC_INLINE PetscErrorCode PetscIntCast(PetscInt64 a,PetscInt *b)
2114: {
2116: #if !defined(PETSC_USE_64BIT_INDICES)
2117:   if ((a) > PETSC_MAX_INT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Integer too long for PetscInt, you may need to ./configure using --with-64-bit-indices");
2118: #endif
2119:   *b =  (PetscInt)(a);
2120:   return(0);
2121: }

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

2127:    Not Collective

2129:    Input Parameter:
2130: .     a - the PetscInt value

2132:    Output Parameter:
2133: .     b - the resulting PetscBLASInt value

2135:    Level: advanced

2137:    Notes:
2138:       Not available from Fortran
2139:       Errors if the integer is negative since PETSc calls to BLAS/LAPACK never need to cast negative integer inputs

2141: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscIntCast()
2142: @*/
2143: PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
2144: {
2146: #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
2147:   *b = 0;
2148:   if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array is too long 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");
2149: #endif
2150:   if (a < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Passing negative integer to BLAS/LAPACK routine");
2151:   *b =  (PetscBLASInt)(a);
2152:   return(0);
2153: }

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

2159:    Not Collective

2161:    Input Parameter:
2162: .     a - the PetscInt value

2164:    Output Parameter:
2165: .     b - the resulting PetscMPIInt value

2167:    Level: advanced

2169:    Not available from Fortran

2171: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntCast()
2172: @*/
2173: PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
2174: {
2176: #if defined(PETSC_USE_64BIT_INDICES)
2177:   *b = 0;
2178:   if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for MPI, which is restricted to 32 bit integers");
2179: #endif
2180:   *b =  (PetscMPIInt)(a);
2181:   return(0);
2182: }

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

2186: /*@C

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

2190:    Not Collective

2192:    Input Parameter:
2193: +     a - the PetscReal value
2194: -     b - the second value

2196:    Returns:
2197:       the result as a PetscInt value

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

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

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

2208:    Not available from Fortran

2210:    Level: advanced

2212: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError(), PetscIntSumError()
2213: @*/
2214: PETSC_STATIC_INLINE PetscInt PetscRealIntMultTruncate(PetscReal a,PetscInt b)
2215: {
2216:   PetscInt64 r;

2218:   r  =  (PetscInt64) (a*(PetscReal)b);
2219:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2220:   return (PetscInt) r;
2221: }

2223: /*@C

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

2227:    Not Collective

2229:    Input Parameter:
2230: +     a - the PetscInt value
2231: -     b - the second value

2233:    Returns:
2234:       the result as a PetscInt value

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

2240:    Not available from Fortran

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

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

2246:    Level: advanced

2248: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError(), PetscIntSumError()
2249: @*/
2250: PETSC_STATIC_INLINE PetscInt PetscIntMultTruncate(PetscInt a,PetscInt b)
2251: {
2252:   PetscInt64 r;

2254:   r  =  PetscInt64Mult(a,b);
2255:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2256:   return (PetscInt) r;
2257: }

2259: /*@C

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

2263:    Not Collective

2265:    Input Parameter:
2266: +     a - the PetscInt value
2267: -     b - the second value

2269:    Returns:
2270:      the result as a PetscInt value

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

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

2278:    Not available from Fortran

2280:    Level: advanced

2282: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError()
2283: @*/
2284: PETSC_STATIC_INLINE PetscInt PetscIntSumTruncate(PetscInt a,PetscInt b)
2285: {
2286:   PetscInt64 r;

2288:   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2289:   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2290:   return (PetscInt) r;
2291: }

2293: /*@C

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

2297:    Not Collective

2299:    Input Parameter:
2300: +     a - the PetscInt value
2301: -     b - the second value

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

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

2309:    Not available from Fortran

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

2313:    Level: advanced

2315: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64(), PetscIntSumError()
2316: @*/
2317: PETSC_STATIC_INLINE PetscErrorCode PetscIntMultError(PetscInt a,PetscInt b,PetscInt *result)
2318: {
2319:   PetscInt64 r;

2322:   r  =  PetscInt64Mult(a,b);
2323: #if !defined(PETSC_USE_64BIT_INDICES)
2324:   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);
2325: #endif
2326:   if (result) *result = (PetscInt) r;
2327:   return(0);
2328: }

2330: /*@C

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

2334:    Not Collective

2336:    Input Parameter:
2337: +     a - the PetscInt value
2338: -     b - the second value

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

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

2346:    Not available from Fortran

2348:    Level: advanced

2350: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError()
2351: @*/
2352: PETSC_STATIC_INLINE PetscErrorCode PetscIntSumError(PetscInt a,PetscInt b,PetscInt *result)
2353: {
2354:   PetscInt64 r;

2357:   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2358: #if !defined(PETSC_USE_64BIT_INDICES)
2359:   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);
2360: #endif
2361:   if (result) *result = (PetscInt) r;
2362:   return(0);
2363: }

2365: /*
2366:      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2367: */
2368: #if defined(hz)
2369: #  undef hz
2370: #endif

2372: #include <limits.h>

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

2376: #define PETSC_BITS_PER_BYTE CHAR_BIT

2378: #if defined(PETSC_HAVE_SYS_TYPES_H)
2379: #  include <sys/types.h>
2380: #endif

2382: /*MC

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


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

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

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

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

2397:     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),
2398:     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
2399:     version number (x.y.z).

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

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

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

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

2409:     Level: intermediate

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

2413: M*/

2415: /*MC

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

2420: $
2421: $#include "petsc/finclude/petscXXX.h"
2422: $         use petscXXX

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

2426: $    XXX variablename
2427: $    type(tXXX) variablename

2429:     For example,

2431: $#include "petsc/finclude/petscvec.h"
2432: $         use petscvec
2433: $
2434: $    Vec b
2435: $    type(tVec) x

2437:     Level: beginner

2439: M*/

2441: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2442: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2443: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2444: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2445: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2446: PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2447: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2448: PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt*,PetscInt*,PetscInt*,PetscInt*);

2450: PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt,const PetscInt[],PetscBool*);
2451: PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt,const PetscMPIInt[],PetscBool*);
2452: PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt,const PetscReal[],PetscBool*);
2453: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2454: PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt,PetscInt[]);
2455: PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt*,PetscInt[]);
2456: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2458: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2459: PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt*);
2460: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2461: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2462: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2463: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt[],PetscInt[],PetscInt[]);
2464: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2465: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2466: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2467: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt,PetscMPIInt[],PetscInt[]);
2468: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2469: PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt,PetscInt[],void*,size_t,void*);
2470: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2471: PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt,PetscReal[],PetscInt[]);
2472: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2473: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2474: PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal,PetscInt,const PetscReal[],PetscReal,PetscInt*);
2475: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2476: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2477: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2478: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt[],const PetscInt[],PetscInt,const PetscInt[],const PetscInt[],PetscInt*,PetscInt**,PetscInt**);
2479: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt*,PetscInt**);
2480: PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt,const PetscMPIInt[],PetscInt,const PetscMPIInt[],PetscInt*,PetscMPIInt**);
2481: PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *);

2483: PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt,void*,size_t,int(*)(const void*,const void*,void*),void*);
2484: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt,PetscInt[]);
2485: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt,PetscMPIInt[]);
2486: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt,PetscReal[]);
2487: PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt,void*,size_t,void*,size_t,int(*)(const void*,const void*,void*),void*);
2488: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt,PetscInt[],PetscInt[]);
2489: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt,PetscMPIInt[],PetscMPIInt[]);
2490: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt,PetscReal[],PetscInt[]);

2492: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2493: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);

2495: /*J
2496:     PetscRandomType - String with the name of a PETSc randomizer

2498:    Level: beginner

2500:    Notes:
2501:    To use SPRNG or RANDOM123 you must have ./configure PETSc
2502:    with the option --download-sprng or --download-random123

2504: .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2505: J*/
2506: typedef const char* PetscRandomType;
2507: #define PETSCRAND       "rand"
2508: #define PETSCRAND48     "rand48"
2509: #define PETSCSPRNG      "sprng"
2510: #define PETSCRANDER48   "rander48"
2511: #define PETSCRANDOM123  "random123"

2513: /* Logging support */
2514: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;

2516: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);

2518: /* Dynamic creation and loading functions */
2519: PETSC_EXTERN PetscFunctionList PetscRandomList;

2521: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2522: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2523: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2524: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType*);
2525: PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom,PetscObject,const char[]);
2526: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);

2528: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2529: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2530: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2531: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2532: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2533: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2534: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2535: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2536: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);

2538: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2539: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2540: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2541: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2542: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2543: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2544: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);
2545: PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2546: PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2547: PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);

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

2551: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscInt*,PetscDataType);
2552: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscInt*,PetscDataType);
2553: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,const void*,PetscInt,PetscDataType);
2554: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,const void*,PetscInt,PetscDataType);
2555: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2556: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2557: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2558: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2559: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2560: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2561: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2562: #if defined(PETSC_USE_SOCKET_VIEWER)
2563: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);
2564: #endif

2566: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2567: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2568: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);

2570: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2571: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool);
2572: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2573: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2574: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2575: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2576: PETSC_EXTERN PetscErrorCode PetscWaitOnError(void);

2578: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2579: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2580: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2581: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2582: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2583: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt*,const void*,PetscMPIInt*,PetscMPIInt**,void*)
2584:   PetscAttrMPIPointerWithType(6,3);
2585: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2586:                                                     PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2587:                                                     PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2588:   PetscAttrMPIPointerWithType(6,3);
2589: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2590:                                                        MPI_Request**,MPI_Request**,
2591:                                                        PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2592:                                                        PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2593:   PetscAttrMPIPointerWithType(6,3);

2595: PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2596: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2597: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);

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

2601: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);

2603: PETSC_EXTERN const char *const PetscSubcommTypes[];

2605: struct _n_PetscSubcomm {
2606:   MPI_Comm         parent;           /* parent communicator */
2607:   MPI_Comm         dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2608:   MPI_Comm         child;            /* the sub-communicator */
2609:   PetscMPIInt      n;                /* num of subcommunicators under the parent communicator */
2610:   PetscMPIInt      color;            /* color of processors belong to this communicator */
2611:   PetscMPIInt      *subsize;         /* size of subcommunicator[color] */
2612:   PetscSubcommType type;
2613:   char             *subcommprefix;
2614: };

2616: PETSC_STATIC_INLINE MPI_Comm PetscSubcommParent(PetscSubcomm scomm) {return scomm->parent;}
2617: PETSC_STATIC_INLINE MPI_Comm PetscSubcommChild(PetscSubcomm scomm) {return scomm->child;}
2618: PETSC_STATIC_INLINE MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm) {return scomm->dupparent;}
2619: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2620: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2621: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2622: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2623: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2624: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2625: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2626: PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm,const char[]);
2627: PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm,MPI_Comm*);
2628: PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm,MPI_Comm*);
2629: PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm,MPI_Comm*);

2631: PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt,PetscHeap*);
2632: PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap,PetscInt,PetscInt);
2633: PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap,PetscInt*,PetscInt*);
2634: PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap,PetscInt*,PetscInt*);
2635: PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap,PetscInt,PetscInt);
2636: PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2637: PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap*);
2638: PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap,PetscViewer);

2640: PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2641: PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm,PetscShmComm*);
2642: PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2643: PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2644: PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm,MPI_Comm*);

2646: /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2647: PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm,PetscInt,PetscOmpCtrl*);
2648: PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl,MPI_Comm*,MPI_Comm*,PetscBool*);
2649: PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl*);
2650: PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2651: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2652: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);

2654: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2655: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2656: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2657: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2658: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2659: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2660: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2661: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);


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

2669: extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2670: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted*);
2671: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted*);
2672: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted,const char*,const char*,PetscBool*);

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

2678: PETSC_EXTERN PetscSegBuffer PetscCitationsList;

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

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

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

2689:    Level: intermediate

2691:    Not available from Fortran

2693:      Options Database:
2694: .     -citations [filename]   - print out the bibtex entries for the given computation
2695: @*/
2696: PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2697: {
2698:   size_t         len;
2699:   char           *vstring;

2703:   if (set && *set) return(0);
2704:   PetscStrlen(cit,&len);
2705:   PetscSegBufferGet(PetscCitationsList,len,&vstring);
2706:   PetscArraycpy(vstring,cit,len);
2707:   if (set) *set = PETSC_TRUE;
2708:   return(0);
2709: }

2711: PETSC_EXTERN PetscErrorCode PetscURLShorten(const char[],char[],size_t);
2712: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm,char[],char[],size_t);
2713: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm,const char[],char[],size_t);
2714: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm,const char[],const char []);

2716: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm,char[],char[],size_t);
2717: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm,const char[],char[],char[],size_t);

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

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

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


2728: #if defined(PETSC_USE_DEBUG)
2729: /*
2730:    Verify that all processes in the communicator have called this from the same line of code
2731:  */
2732: PETSC_EXTERN PetscErrorCode PetscAllreduceBarrierCheck(MPI_Comm,PetscMPIInt,int,const char*,const char *);

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

2739:    Collective

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

2745:    Input Parameters:
2746: +  indata - pointer to the input data to be reduced
2747: .  count - the number of MPI data items in indata and outdata
2748: .  datatype - the MPI datatype, for example MPI_INT
2749: .  op - the MPI operation, for example MPI_SUM
2750: -  comm - the MPI communicator on which the operation occurs

2752:    Output Parameter:
2753: .  outdata - the reduced values

2755:    Notes:
2756:    In optimized mode this directly calls MPI_Allreduce()

2758:    Level: developer

2760: .seealso: MPI_Allreduce()
2761: M*/
2762: #define MPIU_Allreduce(a,b,c,d,e,fcomm) (PetscAllreduceBarrierCheck(fcomm,c,__LINE__,PETSC_FUNCTION_NAME,__FILE__) || MPI_Allreduce(a,b,c,d,e,fcomm))
2763: #else
2764: #define MPIU_Allreduce(a,b,c,d,e,fcomm) MPI_Allreduce(a,b,c,d,e,fcomm)
2765: #endif

2767: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
2768: PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint,PetscMPIInt,MPI_Info,MPI_Comm,void*,MPI_Win*);
2769: PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win,PetscMPIInt,MPI_Aint*,PetscMPIInt*,void*);
2770: #endif

2772: /*
2773:     Returned from PETSc functions that are called from MPI, such as related to attributes
2774: */
2775: PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CLASS;
2776: PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CODE;

2778: #endif