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