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