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