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