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