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