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