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