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