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