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