xref: /petsc/include/petscsys.h (revision ecf241c801971573cb09fea9ebbae6d2eb887520)
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>
16 
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
20 
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
39 
40 #include <petscsystypes.h>
41 
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
51 
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
61 
62 #if defined(__cplusplus)
63 #  define PETSC_RESTRICT PETSC_CXX_RESTRICT
64 #else
65 #  define PETSC_RESTRICT PETSC_C_RESTRICT
66 #endif
67 
68 #if defined(__cplusplus)
69 #  define PETSC_INLINE PETSC_CXX_INLINE
70 #else
71 #  define PETSC_INLINE PETSC_C_INLINE
72 #endif
73 
74 #define PETSC_STATIC_INLINE static PETSC_INLINE
75 
76 #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES) /* For Win32 shared libraries */
77 #  define PETSC_DLLEXPORT __declspec(dllexport)
78 #  define PETSC_DLLIMPORT __declspec(dllimport)
79 #  define PETSC_VISIBILITY_INTERNAL
80 #elif defined(PETSC_USE_VISIBILITY_CXX) && defined(__cplusplus)
81 #  define PETSC_DLLEXPORT __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 PETSC_DLLEXPORT __attribute__((visibility ("default")))
86 #  define PETSC_DLLIMPORT __attribute__((visibility ("default")))
87 #  define PETSC_VISIBILITY_INTERNAL __attribute__((visibility ("hidden")))
88 #else
89 #  define PETSC_DLLEXPORT
90 #  define PETSC_DLLIMPORT
91 #  define PETSC_VISIBILITY_INTERNAL
92 #endif
93 
94 #if defined(petsc_EXPORTS)      /* CMake defines this when building the shared library */
95 #  define PETSC_VISIBILITY_PUBLIC PETSC_DLLEXPORT
96 #else  /* Win32 users need this to import symbols from petsc.dll */
97 #  define PETSC_VISIBILITY_PUBLIC PETSC_DLLIMPORT
98 #endif
99 
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
114 
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"
117 
118 /* ========================================================================== */
119 
120 /*
121     Defines the interface to MPI allowing the use of all MPI functions.
122 
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
140 
141 /*MC
142   PetscDefined - determine whether a boolean macro is defined
143 
144   Notes:
145   Typical usage is within normal code,
146 
147 $   if (PetscDefined(USE_DEBUG)) { ... }
148 
149   but can also be used in the preprocessor,
150 
151 $   #if PetscDefined(USE_DEBUG)
152 $     ...
153 $   #else
154 
155   Either way evaluates true if PETSC_USE_DEBUG is defined (merely defined or defined to 1) or undefined.  This macro
156   should not be used if its argument may be defined to a non-empty value other than 1.
157 
158   Developer Notes:
159   Getting something that works in C and CPP for an arg that may or may not be defined is tricky.  Here, if we have
160   "#define PETSC_HAVE_BOOGER 1" we match on the placeholder define, insert the "0," for arg1 and generate the triplet
161   (0, 1, 0).  Then the last step cherry picks the 2nd arg (a one).  When PETSC_HAVE_BOOGER is not defined, we generate
162   a (... 1, 0) pair, and when the last step cherry picks the 2nd arg, we get a zero.
163 
164   Our extra expansion via PetscDefined__take_second_expand() is needed with MSVC, which has a nonconforming
165   implementation of variadic macros.
166 
167   Level: developer
168 M*/
169 #if !defined(PETSC_SKIP_VARIADIC_MACROS)
170 #  define PetscDefined_arg_1    shift,
171 #  define PetscDefined_arg_     shift,
172 #  define PetscDefined__take_second_expanded(ignored, val, ...) val
173 #  define PetscDefined__take_second_expand(args) PetscDefined__take_second_expanded args
174 #  define PetscDefined__take_second(...) PetscDefined__take_second_expand((__VA_ARGS__))
175 #  define PetscDefined___(arg1_or_junk) PetscDefined__take_second(arg1_or_junk 1, 0, at_)
176 #  define PetscDefined__(value) PetscDefined___(PetscDefined_arg_ ## value)
177 #  define PetscDefined_(d)      PetscDefined__(d)
178 #  define PetscDefined(d)       PetscDefined_(PETSC_ ## d)
179 #endif
180 
181 /*
182    Perform various sanity checks that the correct mpi.h is being included at compile time.
183    This usually happens because
184       * either an unexpected mpi.h is in the default compiler path (i.e. in /usr/include) or
185       * an extra include path -I/something (which contains the unexpected mpi.h) is being passed to the compiler
186 */
187 #if defined(PETSC_HAVE_MPIUNI)
188 #  if !defined(MPIUNI_H)
189 #    error "PETSc was configured with --with-mpi=0 but now appears to be compiling using a different mpi.h"
190 #  endif
191 #elif defined(PETSC_HAVE_I_MPI_NUMVERSION)
192 #  if !defined(I_MPI_NUMVERSION)
193 #    error "PETSc was configured with I_MPI but now appears to be compiling using a non-I_MPI mpi.h"
194 #  elif I_MPI_NUMVERSION != PETSC_HAVE_I_MPI_NUMVERSION
195 #    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"
196 #  endif
197 #elif defined(PETSC_HAVE_MVAPICH2_NUMVERSION)
198 #  if !defined(MVAPICH2_NUMVERSION)
199 #    error "PETSc was configured with MVAPICH2 but now appears to be compiling using a non-MVAPICH2 mpi.h"
200 #  elif MVAPICH2_NUMVERSION != PETSC_HAVE_MVAPICH2_NUMVERSION
201 #    error "PETSc was configured with one MVAPICH2 mpi.h version but now appears to be compiling using a different MVAPICH2 mpi.h version"
202 #  endif
203 #elif defined(PETSC_HAVE_MPICH_NUMVERSION)
204 #  if !defined(MPICH_NUMVERSION) || defined(MVAPICH2_NUMVERSION) || defined(I_MPI_NUMVERSION)
205 #    error "PETSc was configured with MPICH but now appears to be compiling using a non-MPICH mpi.h"
206 #  elif (MPICH_NUMVERSION/100000 != PETSC_HAVE_MPICH_NUMVERSION/100000) || (MPICH_NUMVERSION%100000/1000 < PETSC_HAVE_MPICH_NUMVERSION%100000/1000)
207 #    error "PETSc was configured with one MPICH mpi.h version but now appears to be compiling using a different MPICH mpi.h version"
208 #  endif
209 #elif defined(PETSC_HAVE_OMPI_MAJOR_VERSION)
210 #  if !defined(OMPI_MAJOR_VERSION)
211 #    error "PETSc was configured with OpenMPI but now appears to be compiling using a non-OpenMPI mpi.h"
212 #  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)
213 #    error "PETSc was configured with one OpenMPI mpi.h version but now appears to be compiling using a different OpenMPI mpi.h version"
214 #  endif
215 #elif defined(PETSC_HAVE_MSMPI_VERSION)
216 #  if !defined(MSMPI_VER)
217 #    error "PETSc was configured with MSMPI but now appears to be compiling using a non-MSMPI mpi.h"
218 #  elif (MSMPI_VER != PETSC_HAVE_MSMPI_VERSION)
219 #    error "PETSc was configured with one MSMPI mpi.h version but now appears to be compiling using a different MSMPI mpi.h version"
220 #  endif
221 #elif defined(OMPI_MAJOR_VERSION) || defined(MPICH_NUMVERSION) || defined(MSMPI_VER)
222 #  error "PETSc was configured with undetermined MPI - but now appears to be compiling using any of OpenMPI, MS-MPI or a MPICH variant"
223 #endif
224 
225 /*
226     Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler
227     see the top of mpicxx.h in the MPICH2 distribution.
228 */
229 #include <stdio.h>
230 
231 /* MSMPI on 32bit windows requires this yukky hack - that breaks MPI standard compliance */
232 #if !defined(MPIAPI)
233 #define MPIAPI
234 #endif
235 
236 /*
237    Support for Clang (>=3.2) matching type tag arguments with void* buffer types.
238    This allows the compiler to detect cases where the MPI datatype argument passed to a MPI routine
239    does not match the actual type of the argument being passed in
240 */
241 #if defined(__has_attribute) && defined(works_with_const_which_is_not_true)
242 #  if __has_attribute(argument_with_type_tag) && __has_attribute(pointer_with_type_tag) && __has_attribute(type_tag_for_datatype)
243 #    define PetscAttrMPIPointerWithType(bufno,typeno) __attribute__((pointer_with_type_tag(MPI,bufno,typeno)))
244 #    define PetscAttrMPITypeTag(type)                 __attribute__((type_tag_for_datatype(MPI,type)))
245 #    define PetscAttrMPITypeTagLayoutCompatible(type) __attribute__((type_tag_for_datatype(MPI,type,layout_compatible)))
246 #  endif
247 #endif
248 #if !defined(PetscAttrMPIPointerWithType)
249 #  define PetscAttrMPIPointerWithType(bufno,typeno)
250 #  define PetscAttrMPITypeTag(type)
251 #  define PetscAttrMPITypeTagLayoutCompatible(type)
252 #endif
253 
254 PETSC_EXTERN MPI_Datatype MPIU_ENUM PetscAttrMPITypeTag(PetscEnum);
255 PETSC_EXTERN MPI_Datatype MPIU_BOOL PetscAttrMPITypeTag(PetscBool);
256 
257 /*MC
258    MPIU_INT - MPI datatype corresponding to PetscInt
259 
260    Notes:
261    In MPI calls that require an MPI datatype that matches a PetscInt or array of PetscInt values, pass this value.
262 
263    Level: beginner
264 
265 .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_SCALAR, MPIU_COMPLEX
266 M*/
267 
268 #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 */
269 #  define MPIU_INT64 MPI_INT64_T
270 #  define PetscInt64_FMT PRId64
271 #elif (PETSC_SIZEOF_LONG_LONG == 8)
272 #  define MPIU_INT64 MPI_LONG_LONG_INT
273 #  define PetscInt64_FMT "lld"
274 #elif defined(PETSC_HAVE___INT64)
275 #  define MPIU_INT64 MPI_INT64_T
276 #  define PetscInt64_FMT "ld"
277 #else
278 #  error "cannot determine PetscInt64 type"
279 #endif
280 
281 PETSC_EXTERN MPI_Datatype MPIU_FORTRANADDR;
282 
283 #if defined(PETSC_USE_64BIT_INDICES)
284 #  define MPIU_INT MPIU_INT64
285 #  define PetscInt_FMT PetscInt64_FMT
286 #else
287 #  define MPIU_INT MPI_INT
288 #  define PetscInt_FMT "d"
289 #endif
290 
291 /*
292     For the rare cases when one needs to send a size_t object with MPI
293 */
294 PETSC_EXTERN MPI_Datatype MPIU_SIZE_T;
295 
296 /*
297       You can use PETSC_STDOUT as a replacement of stdout. You can also change
298     the value of PETSC_STDOUT to redirect all standard output elsewhere
299 */
300 PETSC_EXTERN FILE* PETSC_STDOUT;
301 
302 /*
303       You can use PETSC_STDERR as a replacement of stderr. You can also change
304     the value of PETSC_STDERR to redirect all standard error elsewhere
305 */
306 PETSC_EXTERN FILE* PETSC_STDERR;
307 
308 /*MC
309     PetscUnlikely - hints the compiler that the given condition is usually FALSE
310 
311     Synopsis:
312     #include <petscsys.h>
313     PetscBool  PetscUnlikely(PetscBool  cond)
314 
315     Not Collective
316 
317     Input Parameters:
318 .   cond - condition or expression
319 
320     Notes:
321     This returns the same truth value, it is only a hint to compilers that the resulting
322     branch is unlikely.
323 
324     Level: advanced
325 
326 .seealso: PetscUnlikelyDebug(), PetscLikely(), CHKERRQ
327 M*/
328 
329 /*MC
330     PetscLikely - hints the compiler that the given condition is usually TRUE
331 
332     Synopsis:
333     #include <petscsys.h>
334     PetscBool  PetscLikely(PetscBool  cond)
335 
336     Not Collective
337 
338     Input Parameters:
339 .   cond - condition or expression
340 
341     Notes:
342     This returns the same truth value, it is only a hint to compilers that the resulting
343     branch is likely.
344 
345     Level: advanced
346 
347 .seealso: PetscUnlikely()
348 M*/
349 #if defined(PETSC_HAVE_BUILTIN_EXPECT)
350 #  define PetscUnlikely(cond)   __builtin_expect(!!(cond),0)
351 #  define PetscLikely(cond)     __builtin_expect(!!(cond),1)
352 #else
353 #  define PetscUnlikely(cond)   (cond)
354 #  define PetscLikely(cond)     (cond)
355 #endif
356 
357 /*MC
358     PetscUnlikelyDebug - hints the compiler that the given condition is usually FALSE, eliding the check in optimized mode
359 
360     Synopsis:
361     #include <petscsys.h>
362     PetscBool  PetscUnlikelyDebug(PetscBool  cond)
363 
364     Not Collective
365 
366     Input Parameters:
367 .   cond - condition or expression
368 
369     Notes:
370     This returns the same truth value, it is only a hint to compilers that the resulting
371     branch is unlikely.  When compiled in optimized mode, it always returns false.
372 
373     Level: advanced
374 
375 .seealso: PetscUnlikely(), CHKERRQ, SETERRQ
376 M*/
377 #define PetscUnlikelyDebug(cond) (PetscDefined(USE_DEBUG) && PetscUnlikely(cond))
378 
379 /* PetscPragmaSIMD - from CeedPragmaSIMD */
380 
381 #if defined(__INTEL_COMPILER) && !defined(_WIN32)
382 #  define PetscPragmaSIMD _Pragma("vector")
383 #elif defined(__GNUC__) && __GNUC__ >= 5
384 #  define PetscPragmaSIMD _Pragma("GCC ivdep")
385 #elif defined(_OPENMP) && _OPENMP >= 201307
386 #  define PetscPragmaSIMD _Pragma("omp simd")
387 #elif defined(PETSC_HAVE_CRAY_VECTOR)
388 #  define PetscPragmaSIMD _Pragma("_CRI ivdep")
389 #else
390 #  define PetscPragmaSIMD
391 #endif
392 
393 /*
394     Declare extern C stuff after including external header files
395 */
396 
397 PETSC_EXTERN const char *const PetscBools[];
398 
399 /*
400     Defines elementary mathematics functions and constants.
401 */
402 #include <petscmath.h>
403 
404 PETSC_EXTERN const char *const PetscCopyModes[];
405 
406 /*MC
407     PETSC_IGNORE - same as NULL, means PETSc will ignore this argument
408 
409    Level: beginner
410 
411    Note:
412    Accepted by many PETSc functions to not set a parameter and instead use
413           some default
414 
415    Fortran Notes:
416    This macro does not exist in Fortran; you must use PETSC_NULL_INTEGER,
417           PETSC_NULL_DOUBLE_PRECISION etc
418 
419 .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_DETERMINE
420 
421 M*/
422 #define PETSC_IGNORE NULL
423 
424 /* This is deprecated */
425 #define PETSC_NULL NULL
426 
427 /*MC
428     PETSC_DECIDE - standard way of passing in integer or floating point parameter
429        where you wish PETSc to use the default.
430 
431    Level: beginner
432 
433 .seealso: PETSC_DEFAULT, PETSC_IGNORE, PETSC_DETERMINE
434 
435 M*/
436 #define PETSC_DECIDE -1
437 
438 /*MC
439     PETSC_DETERMINE - standard way of passing in integer or floating point parameter
440        where you wish PETSc to compute the required value.
441 
442    Level: beginner
443 
444    Developer Note:
445    I would like to use const PetscInt PETSC_DETERMINE = PETSC_DECIDE; but for
446      some reason this is not allowed by the standard even though PETSC_DECIDE is a constant value.
447 
448 .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_IGNORE, VecSetSizes()
449 
450 M*/
451 #define PETSC_DETERMINE PETSC_DECIDE
452 
453 /*MC
454     PETSC_DEFAULT - standard way of passing in integer or floating point parameter
455        where you wish PETSc to use the default.
456 
457    Level: beginner
458 
459    Fortran Notes:
460    You need to use PETSC_DEFAULT_INTEGER or PETSC_DEFAULT_REAL.
461 
462 .seealso: PETSC_DECIDE, PETSC_IGNORE, PETSC_DETERMINE
463 
464 M*/
465 #define PETSC_DEFAULT -2
466 
467 /*MC
468     PETSC_COMM_WORLD - the equivalent of the MPI_COMM_WORLD communicator which represents
469            all the processes that PETSc knows about.
470 
471    Level: beginner
472 
473    Notes:
474    By default PETSC_COMM_WORLD and MPI_COMM_WORLD are identical unless you wish to
475           run PETSc on ONLY a subset of MPI_COMM_WORLD. In that case create your new (smaller)
476           communicator, call it, say comm, and set PETSC_COMM_WORLD = comm BEFORE calling
477           PetscInitialize(), but after MPI_Init() has been called.
478 
479           The value of PETSC_COMM_WORLD should never be USED/accessed before PetscInitialize()
480           is called because it may not have a valid value yet.
481 
482 .seealso: PETSC_COMM_SELF
483 
484 M*/
485 PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;
486 
487 /*MC
488     PETSC_COMM_SELF - This is always MPI_COMM_SELF
489 
490    Level: beginner
491 
492    Notes:
493    Do not USE/access or set this variable before PetscInitialize() has been called.
494 
495 .seealso: PETSC_COMM_WORLD
496 
497 M*/
498 #define PETSC_COMM_SELF MPI_COMM_SELF
499 
500 /*MC
501     PETSC_MPI_THREAD_REQUIRED - the required threading support used if PETSc initializes
502            MPI with MPI_Init_thread.
503 
504    Level: beginner
505 
506    Notes:
507    By default PETSC_MPI_THREAD_REQUIRED equals MPI_THREAD_FUNNELED.
508 
509 .seealso: PetscInitialize()
510 
511 M*/
512 PETSC_EXTERN PetscMPIInt PETSC_MPI_THREAD_REQUIRED;
513 
514 PETSC_EXTERN PetscBool PetscBeganMPI;
515 PETSC_EXTERN PetscBool PetscErrorHandlingInitialized;
516 PETSC_EXTERN PetscBool PetscInitializeCalled;
517 PETSC_EXTERN PetscBool PetscFinalizeCalled;
518 PETSC_EXTERN PetscBool PetscViennaCLSynchronize;
519 
520 PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm));
521 PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*);
522 PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm*);
523 
524 #if defined(PETSC_HAVE_CUDA)
525 PETSC_EXTERN PetscBool      PetscCUDASynchronize;
526 PETSC_EXTERN PetscErrorCode PetscCUDAInitialize(MPI_Comm,PetscInt);
527 PETSC_EXTERN PetscErrorCode PetscCUDAInitializeCheck(void);
528 #endif
529 
530 #if defined(PETSC_HAVE_HIP)
531 PETSC_EXTERN PetscBool      PetscHIPSynchronize;
532 PETSC_EXTERN PetscErrorCode PetscHIPInitialize(MPI_Comm,PetscInt);
533 PETSC_EXTERN PetscErrorCode PetscHIPInitializeCheck(void);
534 #endif
535 
536 #if defined(PETSC_HAVE_KOKKOS)
537 PETSC_EXTERN PetscErrorCode PetscKokkosInitializeCheck(void);  /* Initialize Kokkos if not yet. */
538 #endif
539 
540 #if defined(PETSC_HAVE_ELEMENTAL)
541 PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
542 PETSC_EXTERN PetscErrorCode PetscElementalInitialized(PetscBool*);
543 PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
544 #endif
545 
546 /*MC
547    PetscMalloc - Allocates memory, One should use PetscNew(), PetscMalloc1() or PetscCalloc1() usually instead of this
548 
549    Synopsis:
550     #include <petscsys.h>
551    PetscErrorCode PetscMalloc(size_t m,void **result)
552 
553    Not Collective
554 
555    Input Parameter:
556 .  m - number of bytes to allocate
557 
558    Output Parameter:
559 .  result - memory allocated
560 
561    Level: beginner
562 
563    Notes:
564    Memory is always allocated at least double aligned
565 
566    It is safe to allocate size 0 and pass the resulting pointer (which may or may not be NULL) to PetscFree().
567 
568 .seealso: PetscFree(), PetscNew()
569 
570 M*/
571 #define PetscMalloc(a,b)  ((*PetscTrMalloc)((a),PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(void**)(b)))
572 
573 /*MC
574    PetscRealloc - Rellocates memory
575 
576    Synopsis:
577     #include <petscsys.h>
578    PetscErrorCode PetscRealloc(size_t m,void **result)
579 
580    Not Collective
581 
582    Input Parameters:
583 +  m - number of bytes to allocate
584 -  result - previous memory
585 
586    Output Parameter:
587 .  result - new memory allocated
588 
589    Level: developer
590 
591    Notes:
592    Memory is always allocated at least double aligned
593 
594 .seealso: PetscMalloc(), PetscFree(), PetscNew()
595 
596 M*/
597 #define PetscRealloc(a,b)  ((*PetscTrRealloc)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__,(void**)(b)))
598 
599 /*MC
600    PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment
601 
602    Synopsis:
603     #include <petscsys.h>
604    void *PetscAddrAlign(void *addr)
605 
606    Not Collective
607 
608    Input Parameters:
609 .  addr - address to align (any pointer type)
610 
611    Level: developer
612 
613 .seealso: PetscMallocAlign()
614 
615 M*/
616 #define PetscAddrAlign(a) (void*)((((PETSC_UINTPTR_T)(a))+(PETSC_MEMALIGN-1)) & ~(PETSC_MEMALIGN-1))
617 
618 /*MC
619    PetscMalloc1 - Allocates an array of memory aligned to PETSC_MEMALIGN
620 
621    Synopsis:
622     #include <petscsys.h>
623    PetscErrorCode PetscMalloc1(size_t m1,type **r1)
624 
625    Not Collective
626 
627    Input Parameter:
628 .  m1 - number of elements to allocate  (may be zero)
629 
630    Output Parameter:
631 .  r1 - memory allocated
632 
633    Note:
634    This uses the sizeof() of the memory type requested to determine the total memory to be allocated, therefore you should not
635          multiply the number of elements requested by the sizeof() the type. For example use
636 $  PetscInt *id;
637 $  PetscMalloc1(10,&id);
638         not
639 $  PetscInt *id;
640 $  PetscMalloc1(10*sizeof(PetscInt),&id);
641 
642         Does not zero the memory allocated, use PetscCalloc1() to obtain memory that has been zeroed.
643 
644    Level: beginner
645 
646 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc1(), PetscMalloc2()
647 
648 M*/
649 #define PetscMalloc1(m1,r1) PetscMallocA(1,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1))
650 
651 /*MC
652    PetscCalloc1 - Allocates a cleared (zeroed) array of memory aligned to PETSC_MEMALIGN
653 
654    Synopsis:
655     #include <petscsys.h>
656    PetscErrorCode PetscCalloc1(size_t m1,type **r1)
657 
658    Not Collective
659 
660    Input Parameter:
661 .  m1 - number of elements to allocate in 1st chunk  (may be zero)
662 
663    Output Parameter:
664 .  r1 - memory allocated
665 
666    Notes:
667    See PetsMalloc1() for more details on usage.
668 
669    Level: beginner
670 
671 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc1(), PetscCalloc2()
672 
673 M*/
674 #define PetscCalloc1(m1,r1) PetscMallocA(1,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1))
675 
676 /*MC
677    PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN
678 
679    Synopsis:
680     #include <petscsys.h>
681    PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)
682 
683    Not Collective
684 
685    Input Parameter:
686 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
687 -  m2 - number of elements to allocate in 2nd chunk  (may be zero)
688 
689    Output Parameter:
690 +  r1 - memory allocated in first chunk
691 -  r2 - memory allocated in second chunk
692 
693    Level: developer
694 
695 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc1(), PetscCalloc2()
696 
697 M*/
698 #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))
699 
700 /*MC
701    PetscCalloc2 - Allocates 2 cleared (zeroed) arrays of memory both aligned to PETSC_MEMALIGN
702 
703    Synopsis:
704     #include <petscsys.h>
705    PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)
706 
707    Not Collective
708 
709    Input Parameter:
710 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
711 -  m2 - number of elements to allocate in 2nd chunk  (may be zero)
712 
713    Output Parameter:
714 +  r1 - memory allocated in first chunk
715 -  r2 - memory allocated in second chunk
716 
717    Level: developer
718 
719 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc1(), PetscMalloc2()
720 
721 M*/
722 #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))
723 
724 /*MC
725    PetscMalloc3 - Allocates 3 arrays of memory, all aligned to PETSC_MEMALIGN
726 
727    Synopsis:
728     #include <petscsys.h>
729    PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
730 
731    Not Collective
732 
733    Input Parameter:
734 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
735 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
736 -  m3 - number of elements to allocate in 3rd chunk  (may be zero)
737 
738    Output Parameter:
739 +  r1 - memory allocated in first chunk
740 .  r2 - memory allocated in second chunk
741 -  r3 - memory allocated in third chunk
742 
743    Level: developer
744 
745 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc3(), PetscFree3()
746 
747 M*/
748 #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))
749 
750 /*MC
751    PetscCalloc3 - Allocates 3 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN
752 
753    Synopsis:
754     #include <petscsys.h>
755    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
756 
757    Not Collective
758 
759    Input Parameter:
760 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
761 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
762 -  m3 - number of elements to allocate in 3rd chunk  (may be zero)
763 
764    Output Parameter:
765 +  r1 - memory allocated in first chunk
766 .  r2 - memory allocated in second chunk
767 -  r3 - memory allocated in third chunk
768 
769    Level: developer
770 
771 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc2(), PetscMalloc3(), PetscFree3()
772 
773 M*/
774 #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))
775 
776 /*MC
777    PetscMalloc4 - Allocates 4 arrays of memory, all aligned to PETSC_MEMALIGN
778 
779    Synopsis:
780     #include <petscsys.h>
781    PetscErrorCode PetscMalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
782 
783    Not Collective
784 
785    Input Parameter:
786 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
787 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
788 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
789 -  m4 - number of elements to allocate in 4th chunk  (may be zero)
790 
791    Output Parameter:
792 +  r1 - memory allocated in first chunk
793 .  r2 - memory allocated in second chunk
794 .  r3 - memory allocated in third chunk
795 -  r4 - memory allocated in fourth chunk
796 
797    Level: developer
798 
799 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc4(), PetscFree4()
800 
801 M*/
802 #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))
803 
804 /*MC
805    PetscCalloc4 - Allocates 4 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN
806 
807    Synopsis:
808     #include <petscsys.h>
809    PetscErrorCode PetscCalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
810 
811    Not Collective
812 
813    Input Parameters:
814 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
815 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
816 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
817 -  m4 - number of elements to allocate in 4th chunk  (may be zero)
818 
819    Output Parameters:
820 +  r1 - memory allocated in first chunk
821 .  r2 - memory allocated in second chunk
822 .  r3 - memory allocated in third chunk
823 -  r4 - memory allocated in fourth chunk
824 
825    Level: developer
826 
827 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc4(), PetscFree4()
828 
829 M*/
830 #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))
831 
832 /*MC
833    PetscMalloc5 - Allocates 5 arrays of memory, all aligned to PETSC_MEMALIGN
834 
835    Synopsis:
836     #include <petscsys.h>
837    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)
838 
839    Not Collective
840 
841    Input Parameters:
842 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
843 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
844 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
845 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
846 -  m5 - number of elements to allocate in 5th chunk  (may be zero)
847 
848    Output Parameters:
849 +  r1 - memory allocated in first chunk
850 .  r2 - memory allocated in second chunk
851 .  r3 - memory allocated in third chunk
852 .  r4 - memory allocated in fourth chunk
853 -  r5 - memory allocated in fifth chunk
854 
855    Level: developer
856 
857 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc5(), PetscFree5()
858 
859 M*/
860 #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))
861 
862 /*MC
863    PetscCalloc5 - Allocates 5 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN
864 
865    Synopsis:
866     #include <petscsys.h>
867    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)
868 
869    Not Collective
870 
871    Input Parameters:
872 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
873 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
874 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
875 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
876 -  m5 - number of elements to allocate in 5th chunk  (may be zero)
877 
878    Output Parameters:
879 +  r1 - memory allocated in first chunk
880 .  r2 - memory allocated in second chunk
881 .  r3 - memory allocated in third chunk
882 .  r4 - memory allocated in fourth chunk
883 -  r5 - memory allocated in fifth chunk
884 
885    Level: developer
886 
887 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc5(), PetscFree5()
888 
889 M*/
890 #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))
891 
892 /*MC
893    PetscMalloc6 - Allocates 6 arrays of memory, all aligned to PETSC_MEMALIGN
894 
895    Synopsis:
896     #include <petscsys.h>
897    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)
898 
899    Not Collective
900 
901    Input Parameters:
902 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
903 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
904 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
905 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
906 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
907 -  m6 - number of elements to allocate in 6th chunk  (may be zero)
908 
909    Output Parameteasr:
910 +  r1 - memory allocated in first chunk
911 .  r2 - memory allocated in second chunk
912 .  r3 - memory allocated in third chunk
913 .  r4 - memory allocated in fourth chunk
914 .  r5 - memory allocated in fifth chunk
915 -  r6 - memory allocated in sixth chunk
916 
917    Level: developer
918 
919 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc6(), PetscFree3(), PetscFree4(), PetscFree5(), PetscFree6()
920 
921 M*/
922 #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))
923 
924 /*MC
925    PetscCalloc6 - Allocates 6 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN
926 
927    Synopsis:
928     #include <petscsys.h>
929    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)
930 
931    Not Collective
932 
933    Input Parameters:
934 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
935 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
936 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
937 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
938 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
939 -  m6 - number of elements to allocate in 6th chunk  (may be zero)
940 
941    Output Parameters:
942 +  r1 - memory allocated in first chunk
943 .  r2 - memory allocated in second chunk
944 .  r3 - memory allocated in third chunk
945 .  r4 - memory allocated in fourth chunk
946 .  r5 - memory allocated in fifth chunk
947 -  r6 - memory allocated in sixth chunk
948 
949    Level: developer
950 
951 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscMalloc6(), PetscFree6()
952 
953 M*/
954 #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))
955 
956 /*MC
957    PetscMalloc7 - Allocates 7 arrays of memory, all aligned to PETSC_MEMALIGN
958 
959    Synopsis:
960     #include <petscsys.h>
961    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)
962 
963    Not Collective
964 
965    Input Parameters:
966 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
967 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
968 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
969 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
970 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
971 .  m6 - number of elements to allocate in 6th chunk  (may be zero)
972 -  m7 - number of elements to allocate in 7th chunk  (may be zero)
973 
974    Output Parameters:
975 +  r1 - memory allocated in first chunk
976 .  r2 - memory allocated in second chunk
977 .  r3 - memory allocated in third chunk
978 .  r4 - memory allocated in fourth chunk
979 .  r5 - memory allocated in fifth chunk
980 .  r6 - memory allocated in sixth chunk
981 -  r7 - memory allocated in seventh chunk
982 
983    Level: developer
984 
985 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc7(), PetscFree7()
986 
987 M*/
988 #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))
989 
990 /*MC
991    PetscCalloc7 - Allocates 7 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN
992 
993    Synopsis:
994     #include <petscsys.h>
995    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)
996 
997    Not Collective
998 
999    Input Parameters:
1000 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
1001 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
1002 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
1003 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
1004 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
1005 .  m6 - number of elements to allocate in 6th chunk  (may be zero)
1006 -  m7 - number of elements to allocate in 7th chunk  (may be zero)
1007 
1008    Output Parameters:
1009 +  r1 - memory allocated in first chunk
1010 .  r2 - memory allocated in second chunk
1011 .  r3 - memory allocated in third chunk
1012 .  r4 - memory allocated in fourth chunk
1013 .  r5 - memory allocated in fifth chunk
1014 .  r6 - memory allocated in sixth chunk
1015 -  r7 - memory allocated in seventh chunk
1016 
1017    Level: developer
1018 
1019 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscMalloc7(), PetscFree7()
1020 
1021 M*/
1022 #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))
1023 
1024 /*MC
1025    PetscNew - Allocates memory of a particular type, zeros the memory! Aligned to PETSC_MEMALIGN
1026 
1027    Synopsis:
1028     #include <petscsys.h>
1029    PetscErrorCode PetscNew(type **result)
1030 
1031    Not Collective
1032 
1033    Output Parameter:
1034 .  result - memory allocated, sized to match pointer type
1035 
1036    Level: beginner
1037 
1038 .seealso: PetscFree(), PetscMalloc(), PetscNewLog(), PetscCalloc1(), PetscMalloc1()
1039 
1040 M*/
1041 #define PetscNew(b)      PetscCalloc1(1,(b))
1042 
1043 /*MC
1044    PetscNewLog - Allocates memory of a type matching pointer, zeros the memory! Aligned to PETSC_MEMALIGN. Associates the memory allocated
1045          with the given object using PetscLogObjectMemory().
1046 
1047    Synopsis:
1048     #include <petscsys.h>
1049    PetscErrorCode PetscNewLog(PetscObject obj,type **result)
1050 
1051    Not Collective
1052 
1053    Input Parameter:
1054 .  obj - object memory is logged to
1055 
1056    Output Parameter:
1057 .  result - memory allocated, sized to match pointer type
1058 
1059    Level: developer
1060 
1061 .seealso: PetscFree(), PetscMalloc(), PetscNew(), PetscLogObjectMemory(), PetscCalloc1(), PetscMalloc1()
1062 
1063 M*/
1064 #define PetscNewLog(o,b) (PetscNew((b)) || PetscLogObjectMemory((PetscObject)o,sizeof(**(b))))
1065 
1066 /*MC
1067    PetscFree - Frees memory
1068 
1069    Synopsis:
1070     #include <petscsys.h>
1071    PetscErrorCode PetscFree(void *memory)
1072 
1073    Not Collective
1074 
1075    Input Parameter:
1076 .   memory - memory to free (the pointer is ALWAYS set to NULL upon success)
1077 
1078    Level: beginner
1079 
1080    Notes:
1081    Do not free memory obtained with PetscMalloc2(), PetscCalloc2() etc, they must be freed with PetscFree2() etc.
1082 
1083    It is safe to call PetscFree() on a NULL pointer.
1084 
1085 .seealso: PetscNew(), PetscMalloc(), PetscNewLog(), PetscMalloc1(), PetscCalloc1()
1086 
1087 M*/
1088 #define PetscFree(a)   ((*PetscTrFree)((void*)(a),__LINE__,PETSC_FUNCTION_NAME,__FILE__) || ((a) = NULL,0))
1089 
1090 /*MC
1091    PetscFree2 - Frees 2 chunks of memory obtained with PetscMalloc2()
1092 
1093    Synopsis:
1094     #include <petscsys.h>
1095    PetscErrorCode PetscFree2(void *memory1,void *memory2)
1096 
1097    Not Collective
1098 
1099    Input Parameters:
1100 +   memory1 - memory to free
1101 -   memory2 - 2nd memory to free
1102 
1103    Level: developer
1104 
1105    Notes:
1106     Memory must have been obtained with PetscMalloc2()
1107 
1108 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree()
1109 
1110 M*/
1111 #define PetscFree2(m1,m2)   PetscFreeA(2,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2))
1112 
1113 /*MC
1114    PetscFree3 - Frees 3 chunks of memory obtained with PetscMalloc3()
1115 
1116    Synopsis:
1117     #include <petscsys.h>
1118    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)
1119 
1120    Not Collective
1121 
1122    Input Parameters:
1123 +   memory1 - memory to free
1124 .   memory2 - 2nd memory to free
1125 -   memory3 - 3rd memory to free
1126 
1127    Level: developer
1128 
1129    Notes:
1130     Memory must have been obtained with PetscMalloc3()
1131 
1132 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3()
1133 
1134 M*/
1135 #define PetscFree3(m1,m2,m3)   PetscFreeA(3,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3))
1136 
1137 /*MC
1138    PetscFree4 - Frees 4 chunks of memory obtained with PetscMalloc4()
1139 
1140    Synopsis:
1141     #include <petscsys.h>
1142    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)
1143 
1144    Not Collective
1145 
1146    Input Parameters:
1147 +   m1 - memory to free
1148 .   m2 - 2nd memory to free
1149 .   m3 - 3rd memory to free
1150 -   m4 - 4th memory to free
1151 
1152    Level: developer
1153 
1154    Notes:
1155     Memory must have been obtained with PetscMalloc4()
1156 
1157 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4()
1158 
1159 M*/
1160 #define PetscFree4(m1,m2,m3,m4)   PetscFreeA(4,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4))
1161 
1162 /*MC
1163    PetscFree5 - Frees 5 chunks of memory obtained with PetscMalloc5()
1164 
1165    Synopsis:
1166     #include <petscsys.h>
1167    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)
1168 
1169    Not Collective
1170 
1171    Input Parameters:
1172 +   m1 - memory to free
1173 .   m2 - 2nd memory to free
1174 .   m3 - 3rd memory to free
1175 .   m4 - 4th memory to free
1176 -   m5 - 5th memory to free
1177 
1178    Level: developer
1179 
1180    Notes:
1181     Memory must have been obtained with PetscMalloc5()
1182 
1183 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5()
1184 
1185 M*/
1186 #define PetscFree5(m1,m2,m3,m4,m5)   PetscFreeA(5,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5))
1187 
1188 /*MC
1189    PetscFree6 - Frees 6 chunks of memory obtained with PetscMalloc6()
1190 
1191    Synopsis:
1192     #include <petscsys.h>
1193    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)
1194 
1195    Not Collective
1196 
1197    Input Parameters:
1198 +   m1 - memory to free
1199 .   m2 - 2nd memory to free
1200 .   m3 - 3rd memory to free
1201 .   m4 - 4th memory to free
1202 .   m5 - 5th memory to free
1203 -   m6 - 6th memory to free
1204 
1205 
1206    Level: developer
1207 
1208    Notes:
1209     Memory must have been obtained with PetscMalloc6()
1210 
1211 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6()
1212 
1213 M*/
1214 #define PetscFree6(m1,m2,m3,m4,m5,m6)   PetscFreeA(6,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5),&(m6))
1215 
1216 /*MC
1217    PetscFree7 - Frees 7 chunks of memory obtained with PetscMalloc7()
1218 
1219    Synopsis:
1220     #include <petscsys.h>
1221    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)
1222 
1223    Not Collective
1224 
1225    Input Parameters:
1226 +   m1 - memory to free
1227 .   m2 - 2nd memory to free
1228 .   m3 - 3rd memory to free
1229 .   m4 - 4th memory to free
1230 .   m5 - 5th memory to free
1231 .   m6 - 6th memory to free
1232 -   m7 - 7th memory to free
1233 
1234 
1235    Level: developer
1236 
1237    Notes:
1238     Memory must have been obtained with PetscMalloc7()
1239 
1240 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1241           PetscMalloc7()
1242 
1243 M*/
1244 #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   PetscFreeA(7,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5),&(m6),&(m7))
1245 
1246 PETSC_EXTERN PetscErrorCode PetscMallocA(int,PetscBool,int,const char *,const char *,size_t,void *,...);
1247 PETSC_EXTERN PetscErrorCode PetscFreeA(int,int,const char *,const char *,void *,...);
1248 PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,PetscBool,int,const char[],const char[],void**);
1249 PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1250 PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t,int,const char[],const char[],void**);
1251 PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1252 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 **));
1253 PETSC_EXTERN PetscErrorCode PetscMallocClear(void);
1254 
1255 /*
1256   Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1257 */
1258 PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1259 PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);
1260 #if defined(PETSC_HAVE_CUDA)
1261 PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void);
1262 PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void);
1263 #endif
1264 
1265 #define MPIU_PETSCLOGDOUBLE  MPI_DOUBLE
1266 #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION
1267 
1268 /*
1269    Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1270 */
1271 PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1272 PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1273 PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1274 PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1275 PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1276 PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int,PetscLogDouble*);
1277 PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool,PetscBool);
1278 PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*,PetscBool*,PetscBool*);
1279 PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1280 PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1281 PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool*);
1282 PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool);
1283 PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool*);
1284 
1285 PETSC_EXTERN const char *const PetscDataTypes[];
1286 PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1287 PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1288 PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1289 PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);
1290 
1291 /*
1292     Basic memory and string operations. These are usually simple wrappers
1293    around the basic Unix system calls, but a few of them have additional
1294    functionality and/or error checking.
1295 */
1296 PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1297 PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1298 PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***);
1299 PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1300 PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool  *);
1301 PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool  *);
1302 PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1303 PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1304 PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1305 PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1306 PETSC_EXTERN PetscErrorCode PetscStrlcat(char[],const char[],size_t);
1307 PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1308 PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1309 PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1310 PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1311 PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1312 PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1313 PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1314 PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1315 PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*);
1316 PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1317 PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1318 PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1319 PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1320 PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt,const char *const*,char***);
1321 PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt,char***);
1322 PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);
1323 
1324 PETSC_EXTERN void PetscStrcmpNoError(const char[],const char[],PetscBool  *);
1325 
1326 PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1327 PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1328 PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);
1329 
1330 PETSC_EXTERN PetscErrorCode PetscStrInList(const char[],const char[],char,PetscBool*);
1331 PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1332 PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);
1333 
1334 /*
1335    These are MPI operations for MPI_Allreduce() etc
1336 */
1337 PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1338 #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1339 PETSC_EXTERN MPI_Op MPIU_SUM;
1340 #else
1341 #define MPIU_SUM MPI_SUM
1342 #endif
1343 #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1344 PETSC_EXTERN MPI_Op MPIU_MAX;
1345 PETSC_EXTERN MPI_Op MPIU_MIN;
1346 #else
1347 #define MPIU_MAX MPI_MAX
1348 #define MPIU_MIN MPI_MIN
1349 #endif
1350 PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);
1351 
1352 PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1353 PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1354 
1355 PETSC_EXTERN const char *const PetscFileModes[];
1356 
1357 /*
1358     Defines PETSc error handling.
1359 */
1360 #include <petscerror.h>
1361 
1362 #define PETSC_SMALLEST_CLASSID  1211211
1363 PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1364 PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1365 PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);
1366 PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject,PetscObjectId*);
1367 PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject,PetscObjectId,PetscBool*);
1368 
1369 
1370 /*
1371    Routines that get memory usage information from the OS
1372 */
1373 PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1374 PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1375 PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1376 PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);
1377 
1378 PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);
1379 
1380 /*
1381    Initialization of PETSc
1382 */
1383 PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1384 PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1385 PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1386 PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1387 PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1388 PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1389 PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1390 PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1391 PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1392 PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);
1393 
1394 PETSC_EXTERN PetscErrorCode PetscEnd(void);
1395 PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);
1396 
1397 PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1398 PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1399 PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1400 PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);
1401 
1402 PETSC_EXTERN PetscErrorCode PetscMonitorCompare(PetscErrorCode (*)(void),void *,PetscErrorCode (*)(void**),PetscErrorCode (*)(void),void *,PetscErrorCode (*)(void**),PetscBool *);
1403 
1404 
1405 /*
1406      These are so that in extern C code we can caste function pointers to non-extern C
1407    function pointers. Since the regular C++ code expects its function pointers to be C++
1408 */
1409 PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1410 PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1411 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);
1412 
1413 /*
1414     Functions that can act on any PETSc object.
1415 */
1416 PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1417 PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1418 PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1419 PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1420 PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1421 PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1422 PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1423 PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1424 PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1425 PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1426 PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1427 PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1428 PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1429 PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1430 PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt*);
1431 PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1432 PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1433 PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject*);
1434 PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void));
1435 #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d))
1436 PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1437 PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1438 PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1439 PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject,PetscObject);
1440 PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt*);
1441 
1442 #include <petscviewertypes.h>
1443 #include <petscoptions.h>
1444 
1445 PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer,PetscBool,PetscLogDouble);
1446 PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool*);
1447 
1448 PETSC_EXTERN PetscErrorCode PetscObjectsListGetGlobalNumbering(MPI_Comm,PetscInt,PetscObject*,PetscInt*,PetscInt*);
1449 
1450 PETSC_EXTERN PetscErrorCode PetscMemoryShowUsage(PetscViewer,const char[]);
1451 PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer,const char[]);
1452 PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1453 PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1454 #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1455 PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1456 PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1457 PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1458 PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1459 PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1460 PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1461 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1462 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1463 PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,PetscObject,const char[]);
1464 PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1465 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1466 PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject,const char[],PetscBool *);
1467 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1468 PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1469 PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1470 PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);
1471 
1472 #if defined(PETSC_HAVE_SAWS)
1473 PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1474 PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1475 PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1476 PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1477 PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1478 PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1479 PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1480 PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1481 PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1482 PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);
1483 
1484 #else
1485 #define PetscSAWsBlock()                        0
1486 #define PetscObjectSAWsViewOff(obj)             0
1487 #define PetscObjectSAWsSetBlock(obj,flg)        0
1488 #define PetscObjectSAWsBlock(obj)               0
1489 #define PetscObjectSAWsGrantAccess(obj)         0
1490 #define PetscObjectSAWsTakeAccess(obj)          0
1491 #define PetscStackViewSAWs()                    0
1492 #define PetscStackSAWsViewOff()                 0
1493 #define PetscStackSAWsTakeAccess()
1494 #define PetscStackSAWsGrantAccess()
1495 
1496 #endif
1497 
1498 PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1499 PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1500 PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);
1501 
1502 #if defined(PETSC_USE_DEBUG)
1503 PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1504 #endif
1505 PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);
1506 
1507 PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1508 PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1509 PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1510 PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1511 PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1512 PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);
1513 
1514 /*
1515     Dynamic library lists. Lists of names of routines in objects or in dynamic
1516   link libraries that will be loaded as needed.
1517 */
1518 
1519 #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1520 PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1521 PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1522 #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1523 PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1524 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[],const char[]);
1525 PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1526 PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1527 PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);
1528 
1529 PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1530 PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1531 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1532 PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1533 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1534 PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1535 PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1536 PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);
1537 
1538 /*
1539      Useful utility routines
1540 */
1541 PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1542 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1543 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm,PetscInt*,PetscInt*);
1544 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1545 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1546 PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1547 PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);
1548 PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm,PetscInt[],PetscInt[]);
1549 PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm,PetscReal[],PetscReal[]);
1550 
1551 /*MC
1552     PetscNot - negates a logical type value and returns result as a PetscBool
1553 
1554     Notes:
1555     This is useful in cases like
1556 $     int        *a;
1557 $     PetscBool  flag = PetscNot(a)
1558      where !a would not return a PetscBool because we cannot provide a cast from int to PetscBool in C.
1559 
1560     Level: beginner
1561 
1562     .seealso : PetscBool, PETSC_TRUE, PETSC_FALSE
1563 M*/
1564 #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)
1565 
1566 /*MC
1567     PetscHelpPrintf - Prints help messages.
1568 
1569    Synopsis:
1570     #include <petscsys.h>
1571      PetscErrorCode (*PetscHelpPrintf)(const char format[],...);
1572 
1573     Not Collective
1574 
1575     Input Parameters:
1576 .   format - the usual printf() format string
1577 
1578    Level: developer
1579 
1580     Fortran Note:
1581     This routine is not supported in Fortran.
1582 
1583 
1584 .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1585 M*/
1586 PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);
1587 
1588 /*
1589      Defines PETSc profiling.
1590 */
1591 #include <petsclog.h>
1592 
1593 /*
1594       Simple PETSc parallel IO for ASCII printing
1595 */
1596 PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1597 PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1598 PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1599 PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1600 PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1601 PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1602 PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);
1603 PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[],size_t,const char*,PetscInt,const PetscReal[]);
1604 
1605 PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1606 PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1607 PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);
1608 
1609 PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char*,size_t*);
1610 PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char*,char *);
1611 
1612 #if defined(PETSC_HAVE_POPEN)
1613 PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1614 PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*);
1615 PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1616 #endif
1617 
1618 PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1619 PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1620 PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*);
1621 PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1622 PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1623 PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1624 PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);
1625 
1626 PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1627 PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void**);
1628 PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void*);
1629 PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1630 PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer*);
1631 PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer,PetscErrorCode (*)(void*));
1632 PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void*);
1633 
1634 /*
1635    For use in debuggers
1636 */
1637 PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1638 PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1639 PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1640 PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1641 PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);
1642 
1643 #include <stddef.h>
1644 #include <string.h>             /* for memcpy, memset */
1645 #include <stdlib.h>
1646 
1647 #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1648 #include <xmmintrin.h>
1649 #endif
1650 
1651 /*@C
1652    PetscMemmove - Copies n bytes, beginning at location b, to the space
1653    beginning at location a. Copying  between regions that overlap will
1654    take place correctly. Use PetscMemcpy() if the locations do not overlap
1655 
1656    Not Collective
1657 
1658    Input Parameters:
1659 +  b - pointer to initial memory space
1660 .  a - pointer to copy space
1661 -  n - length (in bytes) of space to copy
1662 
1663    Level: intermediate
1664 
1665    Note:
1666    PetscArraymove() is preferred
1667    This routine is analogous to memmove().
1668 
1669    Developers Note: This is inlined for performance
1670 
1671 .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscStrallocpy(),
1672           PetscArraymove()
1673 @*/
1674 PETSC_STATIC_INLINE PetscErrorCode PetscMemmove(void *a,const void *b,size_t n)
1675 {
1676   PetscFunctionBegin;
1677   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to null pointer");
1678   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1679 #if !defined(PETSC_HAVE_MEMMOVE)
1680   if (a < b) {
1681     if (a <= b - n) memcpy(a,b,n);
1682     else {
1683       memcpy(a,b,(int)(b - a));
1684       PetscMemmove(b,b + (int)(b - a),n - (int)(b - a));
1685     }
1686   } else {
1687     if (b <= a - n) memcpy(a,b,n);
1688     else {
1689       memcpy(b + n,b + (n - (int)(a - b)),(int)(a - b));
1690       PetscMemmove(a,b,n - (int)(a - b));
1691     }
1692   }
1693 #else
1694   memmove((char*)(a),(char*)(b),n);
1695 #endif
1696   PetscFunctionReturn(0);
1697 }
1698 
1699 /*@C
1700    PetscMemcpy - Copies n bytes, beginning at location b, to the space
1701    beginning at location a. The two memory regions CANNOT overlap, use
1702    PetscMemmove() in that case.
1703 
1704    Not Collective
1705 
1706    Input Parameters:
1707 +  b - pointer to initial memory space
1708 -  n - length (in bytes) of space to copy
1709 
1710    Output Parameter:
1711 .  a - pointer to copy space
1712 
1713    Level: intermediate
1714 
1715    Compile Option:
1716     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1717                                   for memory copies on double precision values.
1718     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1719                                   for memory copies on double precision values.
1720     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1721                                   for memory copies on double precision values.
1722 
1723    Note:
1724    Prefer PetscArraycpy()
1725    This routine is analogous to memcpy().
1726    Not available from Fortran
1727 
1728    Developer Note: this is inlined for fastest performance
1729 
1730 .seealso: PetscMemzero(), PetscMemcmp(), PetscArrayzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy()
1731 
1732 @*/
1733 PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n)
1734 {
1735 #if defined(PETSC_USE_DEBUG)
1736   size_t al = (size_t) a,bl = (size_t) b;
1737   size_t nl = (size_t) n;
1738   PetscFunctionBegin;
1739   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1740   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1741 #else
1742   PetscFunctionBegin;
1743 #endif
1744   if (a != b && n > 0) {
1745 #if defined(PETSC_USE_DEBUG)
1746     if ((al > bl && (al - bl) < nl) || (bl - al) < nl)  SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1747               or make sure your copy regions and lengths are correct. \n\
1748               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1749 #endif
1750 #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1751    if (!(a % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1752       size_t len = n/sizeof(PetscScalar);
1753 #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1754       PetscBLASInt   one = 1,blen;
1755       PetscErrorCode ierr;
1756       ierr = PetscBLASIntCast(len,&blen);CHKERRQ(ierr);
1757       PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one));
1758 #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1759       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1760 #else
1761       size_t      i;
1762       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1763       for (i=0; i<len; i++) y[i] = x[i];
1764 #endif
1765     } else {
1766       memcpy((char*)(a),(char*)(b),n);
1767     }
1768 #else
1769     memcpy((char*)(a),(char*)(b),n);
1770 #endif
1771   }
1772   PetscFunctionReturn(0);
1773 }
1774 
1775 /*@C
1776    PetscMemzero - Zeros the specified memory.
1777 
1778    Not Collective
1779 
1780    Input Parameters:
1781 +  a - pointer to beginning memory location
1782 -  n - length (in bytes) of memory to initialize
1783 
1784    Level: intermediate
1785 
1786    Compile Option:
1787    PETSC_PREFER_BZERO - on certain machines (the IBM RS6000) the bzero() routine happens
1788   to be faster than the memset() routine. This flag causes the bzero() routine to be used.
1789 
1790    Not available from Fortran
1791    Prefer PetscArrayzero()
1792 
1793    Developer Note: this is inlined for fastest performance
1794 
1795 .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy()
1796 @*/
1797 PETSC_STATIC_INLINE PetscErrorCode PetscMemzero(void *a,size_t n)
1798 {
1799   if (n > 0) {
1800 #if defined(PETSC_USE_DEBUG)
1801     if (!a) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer with %zu bytes",n);
1802 #endif
1803 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1804     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1805       size_t      i,len = n/sizeof(PetscScalar);
1806       PetscScalar *x = (PetscScalar*)a;
1807       for (i=0; i<len; i++) x[i] = 0.0;
1808     } else {
1809 #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1810     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1811       PetscInt len = n/sizeof(PetscScalar);
1812       fortranzero_(&len,(PetscScalar*)a);
1813     } else {
1814 #endif
1815 #if defined(PETSC_PREFER_BZERO)
1816       bzero((char *)a,n);
1817 #else
1818       memset((char*)a,0,n);
1819 #endif
1820 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1821     }
1822 #endif
1823   }
1824   return 0;
1825 }
1826 
1827 /*MC
1828    PetscArraycmp - Compares two arrays in memory.
1829 
1830    Synopsis:
1831     #include <petscsys.h>
1832     PetscErrorCode PetscArraycmp(const anytype *str1,const anytype *str2,size_t cnt,PetscBool *e)
1833 
1834    Not Collective
1835 
1836    Input Parameters:
1837 +  str1 - First array
1838 .  str2 - Second array
1839 -  cnt  - Count of the array, not in bytes, but number of entries in the arrays
1840 
1841    Output Parameters:
1842 .   e - PETSC_TRUE if equal else PETSC_FALSE.
1843 
1844    Level: intermediate
1845 
1846    Note:
1847    This routine is a preferred replacement to PetscMemcmp()
1848    The arrays must be of the same type
1849 
1850 .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy(),
1851           PetscArraymove()
1852 M*/
1853 #define  PetscArraycmp(str1,str2,cnt,e) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcmp(str1,str2,(size_t)(cnt)*sizeof(*(str1)),e))
1854 
1855 /*MC
1856    PetscArraymove - Copies from one array in memory to another, the arrays may overlap. Use PetscArraycpy() when the arrays
1857                     do not overlap
1858 
1859    Synopsis:
1860     #include <petscsys.h>
1861     PetscErrorCode PetscArraymove(anytype *str1,const anytype *str2,size_t cnt)
1862 
1863    Not Collective
1864 
1865    Input Parameters:
1866 +  str1 - First array
1867 .  str2 - Second array
1868 -  cnt  - Count of the array, not in bytes, but number of entries in the arrays
1869 
1870    Level: intermediate
1871 
1872    Note:
1873    This routine is a preferred replacement to PetscMemmove()
1874    The arrays must be of the same type
1875 
1876 .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycpy(), PetscMemmove(), PetscArraycmp(), PetscStrallocpy()
1877 M*/
1878 #define  PetscArraymove(str1,str2,cnt) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemmove(str1,str2,(size_t)(cnt)*sizeof(*(str1))))
1879 
1880 /*MC
1881    PetscArraycpy - Copies from one array in memory to another
1882 
1883    Synopsis:
1884     #include <petscsys.h>
1885     PetscErrorCode PetscArraycpy(anytype *str1,const anytype *str2,size_t cnt)
1886 
1887    Not Collective
1888 
1889    Input Parameters:
1890 +  str1 - First array (destination)
1891 .  str2 - Second array (source)
1892 -  cnt  - Count of the array, not in bytes, but number of entries in the arrays
1893 
1894    Level: intermediate
1895 
1896    Note:
1897    This routine is a preferred replacement to PetscMemcpy()
1898    The arrays must be of the same type
1899 
1900 .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraymove(), PetscMemmove(), PetscArraycmp(), PetscStrallocpy()
1901 M*/
1902 #define  PetscArraycpy(str1,str2,cnt) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcpy(str1,str2,(size_t)(cnt)*sizeof(*(str1))))
1903 
1904 /*MC
1905    PetscArrayzero - Zeros an array in memory.
1906 
1907    Synopsis:
1908     #include <petscsys.h>
1909     PetscErrorCode PetscArrayzero(anytype *str1,size_t cnt)
1910 
1911    Not Collective
1912 
1913    Input Parameters:
1914 +  str1 - array
1915 -  cnt  - Count of the array, not in bytes, but number of entries in the array
1916 
1917    Level: intermediate
1918 
1919    Note:
1920    This routine is a preferred replacement to PetscMemzero()
1921 
1922 .seealso: PetscMemcpy(), PetscMemcmp(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy(), PetscArraymove()
1923 M*/
1924 #define  PetscArrayzero(str1,cnt) PetscMemzero(str1,(size_t)(cnt)*sizeof(*(str1)))
1925 
1926 /*MC
1927    PetscPrefetchBlock - Prefetches a block of memory
1928 
1929    Synopsis:
1930     #include <petscsys.h>
1931     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)
1932 
1933    Not Collective
1934 
1935    Input Parameters:
1936 +  a - pointer to first element to fetch (any type but usually PetscInt or PetscScalar)
1937 .  n - number of elements to fetch
1938 .  rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors)
1939 -  t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note
1940 
1941    Level: developer
1942 
1943    Notes:
1944    The last two arguments (rw and t) must be compile-time constants.
1945 
1946    Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality.  Not all architectures offer
1947    equivalent locality hints, but the following macros are always defined to their closest analogue.
1948 +  PETSC_PREFETCH_HINT_NTA - Non-temporal.  Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched).
1949 .  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.
1950 .  PETSC_PREFETCH_HINT_T1 - Fetch to level 2 and higher (not L1).
1951 -  PETSC_PREFETCH_HINT_T2 - Fetch to high-level cache only.  (On many systems, T0 and T1 are equivalent.)
1952 
1953    This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid
1954    address).
1955 
1956 M*/
1957 #define PetscPrefetchBlock(a,n,rw,t) do {                               \
1958     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
1959     for (; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
1960   } while (0)
1961 
1962 /*
1963       Determine if some of the kernel computation routines use
1964    Fortran (rather than C) for the numerical calculations. On some machines
1965    and compilers (like complex numbers) the Fortran version of the routines
1966    is faster than the C/C++ versions. The flag --with-fortran-kernels
1967    should be used with ./configure to turn these on.
1968 */
1969 #if defined(PETSC_USE_FORTRAN_KERNELS)
1970 
1971 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1972 #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
1973 #endif
1974 
1975 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
1976 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
1977 #endif
1978 
1979 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1980 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
1981 #endif
1982 
1983 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1984 #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
1985 #endif
1986 
1987 #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
1988 #define PETSC_USE_FORTRAN_KERNEL_NORM
1989 #endif
1990 
1991 #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
1992 #define PETSC_USE_FORTRAN_KERNEL_MAXPY
1993 #endif
1994 
1995 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1996 #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
1997 #endif
1998 
1999 #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
2000 #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
2001 #endif
2002 
2003 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
2004 #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
2005 #endif
2006 
2007 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
2008 #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
2009 #endif
2010 
2011 #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
2012 #define PETSC_USE_FORTRAN_KERNEL_MDOT
2013 #endif
2014 
2015 #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
2016 #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
2017 #endif
2018 
2019 #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
2020 #define PETSC_USE_FORTRAN_KERNEL_AYPX
2021 #endif
2022 
2023 #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
2024 #define PETSC_USE_FORTRAN_KERNEL_WAXPY
2025 #endif
2026 
2027 #endif
2028 
2029 /*
2030     Macros for indicating code that should be compiled with a C interface,
2031    rather than a C++ interface. Any routines that are dynamically loaded
2032    (such as the PCCreate_XXX() routines) must be wrapped so that the name
2033    mangler does not change the functions symbol name. This just hides the
2034    ugly extern "C" {} wrappers.
2035 */
2036 #if defined(__cplusplus)
2037 #  define EXTERN_C_BEGIN extern "C" {
2038 #  define EXTERN_C_END }
2039 #else
2040 #  define EXTERN_C_BEGIN
2041 #  define EXTERN_C_END
2042 #endif
2043 
2044 /* --------------------------------------------------------------------*/
2045 
2046 /*MC
2047     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
2048         communication
2049 
2050    Level: beginner
2051 
2052    Note: This manual page is a place-holder because MPICH does not have a manual page for MPI_Comm
2053 
2054 .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
2055 M*/
2056 
2057 #if defined(PETSC_HAVE_MPIIO)
2058 PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2059 PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2060 PETSC_EXTERN PetscErrorCode MPIU_File_write_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2061 PETSC_EXTERN PetscErrorCode MPIU_File_read_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2062 PETSC_EXTERN PetscErrorCode MPIU_File_write_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2063 PETSC_EXTERN PetscErrorCode MPIU_File_read_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2064 #endif
2065 
2066 /* the following petsc_static_inline require petscerror.h */
2067 
2068 /* Limit MPI to 32-bits */
2069 #define PETSC_MPI_INT_MAX  2147483647
2070 #define PETSC_MPI_INT_MIN -2147483647
2071 /* Limit BLAS to 32-bits */
2072 #define PETSC_BLAS_INT_MAX  2147483647
2073 #define PETSC_BLAS_INT_MIN -2147483647
2074 
2075 /*@C
2076     PetscIntCast - casts a PetscInt64 (which is 64 bits in size) to a PetscInt (which may be 32 bits in size), generates an
2077          error if the PetscInt is not large enough to hold the number.
2078 
2079    Not Collective
2080 
2081    Input Parameter:
2082 .     a - the PetscInt64 value
2083 
2084    Output Parameter:
2085 .     b - the resulting PetscInt value
2086 
2087    Level: advanced
2088 
2089    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
2090 
2091    Not available from Fortran
2092 
2093 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscBLASIntCast()
2094 @*/
2095 PETSC_STATIC_INLINE PetscErrorCode PetscIntCast(PetscInt64 a,PetscInt *b)
2096 {
2097   PetscFunctionBegin;
2098 #if !defined(PETSC_USE_64BIT_INDICES)
2099   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");
2100 #endif
2101   *b =  (PetscInt)(a);
2102   PetscFunctionReturn(0);
2103 }
2104 
2105 /*@C
2106     PetscBLASIntCast - casts a PetscInt (which may be 64 bits in size) to a PetscBLASInt (which may be 32 bits in size), generates an
2107          error if the PetscBLASInt is not large enough to hold the number.
2108 
2109    Not Collective
2110 
2111    Input Parameter:
2112 .     a - the PetscInt value
2113 
2114    Output Parameter:
2115 .     b - the resulting PetscBLASInt value
2116 
2117    Level: advanced
2118 
2119    Notes:
2120       Not available from Fortran
2121       Errors if the integer is negative since PETSc calls to BLAS/LAPACK never need to cast negative integer inputs
2122 
2123 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscIntCast()
2124 @*/
2125 PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
2126 {
2127   PetscFunctionBegin;
2128 #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
2129   *b = 0;
2130   if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK");
2131 #endif
2132   if (a < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Passing negative integer to BLAS/LAPACK routine");
2133   *b =  (PetscBLASInt)(a);
2134   PetscFunctionReturn(0);
2135 }
2136 
2137 /*@C
2138     PetscMPIIntCast - casts a PetscInt (which may be 64 bits in size) to a PetscMPIInt (which may be 32 bits in size), generates an
2139          error if the PetscMPIInt is not large enough to hold the number.
2140 
2141    Not Collective
2142 
2143    Input Parameter:
2144 .     a - the PetscInt value
2145 
2146    Output Parameter:
2147 .     b - the resulting PetscMPIInt value
2148 
2149    Level: advanced
2150 
2151    Not available from Fortran
2152 
2153 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntCast()
2154 @*/
2155 PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
2156 {
2157   PetscFunctionBegin;
2158 #if defined(PETSC_USE_64BIT_INDICES)
2159   *b = 0;
2160   if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for MPI");
2161 #endif
2162   *b =  (PetscMPIInt)(a);
2163   PetscFunctionReturn(0);
2164 }
2165 
2166 #define PetscInt64Mult(a,b)   ((PetscInt64)(a))*((PetscInt64)(b))
2167 
2168 /*@C
2169 
2170    PetscRealIntMultTruncate - Computes the product of a positive PetscReal and a positive PetscInt and truncates the value to slightly less than the maximal possible value
2171 
2172    Not Collective
2173 
2174    Input Parameter:
2175 +     a - the PetscReal value
2176 -     b - the second value
2177 
2178    Returns:
2179       the result as a PetscInt value
2180 
2181    Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2182    Use PetscIntMultTruncate() to compute the product of two positive PetscInt and truncate to fit a PetscInt
2183    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
2184 
2185    Developers Note:
2186    We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks.
2187 
2188    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2189 
2190    Not available from Fortran
2191 
2192    Level: advanced
2193 
2194 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2195 @*/
2196 PETSC_STATIC_INLINE PetscInt PetscRealIntMultTruncate(PetscReal a,PetscInt b)
2197 {
2198   PetscInt64 r;
2199 
2200   r  =  (PetscInt64) (a*(PetscReal)b);
2201   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2202   return (PetscInt) r;
2203 }
2204 
2205 /*@C
2206 
2207    PetscIntMultTruncate - Computes the product of two positive PetscInt and truncates the value to slightly less than the maximal possible value
2208 
2209    Not Collective
2210 
2211    Input Parameter:
2212 +     a - the PetscInt value
2213 -     b - the second value
2214 
2215    Returns:
2216       the result as a PetscInt value
2217 
2218    Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2219    Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2220    Use PetscIntMultError() to compute the product of two PetscInt if you wish to generate an error if the result will not fit in a PetscInt
2221 
2222    Not available from Fortran
2223 
2224    Developers Note: We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks.
2225 
2226    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2227 
2228    Level: advanced
2229 
2230 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2231 @*/
2232 PETSC_STATIC_INLINE PetscInt PetscIntMultTruncate(PetscInt a,PetscInt b)
2233 {
2234   PetscInt64 r;
2235 
2236   r  =  PetscInt64Mult(a,b);
2237   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2238   return (PetscInt) r;
2239 }
2240 
2241 /*@C
2242 
2243    PetscIntSumTruncate - Computes the sum of two positive PetscInt and truncates the value to slightly less than the maximal possible value
2244 
2245    Not Collective
2246 
2247    Input Parameter:
2248 +     a - the PetscInt value
2249 -     b - the second value
2250 
2251    Returns:
2252      the result as a PetscInt value
2253 
2254    Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2255    Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2256    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
2257 
2258    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2259 
2260    Not available from Fortran
2261 
2262    Level: advanced
2263 
2264 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2265 @*/
2266 PETSC_STATIC_INLINE PetscInt PetscIntSumTruncate(PetscInt a,PetscInt b)
2267 {
2268   PetscInt64 r;
2269 
2270   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2271   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2272   return (PetscInt) r;
2273 }
2274 
2275 /*@C
2276 
2277    PetscIntMultError - Computes the product of two positive PetscInt and generates an error with overflow.
2278 
2279    Not Collective
2280 
2281    Input Parameter:
2282 +     a - the PetscInt value
2283 -     b - the second value
2284 
2285    Output Parameter:ma
2286 .     result - the result as a PetscInt value, or NULL if you do not want the result, you just want to check if it overflows
2287 
2288    Use PetscInt64Mult() to compute the product of two 32 bit PetscInt and store in a PetscInt64
2289    Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt
2290 
2291    Not available from Fortran
2292 
2293    Developers Note: We currently assume that PetscInt addition does not overflow, this is obviously wrong but requires many more checks.
2294 
2295    Level: advanced
2296 
2297 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64(), PetscIntSumError()
2298 @*/
2299 PETSC_STATIC_INLINE PetscErrorCode PetscIntMultError(PetscInt a,PetscInt b,PetscInt *result)
2300 {
2301   PetscInt64 r;
2302 
2303   PetscFunctionBegin;
2304   r  =  PetscInt64Mult(a,b);
2305 #if !defined(PETSC_USE_64BIT_INDICES)
2306   if (r > PETSC_MAX_INT) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Product of two integer %d %d overflow, you must ./configure PETSc with --with-64-bit-indices for the case you are running",a,b);
2307 #endif
2308   if (result) *result = (PetscInt) r;
2309   PetscFunctionReturn(0);
2310 }
2311 
2312 /*@C
2313 
2314    PetscIntSumError - Computes the sum of two positive PetscInt and generates an error with overflow.
2315 
2316    Not Collective
2317 
2318    Input Parameter:
2319 +     a - the PetscInt value
2320 -     b - the second value
2321 
2322    Output Parameter:ma
2323 .     c - the result as a PetscInt value,  or NULL if you do not want the result, you just want to check if it overflows
2324 
2325    Use PetscInt64Mult() to compute the product of two 32 bit PetscInt and store in a PetscInt64
2326    Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt
2327 
2328    Not available from Fortran
2329 
2330    Level: advanced
2331 
2332 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2333 @*/
2334 PETSC_STATIC_INLINE PetscErrorCode PetscIntSumError(PetscInt a,PetscInt b,PetscInt *result)
2335 {
2336   PetscInt64 r;
2337 
2338   PetscFunctionBegin;
2339   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2340 #if !defined(PETSC_USE_64BIT_INDICES)
2341   if (r > PETSC_MAX_INT) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Sum of two integer %d %d overflow, you must ./configure PETSc with --with-64-bit-indices for the case you are running",a,b);
2342 #endif
2343   if (result) *result = (PetscInt) r;
2344   PetscFunctionReturn(0);
2345 }
2346 
2347 /*
2348      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2349 */
2350 #if defined(hz)
2351 #  undef hz
2352 #endif
2353 
2354 #include <limits.h>
2355 
2356 /* The number of bits in a byte */
2357 
2358 #define PETSC_BITS_PER_BYTE CHAR_BIT
2359 
2360 #if defined(PETSC_HAVE_SYS_TYPES_H)
2361 #  include <sys/types.h>
2362 #endif
2363 
2364 /*MC
2365 
2366     PETSC_VERSION - This manual page provides information about how PETSc documents and uses its version information. This information is available to both C/C++
2367                     and Fortran compilers when petscsys.h is included.
2368 
2369 
2370     The current PETSc version and the API for accessing it are defined in petscversion.h
2371 
2372     The complete version number is given as the triple  PETSC_VERSION_MAJOR.PETSC_VERSION_MINOR.PETSC_VERSION_SUBMINOR (in short hand x.y.z)
2373 
2374     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
2375     where only a change in the major version number (x) indicates a change in the API.
2376 
2377     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
2378 
2379     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),
2380     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
2381     version number (x.y.z).
2382 
2383     PETSC_RELEASE_DATE is the date the x.y version was released (i.e. the version before any patch releases)
2384 
2385     PETSC_VERSION_DATE is the date the x.y.z version was released
2386 
2387     PETSC_VERSION_GIT is the last git commit to the repository given in the form vx.y.z-wwwww
2388 
2389     PETSC_VERSION_DATE_GIT is the date of the last git commit to the repository
2390 
2391     Level: intermediate
2392 
2393     PETSC_VERSION_() and PETSC_VERSION_PATCH are deprecated and will eventually be removed. For several releases PETSC_VERSION_PATCH is always 0
2394 
2395 M*/
2396 
2397 /*MC
2398 
2399     UsingFortran - To use PETSc with Fortran you must use both PETSc include files and modules. At the beginning
2400       of every function and module definition you need something like
2401 
2402 $
2403 $#include "petsc/finclude/petscXXX.h"
2404 $         use petscXXX
2405 
2406      You can declare PETSc variables using either of the following.
2407 
2408 $    XXX variablename
2409 $    type(tXXX) variablename
2410 
2411     For example,
2412 
2413 $#include "petsc/finclude/petscvec.h"
2414 $         use petscvec
2415 $
2416 $    Vec b
2417 $    type(tVec) x
2418 
2419     Level: beginner
2420 
2421 M*/
2422 
2423 PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2424 PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2425 PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2426 PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2427 PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2428 PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2429 PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2430 PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt*,PetscInt*,PetscInt*,PetscInt*);
2431 
2432 PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt,const PetscInt[],PetscBool*);
2433 PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt,const PetscMPIInt[],PetscBool*);
2434 PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt,const PetscReal[],PetscBool*);
2435 PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2436 PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt,PetscInt[]);
2437 PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt*,PetscInt[]);
2438 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2439 PETSC_EXTERN PetscErrorCode PetscCheckDupsInt(PetscInt,const PetscInt[],PetscBool*);
2440 PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2441 PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt*);
2442 PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2443 PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2444 PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2445 PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt[],PetscInt[],PetscInt[]);
2446 PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2447 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2448 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2449 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt,PetscMPIInt[],PetscInt[]);
2450 PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2451 PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt,PetscInt[],void*,size_t,void*);
2452 PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2453 PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt,PetscReal[],PetscInt[]);
2454 PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2455 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2456 PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal,PetscInt,const PetscReal[],PetscReal,PetscInt*);
2457 PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2458 PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2459 PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2460 PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt[],const PetscInt[],PetscInt,const PetscInt[],const PetscInt[],PetscInt*,PetscInt**,PetscInt**);
2461 PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt*,PetscInt**);
2462 PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt,const PetscMPIInt[],PetscInt,const PetscMPIInt[],PetscInt*,PetscMPIInt**);
2463 PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *);
2464 
2465 PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt,void*,size_t,int(*)(const void*,const void*,void*),void*);
2466 PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt,PetscInt[]);
2467 PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt,PetscMPIInt[]);
2468 PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt,PetscReal[]);
2469 PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt,void*,size_t,void*,size_t,int(*)(const void*,const void*,void*),void*);
2470 PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt,PetscInt[],PetscInt[]);
2471 PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt,PetscMPIInt[],PetscMPIInt[]);
2472 PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt,PetscReal[],PetscInt[]);
2473 
2474 PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2475 PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);
2476 
2477 /*J
2478     PetscRandomType - String with the name of a PETSc randomizer
2479 
2480    Level: beginner
2481 
2482    Notes:
2483    To use SPRNG or RANDOM123 you must have ./configure PETSc
2484    with the option --download-sprng or --download-random123
2485 
2486 .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2487 J*/
2488 typedef const char* PetscRandomType;
2489 #define PETSCRAND       "rand"
2490 #define PETSCRAND48     "rand48"
2491 #define PETSCSPRNG      "sprng"
2492 #define PETSCRANDER48   "rander48"
2493 #define PETSCRANDOM123  "random123"
2494 
2495 /* Logging support */
2496 PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;
2497 
2498 PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);
2499 
2500 /* Dynamic creation and loading functions */
2501 PETSC_EXTERN PetscFunctionList PetscRandomList;
2502 
2503 PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2504 PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2505 PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2506 PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType*);
2507 PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom,PetscObject,const char[]);
2508 PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);
2509 
2510 PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2511 PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2512 PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2513 PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2514 PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2515 PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2516 PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2517 PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2518 PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);
2519 
2520 PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2521 PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2522 PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2523 PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2524 PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2525 PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2526 PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);
2527 PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2528 PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2529 PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);
2530 
2531 PETSC_STATIC_INLINE PetscBool PetscBinaryBigEndian(void) {long _petsc_v = 1; return ((char*)&_petsc_v)[0] ? PETSC_FALSE : PETSC_TRUE;}
2532 
2533 PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscInt*,PetscDataType);
2534 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscInt*,PetscDataType);
2535 PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,const void*,PetscInt,PetscDataType);
2536 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,const void*,PetscInt,PetscDataType);
2537 PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2538 PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2539 PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2540 PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2541 PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2542 PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2543 PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2544 #if defined(PETSC_USE_SOCKET_VIEWER)
2545 PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);
2546 #endif
2547 
2548 PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2549 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2550 PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);
2551 
2552 PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2553 PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool);
2554 PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2555 PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2556 PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2557 PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2558 PETSC_EXTERN PetscErrorCode PetscWaitOnError(void);
2559 
2560 PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2561 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2562 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2563 PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2564 PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2565 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt*,const void*,PetscMPIInt*,PetscMPIInt**,void*)
2566   PetscAttrMPIPointerWithType(6,3);
2567 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2568                                                     PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2569                                                     PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2570   PetscAttrMPIPointerWithType(6,3);
2571 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2572                                                        MPI_Request**,MPI_Request**,
2573                                                        PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2574                                                        PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2575   PetscAttrMPIPointerWithType(6,3);
2576 
2577 PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2578 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2579 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);
2580 
2581 PETSC_EXTERN PetscErrorCode PetscSSEIsEnabled(MPI_Comm,PetscBool*,PetscBool*);
2582 
2583 PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);
2584 
2585 PETSC_EXTERN const char *const PetscSubcommTypes[];
2586 
2587 struct _n_PetscSubcomm {
2588   MPI_Comm         parent;           /* parent communicator */
2589   MPI_Comm         dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2590   MPI_Comm         child;            /* the sub-communicator */
2591   PetscMPIInt      n;                /* num of subcommunicators under the parent communicator */
2592   PetscMPIInt      color;            /* color of processors belong to this communicator */
2593   PetscMPIInt      *subsize;         /* size of subcommunicator[color] */
2594   PetscSubcommType type;
2595   char             *subcommprefix;
2596 };
2597 
2598 PETSC_STATIC_INLINE MPI_Comm PetscSubcommParent(PetscSubcomm scomm) {return scomm->parent;}
2599 PETSC_STATIC_INLINE MPI_Comm PetscSubcommChild(PetscSubcomm scomm) {return scomm->child;}
2600 PETSC_STATIC_INLINE MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm) {return scomm->dupparent;}
2601 PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2602 PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2603 PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2604 PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2605 PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2606 PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2607 PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2608 PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm,const char[]);
2609 PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm,MPI_Comm*);
2610 PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm,MPI_Comm*);
2611 PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm,MPI_Comm*);
2612 
2613 PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt,PetscHeap*);
2614 PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap,PetscInt,PetscInt);
2615 PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap,PetscInt*,PetscInt*);
2616 PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap,PetscInt*,PetscInt*);
2617 PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap,PetscInt,PetscInt);
2618 PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2619 PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap*);
2620 PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap,PetscViewer);
2621 
2622 PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2623 PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm,PetscShmComm*);
2624 PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2625 PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2626 PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm,MPI_Comm*);
2627 
2628 /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2629 PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm,PetscInt,PetscOmpCtrl*);
2630 PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl,MPI_Comm*,MPI_Comm*,PetscBool*);
2631 PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl*);
2632 PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2633 PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2634 PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);
2635 
2636 PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2637 PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2638 PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2639 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2640 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2641 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2642 PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2643 PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);
2644 
2645 
2646 /* Type-safe wrapper to encourage use of PETSC_RESTRICT. Does not use PetscFunctionBegin because the error handling
2647  * prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as
2648  * possible. */
2649 PETSC_STATIC_INLINE PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg,size_t count,PetscInt *PETSC_RESTRICT *slot) {return PetscSegBufferGet(seg,count,(void**)slot);}
2650 
2651 extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2652 PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted*);
2653 PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted*);
2654 PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted,const char*,const char*,PetscBool*);
2655 
2656 #include <stdarg.h>
2657 PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list);
2658 PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE*,const char[],va_list);
2659 
2660 PETSC_EXTERN PetscSegBuffer PetscCitationsList;
2661 
2662 /*@C
2663       PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code.
2664 
2665      Not Collective - only what is registered on rank 0 of PETSC_COMM_WORLD will be printed
2666 
2667      Input Parameters:
2668 +      cite - the bibtex item, formated to displayed on multiple lines nicely
2669 -      set - a boolean variable initially set to PETSC_FALSE; this is used to insure only a single registration of the citation
2670 
2671    Level: intermediate
2672 
2673    Not available from Fortran
2674 
2675      Options Database:
2676 .     -citations [filename]   - print out the bibtex entries for the given computation
2677 @*/
2678 PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2679 {
2680   size_t         len;
2681   char           *vstring;
2682   PetscErrorCode ierr;
2683 
2684   PetscFunctionBegin;
2685   if (set && *set) PetscFunctionReturn(0);
2686   ierr = PetscStrlen(cit,&len);CHKERRQ(ierr);
2687   ierr = PetscSegBufferGet(PetscCitationsList,len,&vstring);CHKERRQ(ierr);
2688   ierr = PetscArraycpy(vstring,cit,len);CHKERRQ(ierr);
2689   if (set) *set = PETSC_TRUE;
2690   PetscFunctionReturn(0);
2691 }
2692 
2693 PETSC_EXTERN PetscErrorCode PetscURLShorten(const char[],char[],size_t);
2694 PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm,char[],char[],size_t);
2695 PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm,const char[],char[],size_t);
2696 PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm,const char[],const char []);
2697 
2698 PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm,char[],char[],size_t);
2699 PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm,const char[],char[],char[],size_t);
2700 
2701 PETSC_EXTERN PetscErrorCode PetscGlobusGetTransfers(MPI_Comm,const char[],char[],size_t);
2702 
2703 PETSC_EXTERN PetscErrorCode PetscTextBelt(MPI_Comm,const char[],const char[],PetscBool*);
2704 PETSC_EXTERN PetscErrorCode PetscTellMyCell(MPI_Comm,const char[],const char[],PetscBool*);
2705 
2706 PETSC_EXTERN PetscErrorCode PetscPullJSONValue(const char[],const char[],char[],size_t,PetscBool*);
2707 PETSC_EXTERN PetscErrorCode PetscPushJSONValue(char[],const char[],const char[],size_t);
2708 
2709 
2710 #if defined(PETSC_USE_DEBUG)
2711 /*
2712    Verify that all processes in the communicator have called this from the same line of code
2713  */
2714 PETSC_EXTERN PetscErrorCode PetscAllreduceBarrierCheck(MPI_Comm,PetscMPIInt,int,const char*,const char *);
2715 
2716 /*MC
2717    MPIU_Allreduce - a PETSc replacement for MPI_Allreduce() that tries to determine if the call from all the MPI processes occur from the
2718                     same place in the PETSc code. This helps to detect bugs where different MPI processes follow different code paths
2719                     resulting in inconsistent and incorrect calls to MPI_Allreduce().
2720 
2721    Collective
2722 
2723    Synopsis:
2724      #include <petscsys.h>
2725      PetscErrorCode MPIU_Allreduce(void *indata,void *outdata,PetscMPIInt count,MPI_Datatype datatype, MPI_Op op, MPI_Comm comm);
2726 
2727    Input Parameters:
2728 +  indata - pointer to the input data to be reduced
2729 .  count - the number of MPI data items in indata and outdata
2730 .  datatype - the MPI datatype, for example MPI_INT
2731 .  op - the MPI operation, for example MPI_SUM
2732 -  comm - the MPI communicator on which the operation occurs
2733 
2734    Output Parameter:
2735 .  outdata - the reduced values
2736 
2737    Notes:
2738    In optimized mode this directly calls MPI_Allreduce()
2739 
2740    Level: developer
2741 
2742 .seealso: MPI_Allreduce()
2743 M*/
2744 #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))
2745 #else
2746 #define MPIU_Allreduce(a,b,c,d,e,fcomm) MPI_Allreduce(a,b,c,d,e,fcomm)
2747 #endif
2748 
2749 #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
2750 PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint,PetscMPIInt,MPI_Info,MPI_Comm,void*,MPI_Win*);
2751 PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win,PetscMPIInt,MPI_Aint*,PetscMPIInt*,void*);
2752 #endif
2753 
2754 /*
2755     Returned from PETSc functions that are called from MPI, such as related to attributes
2756 */
2757 PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CLASS;
2758 PETSC_EXTERN PetscMPIInt PETSC_MPI_ERROR_CODE;
2759 
2760 #endif
2761