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