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