xref: /petsc/include/petscsys.h (revision e854a5c95e8b401f0c8cdf7d9df27f3a3aee411b)
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 
1249    Level: developer
1250 
1251    Notes:
1252     Memory must have been obtained with PetscMalloc6()
1253 
1254 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6()
1255 
1256 M*/
1257 #define PetscFree6(m1,m2,m3,m4,m5,m6)   PetscFreeA(6,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5),&(m6))
1258 
1259 /*MC
1260    PetscFree7 - Frees 7 chunks of memory obtained with PetscMalloc7()
1261 
1262    Synopsis:
1263     #include <petscsys.h>
1264    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)
1265 
1266    Not Collective
1267 
1268    Input Parameters:
1269 +   m1 - memory to free
1270 .   m2 - 2nd memory to free
1271 .   m3 - 3rd memory to free
1272 .   m4 - 4th memory to free
1273 .   m5 - 5th memory to free
1274 .   m6 - 6th memory to free
1275 -   m7 - 7th memory to free
1276 
1277 
1278    Level: developer
1279 
1280    Notes:
1281     Memory must have been obtained with PetscMalloc7()
1282 
1283 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1284           PetscMalloc7()
1285 
1286 M*/
1287 #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   PetscFreeA(7,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5),&(m6),&(m7))
1288 
1289 PETSC_EXTERN PetscErrorCode PetscMallocA(int,PetscBool,int,const char *,const char *,size_t,void *,...);
1290 PETSC_EXTERN PetscErrorCode PetscFreeA(int,int,const char *,const char *,void *,...);
1291 PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,PetscBool,int,const char[],const char[],void**);
1292 PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1293 PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t,int,const char[],const char[],void**);
1294 PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1295 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 **));
1296 PETSC_EXTERN PetscErrorCode PetscMallocClear(void);
1297 
1298 /*
1299   Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1300 */
1301 PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1302 PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);
1303 #if defined(PETSC_HAVE_CUDA)
1304 PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void);
1305 PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void);
1306 #endif
1307 #if defined(PETSC_HAVE_HIP)
1308 PETSC_EXTERN PetscErrorCode PetscMallocSetHIPHost(void);
1309 PETSC_EXTERN PetscErrorCode PetscMallocResetHIPHost(void);
1310 #endif
1311 
1312 
1313 #define MPIU_PETSCLOGDOUBLE  MPI_DOUBLE
1314 #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION
1315 
1316 /*
1317    Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1318 */
1319 PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1320 PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1321 PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1322 PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1323 PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1324 PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int,PetscLogDouble*);
1325 PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool,PetscBool);
1326 PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*,PetscBool*,PetscBool*);
1327 PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1328 PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1329 PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool*);
1330 PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool);
1331 PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool*);
1332 
1333 PETSC_EXTERN const char *const PetscDataTypes[];
1334 PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1335 PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1336 PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1337 PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);
1338 
1339 /*
1340     Basic memory and string operations. These are usually simple wrappers
1341    around the basic Unix system calls, but a few of them have additional
1342    functionality and/or error checking.
1343 */
1344 PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1345 PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1346 PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***);
1347 PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1348 PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool  *);
1349 PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool  *);
1350 PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1351 PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1352 PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1353 PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1354 PETSC_EXTERN PetscErrorCode PetscStrlcat(char[],const char[],size_t);
1355 PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1356 PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1357 PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1358 PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1359 PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1360 PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1361 PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1362 PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1363 PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*);
1364 PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1365 PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1366 PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1367 PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1368 PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt,const char *const*,char***);
1369 PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt,char***);
1370 PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);
1371 
1372 PETSC_EXTERN void PetscStrcmpNoError(const char[],const char[],PetscBool  *);
1373 
1374 PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1375 PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1376 PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);
1377 
1378 PETSC_EXTERN PetscErrorCode PetscStrInList(const char[],const char[],char,PetscBool*);
1379 PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1380 PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);
1381 
1382 /*
1383    These are MPI operations for MPI_Allreduce() etc
1384 */
1385 PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1386 #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1387 PETSC_EXTERN MPI_Op MPIU_SUM;
1388 PETSC_EXTERN MPI_Op MPIU_MAX;
1389 PETSC_EXTERN MPI_Op MPIU_MIN;
1390 #else
1391 #define MPIU_SUM MPI_SUM
1392 #define MPIU_MAX MPI_MAX
1393 #define MPIU_MIN MPI_MIN
1394 #endif
1395 PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);
1396 
1397 PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1398 PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1399 
1400 PETSC_EXTERN const char *const PetscFileModes[];
1401 
1402 /*
1403     Defines PETSc error handling.
1404 */
1405 #include <petscerror.h>
1406 
1407 #define PETSC_SMALLEST_CLASSID  1211211
1408 PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1409 PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1410 PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);
1411 PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject,PetscObjectId*);
1412 PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject,PetscObjectId,PetscBool*);
1413 
1414 
1415 /*
1416    Routines that get memory usage information from the OS
1417 */
1418 PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1419 PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1420 PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1421 PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);
1422 
1423 PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);
1424 
1425 /*
1426    Initialization of PETSc
1427 */
1428 PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1429 PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1430 PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1431 PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1432 PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1433 PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1434 PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1435 PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1436 PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1437 PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);
1438 
1439 PETSC_EXTERN PetscErrorCode PetscEnd(void);
1440 PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);
1441 
1442 PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1443 PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1444 PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1445 PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);
1446 
1447 PETSC_EXTERN PetscErrorCode PetscMonitorCompare(PetscErrorCode (*)(void),void *,PetscErrorCode (*)(void**),PetscErrorCode (*)(void),void *,PetscErrorCode (*)(void**),PetscBool *);
1448 
1449 
1450 /*
1451      These are so that in extern C code we can caste function pointers to non-extern C
1452    function pointers. Since the regular C++ code expects its function pointers to be C++
1453 */
1454 PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1455 PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1456 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);
1457 
1458 /*
1459     Functions that can act on any PETSc object.
1460 */
1461 PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1462 PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1463 PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1464 PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1465 PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1466 PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1467 PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1468 PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1469 PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1470 PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1471 PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1472 PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1473 PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1474 PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1475 PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt*);
1476 PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1477 PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1478 PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject*);
1479 PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void));
1480 #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d))
1481 PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1482 PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1483 PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1484 PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject,PetscObject);
1485 PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt*);
1486 
1487 #include <petscviewertypes.h>
1488 #include <petscoptions.h>
1489 
1490 PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer,PetscBool,PetscLogDouble);
1491 PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool*);
1492 
1493 PETSC_EXTERN PetscErrorCode PetscObjectsListGetGlobalNumbering(MPI_Comm,PetscInt,PetscObject*,PetscInt*,PetscInt*);
1494 
1495 PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer,const char[]);
1496 PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1497 PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1498 #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1499 PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1500 PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1501 PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1502 PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1503 PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1504 PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1505 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1506 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1507 PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,PetscObject,const char[]);
1508 PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1509 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1510 PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject,const char[],PetscBool *);
1511 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1512 PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1513 PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1514 PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);
1515 
1516 #if defined(PETSC_HAVE_SAWS)
1517 PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1518 PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1519 PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1520 PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1521 PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1522 PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1523 PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1524 PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1525 PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1526 PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);
1527 
1528 #else
1529 #define PetscSAWsBlock()                        0
1530 #define PetscObjectSAWsViewOff(obj)             0
1531 #define PetscObjectSAWsSetBlock(obj,flg)        0
1532 #define PetscObjectSAWsBlock(obj)               0
1533 #define PetscObjectSAWsGrantAccess(obj)         0
1534 #define PetscObjectSAWsTakeAccess(obj)          0
1535 #define PetscStackViewSAWs()                    0
1536 #define PetscStackSAWsViewOff()                 0
1537 #define PetscStackSAWsTakeAccess()
1538 #define PetscStackSAWsGrantAccess()
1539 
1540 #endif
1541 
1542 PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1543 PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1544 PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);
1545 PETSC_EXTERN PetscErrorCode PetscDLAddr(void (*)(void), char **);
1546 #ifdef PETSC_HAVE_CXX
1547 PETSC_EXTERN PetscErrorCode PetscDemangleSymbol(const char *, char **);
1548 #endif
1549 
1550 #if defined(PETSC_USE_DEBUG)
1551 PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1552 #endif
1553 PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);
1554 
1555 PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1556 PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1557 PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1558 PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1559 PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1560 PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);
1561 
1562 /*
1563     Dynamic library lists. Lists of names of routines in objects or in dynamic
1564   link libraries that will be loaded as needed.
1565 */
1566 
1567 #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1568 PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1569 PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1570 #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1571 PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1572 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[],const char[]);
1573 PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1574 PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1575 PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);
1576 
1577 PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1578 PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1579 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1580 PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1581 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1582 PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1583 PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1584 PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);
1585 
1586 /*
1587      Useful utility routines
1588 */
1589 PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1590 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1591 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm,PetscInt*,PetscInt*);
1592 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1593 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1594 PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1595 PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);
1596 PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm,PetscInt[2],PetscInt[2]);
1597 PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm,PetscReal[2],PetscReal[2]);
1598 
1599 /*MC
1600     PetscNot - negates a logical type value and returns result as a PetscBool
1601 
1602     Notes:
1603     This is useful in cases like
1604 $     int        *a;
1605 $     PetscBool  flag = PetscNot(a)
1606      where !a would not return a PetscBool because we cannot provide a cast from int to PetscBool in C.
1607 
1608     Level: beginner
1609 
1610     .seealso : PetscBool, PETSC_TRUE, PETSC_FALSE
1611 M*/
1612 #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)
1613 
1614 /*MC
1615     PetscHelpPrintf - Prints help messages.
1616 
1617    Synopsis:
1618     #include <petscsys.h>
1619      PetscErrorCode (*PetscHelpPrintf)(MPI_Comm comm, const char format[],args);
1620 
1621     Collective on comm
1622 
1623     Input Parameters:
1624 +  comm - the MPI communicator over which the help message is printed
1625 .  format - the usual printf() format string
1626 -  args - arguments to be printed
1627 
1628    Level: developer
1629 
1630    Fortran Note:
1631      This routine is not supported in Fortran.
1632 
1633    Note:
1634      You can change how help messages are printed by replacing the function pointer with a function that does not simply write to stdout.
1635 
1636       To use, write your own function, for example,
1637 $PetscErrorCode mypetschelpprintf(MPI_Comm comm,const char format[],....)
1638 ${
1639 $ PetscFunctionReturn(0);
1640 $}
1641 then do the assigment
1642 $    PetscHelpPrintf = mypetschelpprintf;
1643    You can do the assignment before PetscInitialize().
1644 
1645   The default routine used is called PetscHelpPrintfDefault().
1646 
1647 .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1648 M*/
1649 PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);
1650 
1651 /*
1652      Defines PETSc profiling.
1653 */
1654 #include <petsclog.h>
1655 
1656 /*
1657       Simple PETSc parallel IO for ASCII printing
1658 */
1659 PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1660 PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1661 PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1662 PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1663 PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1664 PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1665 PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);
1666 PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[],size_t,const char*,PetscInt,const PetscReal[]);
1667 
1668 PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1669 PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1670 PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);
1671 
1672 PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char*,size_t*);
1673 PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char*,char *);
1674 
1675 #if defined(PETSC_HAVE_POPEN)
1676 PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1677 PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*);
1678 PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1679 #endif
1680 
1681 PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1682 PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1683 PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*);
1684 PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1685 PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1686 PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1687 PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);
1688 
1689 PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1690 PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void**);
1691 PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void*);
1692 PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1693 PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer*);
1694 PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer,PetscErrorCode (*)(void*));
1695 PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void*);
1696 
1697 /*
1698    For use in debuggers
1699 */
1700 PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1701 PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1702 PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1703 PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1704 PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);
1705 
1706 #include <stddef.h>
1707 #include <string.h>             /* for memcpy, memset */
1708 #include <stdlib.h>
1709 
1710 #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1711 #include <xmmintrin.h>
1712 #endif
1713 
1714 /*@C
1715    PetscMemmove - Copies n bytes, beginning at location b, to the space
1716    beginning at location a. Copying  between regions that overlap will
1717    take place correctly. Use PetscMemcpy() if the locations do not overlap
1718 
1719    Not Collective
1720 
1721    Input Parameters:
1722 +  b - pointer to initial memory space
1723 .  a - pointer to copy space
1724 -  n - length (in bytes) of space to copy
1725 
1726    Level: intermediate
1727 
1728    Note:
1729    PetscArraymove() is preferred
1730    This routine is analogous to memmove().
1731 
1732    Developers Note: This is inlined for performance
1733 
1734 .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscStrallocpy(),
1735           PetscArraymove()
1736 @*/
1737 PETSC_STATIC_INLINE PetscErrorCode PetscMemmove(void *a,const void *b,size_t n)
1738 {
1739   PetscFunctionBegin;
1740   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to null pointer");
1741   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1742 #if !defined(PETSC_HAVE_MEMMOVE)
1743   if (a < b) {
1744     if (a <= b - n) memcpy(a,b,n);
1745     else {
1746       memcpy(a,b,(int)(b - a));
1747       PetscMemmove(b,b + (int)(b - a),n - (int)(b - a));
1748     }
1749   } else {
1750     if (b <= a - n) memcpy(a,b,n);
1751     else {
1752       memcpy(b + n,b + (n - (int)(a - b)),(int)(a - b));
1753       PetscMemmove(a,b,n - (int)(a - b));
1754     }
1755   }
1756 #else
1757   memmove((char*)(a),(char*)(b),n);
1758 #endif
1759   PetscFunctionReturn(0);
1760 }
1761 
1762 /*@C
1763    PetscMemcpy - Copies n bytes, beginning at location b, to the space
1764    beginning at location a. The two memory regions CANNOT overlap, use
1765    PetscMemmove() in that case.
1766 
1767    Not Collective
1768 
1769    Input Parameters:
1770 +  b - pointer to initial memory space
1771 -  n - length (in bytes) of space to copy
1772 
1773    Output Parameter:
1774 .  a - pointer to copy space
1775 
1776    Level: intermediate
1777 
1778    Compile Option:
1779     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1780                                   for memory copies on double precision values.
1781     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1782                                   for memory copies on double precision values.
1783     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1784                                   for memory copies on double precision values.
1785 
1786    Note:
1787    Prefer PetscArraycpy()
1788    This routine is analogous to memcpy().
1789    Not available from Fortran
1790 
1791    Developer Note: this is inlined for fastest performance
1792 
1793 .seealso: PetscMemzero(), PetscMemcmp(), PetscArrayzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy()
1794 
1795 @*/
1796 PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n)
1797 {
1798 #if defined(PETSC_USE_DEBUG)
1799   size_t al = (size_t) a,bl = (size_t) b;
1800   size_t nl = (size_t) n;
1801   PetscFunctionBegin;
1802   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1803   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1804 #else
1805   PetscFunctionBegin;
1806 #endif
1807   if (a != b && n > 0) {
1808 #if defined(PETSC_USE_DEBUG)
1809     if ((al > bl && (al - bl) < nl) || (bl - al) < nl)  SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1810               or make sure your copy regions and lengths are correct. \n\
1811               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1812 #endif
1813 #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1814    if (!(a % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1815       size_t len = n/sizeof(PetscScalar);
1816 #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1817       PetscBLASInt   one = 1,blen;
1818       PetscErrorCode ierr;
1819       ierr = PetscBLASIntCast(len,&blen);CHKERRQ(ierr);
1820       PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one));
1821 #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1822       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1823 #else
1824       size_t      i;
1825       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1826       for (i=0; i<len; i++) y[i] = x[i];
1827 #endif
1828     } else {
1829       memcpy((char*)(a),(char*)(b),n);
1830     }
1831 #else
1832     memcpy((char*)(a),(char*)(b),n);
1833 #endif
1834   }
1835   PetscFunctionReturn(0);
1836 }
1837 
1838 /*@C
1839    PetscMemzero - Zeros the specified memory.
1840 
1841    Not Collective
1842 
1843    Input Parameters:
1844 +  a - pointer to beginning memory location
1845 -  n - length (in bytes) of memory to initialize
1846 
1847    Level: intermediate
1848 
1849    Compile Option:
1850    PETSC_PREFER_BZERO - on certain machines (the IBM RS6000) the bzero() routine happens
1851   to be faster than the memset() routine. This flag causes the bzero() routine to be used.
1852 
1853    Not available from Fortran
1854    Prefer PetscArrayzero()
1855 
1856    Developer Note: this is inlined for fastest performance
1857 
1858 .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy()
1859 @*/
1860 PETSC_STATIC_INLINE PetscErrorCode PetscMemzero(void *a,size_t n)
1861 {
1862   if (n > 0) {
1863 #if defined(PETSC_USE_DEBUG)
1864     if (!a) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer with %zu bytes",n);
1865 #endif
1866 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1867     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1868       size_t      i,len = n/sizeof(PetscScalar);
1869       PetscScalar *x = (PetscScalar*)a;
1870       for (i=0; i<len; i++) x[i] = 0.0;
1871     } else {
1872 #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1873     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1874       PetscInt len = n/sizeof(PetscScalar);
1875       fortranzero_(&len,(PetscScalar*)a);
1876     } else {
1877 #endif
1878 #if defined(PETSC_PREFER_BZERO)
1879       bzero((char *)a,n);
1880 #else
1881       memset((char*)a,0,n);
1882 #endif
1883 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1884     }
1885 #endif
1886   }
1887   return 0;
1888 }
1889 
1890 /*MC
1891    PetscArraycmp - Compares two arrays in memory.
1892 
1893    Synopsis:
1894     #include <petscsys.h>
1895     PetscErrorCode PetscArraycmp(const anytype *str1,const anytype *str2,size_t cnt,PetscBool *e)
1896 
1897    Not Collective
1898 
1899    Input Parameters:
1900 +  str1 - First array
1901 .  str2 - Second array
1902 -  cnt  - Count of the array, not in bytes, but number of entries in the arrays
1903 
1904    Output Parameters:
1905 .   e - PETSC_TRUE if equal else PETSC_FALSE.
1906 
1907    Level: intermediate
1908 
1909    Note:
1910    This routine is a preferred replacement to PetscMemcmp()
1911    The arrays must be of the same type
1912 
1913 .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy(),
1914           PetscArraymove()
1915 M*/
1916 #define  PetscArraycmp(str1,str2,cnt,e) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcmp(str1,str2,(size_t)(cnt)*sizeof(*(str1)),e))
1917 
1918 /*MC
1919    PetscArraymove - Copies from one array in memory to another, the arrays may overlap. Use PetscArraycpy() when the arrays
1920                     do not overlap
1921 
1922    Synopsis:
1923     #include <petscsys.h>
1924     PetscErrorCode PetscArraymove(anytype *str1,const anytype *str2,size_t cnt)
1925 
1926    Not Collective
1927 
1928    Input Parameters:
1929 +  str1 - First array
1930 .  str2 - Second array
1931 -  cnt  - Count of the array, not in bytes, but number of entries in the arrays
1932 
1933    Level: intermediate
1934 
1935    Note:
1936    This routine is a preferred replacement to PetscMemmove()
1937    The arrays must be of the same type
1938 
1939 .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycpy(), PetscMemmove(), PetscArraycmp(), PetscStrallocpy()
1940 M*/
1941 #define  PetscArraymove(str1,str2,cnt) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemmove(str1,str2,(size_t)(cnt)*sizeof(*(str1))))
1942 
1943 /*MC
1944    PetscArraycpy - Copies from one array in memory to another
1945 
1946    Synopsis:
1947     #include <petscsys.h>
1948     PetscErrorCode PetscArraycpy(anytype *str1,const anytype *str2,size_t cnt)
1949 
1950    Not Collective
1951 
1952    Input Parameters:
1953 +  str1 - First array (destination)
1954 .  str2 - Second array (source)
1955 -  cnt  - Count of the array, not in bytes, but number of entries in the arrays
1956 
1957    Level: intermediate
1958 
1959    Note:
1960    This routine is a preferred replacement to PetscMemcpy()
1961    The arrays must be of the same type
1962 
1963 .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraymove(), PetscMemmove(), PetscArraycmp(), PetscStrallocpy()
1964 M*/
1965 #define  PetscArraycpy(str1,str2,cnt) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcpy(str1,str2,(size_t)(cnt)*sizeof(*(str1))))
1966 
1967 /*MC
1968    PetscArrayzero - Zeros an array in memory.
1969 
1970    Synopsis:
1971     #include <petscsys.h>
1972     PetscErrorCode PetscArrayzero(anytype *str1,size_t cnt)
1973 
1974    Not Collective
1975 
1976    Input Parameters:
1977 +  str1 - array
1978 -  cnt  - Count of the array, not in bytes, but number of entries in the array
1979 
1980    Level: intermediate
1981 
1982    Note:
1983    This routine is a preferred replacement to PetscMemzero()
1984 
1985 .seealso: PetscMemcpy(), PetscMemcmp(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy(), PetscArraymove()
1986 M*/
1987 #define  PetscArrayzero(str1,cnt) PetscMemzero(str1,(size_t)(cnt)*sizeof(*(str1)))
1988 
1989 /*MC
1990    PetscPrefetchBlock - Prefetches a block of memory
1991 
1992    Synopsis:
1993     #include <petscsys.h>
1994     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)
1995 
1996    Not Collective
1997 
1998    Input Parameters:
1999 +  a - pointer to first element to fetch (any type but usually PetscInt or PetscScalar)
2000 .  n - number of elements to fetch
2001 .  rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors)
2002 -  t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note
2003 
2004    Level: developer
2005 
2006    Notes:
2007    The last two arguments (rw and t) must be compile-time constants.
2008 
2009    Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality.  Not all architectures offer
2010    equivalent locality hints, but the following macros are always defined to their closest analogue.
2011 +  PETSC_PREFETCH_HINT_NTA - Non-temporal.  Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched).
2012 .  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.
2013 .  PETSC_PREFETCH_HINT_T1 - Fetch to level 2 and higher (not L1).
2014 -  PETSC_PREFETCH_HINT_T2 - Fetch to high-level cache only.  (On many systems, T0 and T1 are equivalent.)
2015 
2016    This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid
2017    address).
2018 
2019 M*/
2020 #define PetscPrefetchBlock(a,n,rw,t) do {                               \
2021     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
2022     for (; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
2023   } while (0)
2024 
2025 /*
2026       Determine if some of the kernel computation routines use
2027    Fortran (rather than C) for the numerical calculations. On some machines
2028    and compilers (like complex numbers) the Fortran version of the routines
2029    is faster than the C/C++ versions. The flag --with-fortran-kernels
2030    should be used with ./configure to turn these on.
2031 */
2032 #if defined(PETSC_USE_FORTRAN_KERNELS)
2033 
2034 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
2035 #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
2036 #endif
2037 
2038 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
2039 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
2040 #endif
2041 
2042 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
2043 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
2044 #endif
2045 
2046 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
2047 #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
2048 #endif
2049 
2050 #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
2051 #define PETSC_USE_FORTRAN_KERNEL_NORM
2052 #endif
2053 
2054 #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
2055 #define PETSC_USE_FORTRAN_KERNEL_MAXPY
2056 #endif
2057 
2058 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
2059 #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
2060 #endif
2061 
2062 #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
2063 #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
2064 #endif
2065 
2066 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
2067 #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
2068 #endif
2069 
2070 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
2071 #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
2072 #endif
2073 
2074 #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
2075 #define PETSC_USE_FORTRAN_KERNEL_MDOT
2076 #endif
2077 
2078 #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
2079 #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
2080 #endif
2081 
2082 #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
2083 #define PETSC_USE_FORTRAN_KERNEL_AYPX
2084 #endif
2085 
2086 #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
2087 #define PETSC_USE_FORTRAN_KERNEL_WAXPY
2088 #endif
2089 
2090 #endif
2091 
2092 /*
2093     Macros for indicating code that should be compiled with a C interface,
2094    rather than a C++ interface. Any routines that are dynamically loaded
2095    (such as the PCCreate_XXX() routines) must be wrapped so that the name
2096    mangler does not change the functions symbol name. This just hides the
2097    ugly extern "C" {} wrappers.
2098 */
2099 #if defined(__cplusplus)
2100 #  define EXTERN_C_BEGIN extern "C" {
2101 #  define EXTERN_C_END }
2102 #else
2103 #  define EXTERN_C_BEGIN
2104 #  define EXTERN_C_END
2105 #endif
2106 
2107 /* --------------------------------------------------------------------*/
2108 
2109 /*MC
2110     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
2111         communication
2112 
2113    Level: beginner
2114 
2115    Note: This manual page is a place-holder because MPICH does not have a manual page for MPI_Comm
2116 
2117 .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
2118 M*/
2119 
2120 #if defined(PETSC_HAVE_MPIIO)
2121 PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2122 PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2123 PETSC_EXTERN PetscErrorCode MPIU_File_write_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2124 PETSC_EXTERN PetscErrorCode MPIU_File_read_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2125 PETSC_EXTERN PetscErrorCode MPIU_File_write_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2126 PETSC_EXTERN PetscErrorCode MPIU_File_read_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2127 #endif
2128 
2129 /* the following petsc_static_inline require petscerror.h */
2130 
2131 /* Limit MPI to 32-bits */
2132 #define PETSC_MPI_INT_MAX  2147483647
2133 #define PETSC_MPI_INT_MIN -2147483647
2134 /* Limit BLAS to 32-bits */
2135 #define PETSC_BLAS_INT_MAX  2147483647
2136 #define PETSC_BLAS_INT_MIN -2147483647
2137 #define PETSC_CUBLAS_INT_MAX  2147483647
2138 
2139 /*@C
2140     PetscIntCast - casts a PetscInt64 (which is 64 bits in size) to a PetscInt (which may be 32 bits in size), generates an
2141          error if the PetscInt is not large enough to hold the number.
2142 
2143    Not Collective
2144 
2145    Input Parameter:
2146 .     a - the PetscInt64 value
2147 
2148    Output Parameter:
2149 .     b - the resulting PetscInt value
2150 
2151    Level: advanced
2152 
2153    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
2154 
2155    Not available from Fortran
2156 
2157 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscBLASIntCast(), PetscIntMultError(), PetscIntSumError()
2158 @*/
2159 PETSC_STATIC_INLINE PetscErrorCode PetscIntCast(PetscInt64 a,PetscInt *b)
2160 {
2161   PetscFunctionBegin;
2162 #if !defined(PETSC_USE_64BIT_INDICES)
2163   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");
2164 #endif
2165   *b =  (PetscInt)(a);
2166   PetscFunctionReturn(0);
2167 }
2168 
2169 /*@C
2170     PetscBLASIntCast - casts a PetscInt (which may be 64 bits in size) to a PetscBLASInt (which may be 32 bits in size), generates an
2171          error if the PetscBLASInt is not large enough to hold the number.
2172 
2173    Not Collective
2174 
2175    Input Parameter:
2176 .     a - the PetscInt value
2177 
2178    Output Parameter:
2179 .     b - the resulting PetscBLASInt value
2180 
2181    Level: advanced
2182 
2183    Notes:
2184       Not available from Fortran
2185       Errors if the integer is negative since PETSc calls to BLAS/LAPACK never need to cast negative integer inputs
2186 
2187 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscIntCast()
2188 @*/
2189 PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
2190 {
2191   PetscFunctionBegin;
2192 #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
2193   *b = 0;
2194   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");
2195 #endif
2196   if (a < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Passing negative integer to BLAS/LAPACK routine");
2197   *b =  (PetscBLASInt)(a);
2198   PetscFunctionReturn(0);
2199 }
2200 
2201 /*@C
2202     PetscCuBLASIntCast - like PetscBLASIntCast(), but for PetscCuBLASInt.
2203 
2204    Not Collective
2205 
2206    Input Parameter:
2207 .     a - the PetscInt value
2208 
2209    Output Parameter:
2210 .     b - the resulting PetscCuBLASInt value
2211 
2212    Level: advanced
2213 
2214    Notes:
2215       Errors if the integer is negative since PETSc calls to cuBLAS and friends never need to cast negative integer inputs
2216 
2217 .seealso: PetscCuBLASInt, PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscMPIIntCast(), PetscIntCast()
2218 @*/
2219 PETSC_STATIC_INLINE PetscErrorCode PetscCuBLASIntCast(PetscInt a,PetscBLASInt *b)
2220 {
2221   PetscFunctionBegin;
2222 #if defined(PETSC_USE_64BIT_INDICES)
2223   *b = 0;
2224   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.");
2225 #endif
2226   if (a < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Passing negative integer to cuBLAS routine");
2227   *b =  (PetscCuBLASInt)(a);
2228   PetscFunctionReturn(0);
2229 }
2230 
2231 /*@C
2232     PetscMPIIntCast - casts a PetscInt (which may be 64 bits in size) to a PetscMPIInt (which may be 32 bits in size), generates an
2233          error if the PetscMPIInt is not large enough to hold the number.
2234 
2235    Not Collective
2236 
2237    Input Parameter:
2238 .     a - the PetscInt value
2239 
2240    Output Parameter:
2241 .     b - the resulting PetscMPIInt value
2242 
2243    Level: advanced
2244 
2245    Not available from Fortran
2246 
2247 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntCast()
2248 @*/
2249 PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
2250 {
2251   PetscFunctionBegin;
2252 #if defined(PETSC_USE_64BIT_INDICES)
2253   *b = 0;
2254   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");
2255 #endif
2256   *b =  (PetscMPIInt)(a);
2257   PetscFunctionReturn(0);
2258 }
2259 
2260 #define PetscInt64Mult(a,b)   ((PetscInt64)(a))*((PetscInt64)(b))
2261 
2262 /*@C
2263 
2264    PetscRealIntMultTruncate - Computes the product of a positive PetscReal and a positive PetscInt and truncates the value to slightly less than the maximal possible value
2265 
2266    Not Collective
2267 
2268    Input Parameter:
2269 +     a - the PetscReal value
2270 -     b - the second value
2271 
2272    Returns:
2273       the result as a PetscInt value
2274 
2275    Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2276    Use PetscIntMultTruncate() to compute the product of two positive PetscInt and truncate to fit a PetscInt
2277    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
2278 
2279    Developers Note:
2280    We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks.
2281 
2282    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2283 
2284    Not available from Fortran
2285 
2286    Level: advanced
2287 
2288 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError(), PetscIntSumError()
2289 @*/
2290 PETSC_STATIC_INLINE PetscInt PetscRealIntMultTruncate(PetscReal a,PetscInt b)
2291 {
2292   PetscInt64 r;
2293 
2294   r  =  (PetscInt64) (a*(PetscReal)b);
2295   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2296   return (PetscInt) r;
2297 }
2298 
2299 /*@C
2300 
2301    PetscIntMultTruncate - Computes the product of two positive PetscInt and truncates the value to slightly less than the maximal possible value
2302 
2303    Not Collective
2304 
2305    Input Parameter:
2306 +     a - the PetscInt value
2307 -     b - the second value
2308 
2309    Returns:
2310       the result as a PetscInt value
2311 
2312    Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2313    Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2314    Use PetscIntMultError() to compute the product of two PetscInt if you wish to generate an error if the result will not fit in a PetscInt
2315 
2316    Not available from Fortran
2317 
2318    Developers Note: We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks.
2319 
2320    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2321 
2322    Level: advanced
2323 
2324 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError(), PetscIntSumError()
2325 @*/
2326 PETSC_STATIC_INLINE PetscInt PetscIntMultTruncate(PetscInt a,PetscInt b)
2327 {
2328   PetscInt64 r;
2329 
2330   r  =  PetscInt64Mult(a,b);
2331   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2332   return (PetscInt) r;
2333 }
2334 
2335 /*@C
2336 
2337    PetscIntSumTruncate - Computes the sum of two positive PetscInt and truncates the value to slightly less than the maximal possible value
2338 
2339    Not Collective
2340 
2341    Input Parameter:
2342 +     a - the PetscInt value
2343 -     b - the second value
2344 
2345    Returns:
2346      the result as a PetscInt value
2347 
2348    Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2349    Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2350    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
2351 
2352    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2353 
2354    Not available from Fortran
2355 
2356    Level: advanced
2357 
2358 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError()
2359 @*/
2360 PETSC_STATIC_INLINE PetscInt PetscIntSumTruncate(PetscInt a,PetscInt b)
2361 {
2362   PetscInt64 r;
2363 
2364   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2365   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2366   return (PetscInt) r;
2367 }
2368 
2369 /*@C
2370 
2371    PetscIntMultError - Computes the product of two positive PetscInt and generates an error with overflow.
2372 
2373    Not Collective
2374 
2375    Input Parameter:
2376 +     a - the PetscInt value
2377 -     b - the second value
2378 
2379    Output Parameter:ma
2380 .     result - the result as a PetscInt value, or NULL if you do not want the result, you just want to check if it overflows
2381 
2382    Use PetscInt64Mult() to compute the product of two 32 bit PetscInt and store in a PetscInt64
2383    Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt
2384 
2385    Not available from Fortran
2386 
2387    Developers Note: We currently assume that PetscInt addition does not overflow, this is obviously wrong but requires many more checks.
2388 
2389    Level: advanced
2390 
2391 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64(), PetscIntSumError()
2392 @*/
2393 PETSC_STATIC_INLINE PetscErrorCode PetscIntMultError(PetscInt a,PetscInt b,PetscInt *result)
2394 {
2395   PetscInt64 r;
2396 
2397   PetscFunctionBegin;
2398   r  =  PetscInt64Mult(a,b);
2399 #if !defined(PETSC_USE_64BIT_INDICES)
2400   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);
2401 #endif
2402   if (result) *result = (PetscInt) r;
2403   PetscFunctionReturn(0);
2404 }
2405 
2406 /*@C
2407 
2408    PetscIntSumError - Computes the sum of two positive PetscInt and generates an error with overflow.
2409 
2410    Not Collective
2411 
2412    Input Parameter:
2413 +     a - the PetscInt value
2414 -     b - the second value
2415 
2416    Output Parameter:ma
2417 .     c - the result as a PetscInt value,  or NULL if you do not want the result, you just want to check if it overflows
2418 
2419    Use PetscInt64Mult() to compute the product of two 32 bit PetscInt and store in a PetscInt64
2420    Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt
2421 
2422    Not available from Fortran
2423 
2424    Level: advanced
2425 
2426 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult(), PetscIntMultError()
2427 @*/
2428 PETSC_STATIC_INLINE PetscErrorCode PetscIntSumError(PetscInt a,PetscInt b,PetscInt *result)
2429 {
2430   PetscInt64 r;
2431 
2432   PetscFunctionBegin;
2433   r  =  ((PetscInt64)a) + ((PetscInt64)b);
2434 #if !defined(PETSC_USE_64BIT_INDICES)
2435   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);
2436 #endif
2437   if (result) *result = (PetscInt) r;
2438   PetscFunctionReturn(0);
2439 }
2440 
2441 /*
2442      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2443 */
2444 #if defined(hz)
2445 #  undef hz
2446 #endif
2447 
2448 #include <limits.h>
2449 
2450 /* The number of bits in a byte */
2451 
2452 #define PETSC_BITS_PER_BYTE CHAR_BIT
2453 
2454 #if defined(PETSC_HAVE_SYS_TYPES_H)
2455 #  include <sys/types.h>
2456 #endif
2457 
2458 /*MC
2459 
2460     PETSC_VERSION - This manual page provides information about how PETSc documents and uses its version information. This information is available to both C/C++
2461                     and Fortran compilers when petscsys.h is included.
2462 
2463 
2464     The current PETSc version and the API for accessing it are defined in petscversion.h
2465 
2466     The complete version number is given as the triple  PETSC_VERSION_MAJOR.PETSC_VERSION_MINOR.PETSC_VERSION_SUBMINOR (in short hand x.y.z)
2467 
2468     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
2469     where only a change in the major version number (x) indicates a change in the API.
2470 
2471     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
2472 
2473     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),
2474     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
2475     version number (x.y.z).
2476 
2477     PETSC_RELEASE_DATE is the date the x.y version was released (i.e. the version before any patch releases)
2478 
2479     PETSC_VERSION_DATE is the date the x.y.z version was released
2480 
2481     PETSC_VERSION_GIT is the last git commit to the repository given in the form vx.y.z-wwwww
2482 
2483     PETSC_VERSION_DATE_GIT is the date of the last git commit to the repository
2484 
2485     Level: intermediate
2486 
2487     PETSC_VERSION_() and PETSC_VERSION_PATCH are deprecated and will eventually be removed. For several releases PETSC_VERSION_PATCH is always 0
2488 
2489 M*/
2490 
2491 PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2492 PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2493 PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2494 PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2495 PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2496 PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2497 PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2498 PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt*,PetscInt*,PetscInt*,PetscInt*);
2499 
2500 PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt,const PetscInt[],PetscBool*);
2501 PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt,const PetscMPIInt[],PetscBool*);
2502 PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt,const PetscReal[],PetscBool*);
2503 PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2504 PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt,PetscInt[]);
2505 PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt*,PetscInt[]);
2506 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2507 PETSC_EXTERN PetscErrorCode PetscCheckDupsInt(PetscInt,const PetscInt[],PetscBool*);
2508 PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2509 PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt*);
2510 PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2511 PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2512 PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2513 PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt[],PetscInt[],PetscInt[]);
2514 PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2515 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2516 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2517 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt,PetscMPIInt[],PetscInt[]);
2518 PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2519 PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt,PetscInt[],void*,size_t,void*);
2520 PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2521 PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt,PetscReal[],PetscInt[]);
2522 PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2523 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2524 PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal,PetscInt,const PetscReal[],PetscReal,PetscInt*);
2525 PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2526 PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2527 PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2528 PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt[],const PetscInt[],PetscInt,const PetscInt[],const PetscInt[],PetscInt*,PetscInt**,PetscInt**);
2529 PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt*,PetscInt**);
2530 PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt,const PetscMPIInt[],PetscInt,const PetscMPIInt[],PetscInt*,PetscMPIInt**);
2531 PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *);
2532 
2533 PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt,void*,size_t,int(*)(const void*,const void*,void*),void*);
2534 PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt,PetscInt[]);
2535 PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt,PetscMPIInt[]);
2536 PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt,PetscReal[]);
2537 PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt,void*,size_t,void*,size_t,int(*)(const void*,const void*,void*),void*);
2538 PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt,PetscInt[],PetscInt[]);
2539 PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt,PetscMPIInt[],PetscMPIInt[]);
2540 PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt,PetscReal[],PetscInt[]);
2541 
2542 PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2543 PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);
2544 
2545 /*J
2546     PetscRandomType - String with the name of a PETSc randomizer
2547 
2548    Level: beginner
2549 
2550    Notes:
2551    To use SPRNG or RANDOM123 you must have ./configure PETSc
2552    with the option --download-sprng or --download-random123
2553 
2554 .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2555 J*/
2556 typedef const char* PetscRandomType;
2557 #define PETSCRAND       "rand"
2558 #define PETSCRAND48     "rand48"
2559 #define PETSCSPRNG      "sprng"
2560 #define PETSCRANDER48   "rander48"
2561 #define PETSCRANDOM123  "random123"
2562 #define PETSCCURAND     "curand"
2563 
2564 /* Logging support */
2565 PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;
2566 
2567 PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);
2568 
2569 /* Dynamic creation and loading functions */
2570 PETSC_EXTERN PetscFunctionList PetscRandomList;
2571 
2572 PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2573 PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom,PetscRandomType);
2574 PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2575 PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom,PetscRandomType*);
2576 PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom,PetscObject,const char[]);
2577 PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);
2578 
2579 PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2580 PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2581 PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2582 PETSC_EXTERN PetscErrorCode PetscRandomGetValues(PetscRandom,PetscInt,PetscScalar*);
2583 PETSC_EXTERN PetscErrorCode PetscRandomGetValuesReal(PetscRandom,PetscInt,PetscReal*);
2584 PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2585 PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2586 PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2587 PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2588 PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2589 PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);
2590 
2591 PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2592 PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2593 PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2594 PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2595 PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2596 PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2597 PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);
2598 PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2599 PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2600 PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);
2601 
2602 PETSC_STATIC_INLINE PetscBool PetscBinaryBigEndian(void) {long _petsc_v = 1; return ((char*)&_petsc_v)[0] ? PETSC_FALSE : PETSC_TRUE;}
2603 
2604 PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscInt*,PetscDataType);
2605 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscInt*,PetscDataType);
2606 PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,const void*,PetscInt,PetscDataType);
2607 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,const void*,PetscInt,PetscDataType);
2608 PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2609 PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2610 PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2611 PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2612 PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2613 PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2614 PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2615 #if defined(PETSC_USE_SOCKET_VIEWER)
2616 PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);
2617 #endif
2618 
2619 PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2620 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2621 PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);
2622 
2623 PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2624 PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool);
2625 PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2626 PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2627 PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2628 PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2629 PETSC_EXTERN PetscErrorCode PetscWaitOnError(void);
2630 
2631 PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2632 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2633 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2634 PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2635 PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2636 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt*,const void*,PetscMPIInt*,PetscMPIInt**,void*)
2637   PetscAttrMPIPointerWithType(6,3);
2638 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2639                                                     PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2640                                                     PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2641   PetscAttrMPIPointerWithType(6,3);
2642 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2643                                                        MPI_Request**,MPI_Request**,
2644                                                        PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2645                                                        PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2646   PetscAttrMPIPointerWithType(6,3);
2647 
2648 PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2649 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2650 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);
2651 
2652 PETSC_EXTERN PetscErrorCode PetscSSEIsEnabled(MPI_Comm,PetscBool*,PetscBool*);
2653 
2654 PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);
2655 
2656 PETSC_EXTERN const char *const PetscSubcommTypes[];
2657 
2658 struct _n_PetscSubcomm {
2659   MPI_Comm         parent;           /* parent communicator */
2660   MPI_Comm         dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2661   MPI_Comm         child;            /* the sub-communicator */
2662   PetscMPIInt      n;                /* num of subcommunicators under the parent communicator */
2663   PetscMPIInt      color;            /* color of processors belong to this communicator */
2664   PetscMPIInt      *subsize;         /* size of subcommunicator[color] */
2665   PetscSubcommType type;
2666   char             *subcommprefix;
2667 };
2668 
2669 PETSC_STATIC_INLINE MPI_Comm PetscSubcommParent(PetscSubcomm scomm) {return scomm->parent;}
2670 PETSC_STATIC_INLINE MPI_Comm PetscSubcommChild(PetscSubcomm scomm) {return scomm->child;}
2671 PETSC_STATIC_INLINE MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm) {return scomm->dupparent;}
2672 PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2673 PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2674 PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2675 PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2676 PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2677 PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2678 PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2679 PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm,const char[]);
2680 PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm,MPI_Comm*);
2681 PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm,MPI_Comm*);
2682 PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm,MPI_Comm*);
2683 
2684 PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt,PetscHeap*);
2685 PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap,PetscInt,PetscInt);
2686 PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap,PetscInt*,PetscInt*);
2687 PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap,PetscInt*,PetscInt*);
2688 PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap,PetscInt,PetscInt);
2689 PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2690 PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap*);
2691 PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap,PetscViewer);
2692 
2693 PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2694 PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm,PetscShmComm*);
2695 PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2696 PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2697 PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm,MPI_Comm*);
2698 
2699 /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2700 PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm,PetscInt,PetscOmpCtrl*);
2701 PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl,MPI_Comm*,MPI_Comm*,PetscBool*);
2702 PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl*);
2703 PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2704 PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2705 PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);
2706 
2707 PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2708 PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2709 PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2710 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2711 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2712 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2713 PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2714 PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);
2715 
2716 
2717 /* Type-safe wrapper to encourage use of PETSC_RESTRICT. Does not use PetscFunctionBegin because the error handling
2718  * prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as
2719  * possible. */
2720 PETSC_STATIC_INLINE PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg,size_t count,PetscInt *PETSC_RESTRICT *slot) {return PetscSegBufferGet(seg,count,(void**)slot);}
2721 
2722 extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2723 PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted*);
2724 PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted*);
2725 PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted,const char*,const char*,PetscBool*);
2726 
2727 #include <stdarg.h>
2728 PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list);
2729 PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE*,const char[],va_list);
2730 
2731 PETSC_EXTERN PetscSegBuffer PetscCitationsList;
2732 
2733 /*@C
2734       PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code.
2735 
2736      Not Collective - only what is registered on rank 0 of PETSC_COMM_WORLD will be printed
2737 
2738      Input Parameters:
2739 +      cite - the bibtex item, formated to displayed on multiple lines nicely
2740 -      set - a boolean variable initially set to PETSC_FALSE; this is used to insure only a single registration of the citation
2741 
2742    Level: intermediate
2743 
2744    Not available from Fortran
2745 
2746      Options Database:
2747 .     -citations [filename]   - print out the bibtex entries for the given computation
2748 @*/
2749 PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2750 {
2751   size_t         len;
2752   char           *vstring;
2753   PetscErrorCode ierr;
2754 
2755   PetscFunctionBegin;
2756   if (set && *set) PetscFunctionReturn(0);
2757   ierr = PetscStrlen(cit,&len);CHKERRQ(ierr);
2758   ierr = PetscSegBufferGet(PetscCitationsList,len,&vstring);CHKERRQ(ierr);
2759   ierr = PetscArraycpy(vstring,cit,len);CHKERRQ(ierr);
2760   if (set) *set = PETSC_TRUE;
2761   PetscFunctionReturn(0);
2762 }
2763 
2764 PETSC_EXTERN PetscErrorCode PetscURLShorten(const char[],char[],size_t);
2765 PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm,char[],char[],size_t);
2766 PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm,const char[],char[],size_t);
2767 PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm,const char[],const char []);
2768 
2769 PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm,char[],char[],size_t);
2770 PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm,const char[],char[],char[],size_t);
2771 
2772 PETSC_EXTERN PetscErrorCode PetscGlobusGetTransfers(MPI_Comm,const char[],char[],size_t);
2773 
2774 PETSC_EXTERN PetscErrorCode PetscTextBelt(MPI_Comm,const char[],const char[],PetscBool*);
2775 PETSC_EXTERN PetscErrorCode PetscTellMyCell(MPI_Comm,const char[],const char[],PetscBool*);
2776 
2777 PETSC_EXTERN PetscErrorCode PetscPullJSONValue(const char[],const char[],char[],size_t,PetscBool*);
2778 PETSC_EXTERN PetscErrorCode PetscPushJSONValue(char[],const char[],const char[],size_t);
2779 
2780 
2781 #if defined(PETSC_USE_DEBUG)
2782 PETSC_STATIC_INLINE unsigned int PetscStrHash(const char *str)
2783 {
2784   unsigned int c,hash = 5381;
2785 
2786   while ((c = (unsigned int)*str++)) hash = ((hash << 5) + hash) + c; /* hash * 33 + c */
2787   return hash;
2788 }
2789 
2790 /*MC
2791    MPIU_Allreduce - a PETSc replacement for MPI_Allreduce() that tries to determine if the call from all the MPI processes occur from the
2792                     same place in the PETSc code. This helps to detect bugs where different MPI processes follow different code paths
2793                     resulting in inconsistent and incorrect calls to MPI_Allreduce().
2794 
2795    Collective
2796 
2797    Synopsis:
2798      #include <petscsys.h>
2799      PetscErrorCode MPIU_Allreduce(void *indata,void *outdata,PetscMPIInt count,MPI_Datatype datatype, MPI_Op op, MPI_Comm comm);
2800 
2801    Input Parameters:
2802 +  a - pointer to the input data to be reduced
2803 .  c - the number of MPI data items in a and b
2804 .  d - the MPI datatype, for example MPI_INT
2805 .  e - the MPI operation, for example MPI_SUM
2806 -  fcomm - the MPI communicator on which the operation occurs
2807 
2808    Output Parameter:
2809 .  b - the reduced values
2810 
2811    Notes:
2812      In optimized mode this directly calls MPI_Allreduce()
2813 
2814      This is defined as a macro that can return error codes internally so it cannot be used in a subroutine that returns void.
2815 
2816      The error code this returns should be checked with CHKERRMPI()
2817 
2818    Level: developer
2819 
2820 .seealso: MPI_Allreduce()
2821 M*/
2822 #define MPIU_Allreduce(a,b,c,d,e,fcomm) MPI_SUCCESS; do {\
2823   PetscErrorCode _4_ierr; \
2824   PetscMPIInt a_b1[6],a_b2[6];\
2825   a_b1[0] = -(PetscMPIInt)__LINE__;                          a_b1[1] = -a_b1[0];\
2826   a_b1[2] = -(PetscMPIInt)PetscStrHash(PETSC_FUNCTION_NAME); a_b1[3] = -a_b1[2];\
2827   a_b1[4] = -(PetscMPIInt)(c);                               a_b1[5] = -a_b1[4];\
2828   _4_ierr = MPI_Allreduce(a_b1,a_b2,6,MPI_INT,MPI_MAX,fcomm);CHKERRMPI(_4_ierr);\
2829   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");\
2830   if (-a_b2[2] != a_b2[3]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_Allreduce() called in different locations (functions) on different processors");\
2831   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);\
2832   _4_ierr = MPI_Allreduce((a),(b),(c),d,e,(fcomm));CHKERRMPI(_4_ierr);\
2833   } while (0)
2834 #else
2835 #define MPIU_Allreduce(a,b,c,d,e,fcomm) MPI_Allreduce((a),(b),(c),d,e,(fcomm))
2836 #endif
2837 
2838 #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
2839 PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint,PetscMPIInt,MPI_Info,MPI_Comm,void*,MPI_Win*);
2840 PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win,PetscMPIInt,MPI_Aint*,PetscMPIInt*,void*);
2841 #endif
2842 
2843 /*
2844     List of external packages and queries on it
2845 */
2846 PETSC_EXTERN PetscErrorCode  PetscHasExternalPackage(const char[],PetscBool*);
2847 
2848 #endif
2849