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