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