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