xref: /petsc/include/petscsys.h (revision 82f73ecaae2ba6158892468c437f14f4e82707f0)
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 /* ========================================================================== */
37 /*
38    This facilitates using the C version of PETSc from C++ and the C++ version from C.
39 */
40 #if defined(__cplusplus)
41 #  define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_CXX
42 #else
43 #  define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_C
44 #endif
45 
46 /* ========================================================================== */
47 /*
48    Since PETSc manages its own extern "C" handling users should never include PETSc include
49    files within extern "C". This will generate a compiler error if a user does put the include
50    file within an extern "C".
51 */
52 #if defined(__cplusplus)
53 void assert_never_put_petsc_headers_inside_an_extern_c(int); void assert_never_put_petsc_headers_inside_an_extern_c(double);
54 #endif
55 
56 #if defined(__cplusplus)
57 #  define PETSC_RESTRICT PETSC_CXX_RESTRICT
58 #else
59 #  define PETSC_RESTRICT PETSC_C_RESTRICT
60 #endif
61 
62 #if defined(__cplusplus)
63 #  define PETSC_STATIC_INLINE PETSC_CXX_STATIC_INLINE
64 #else
65 #  define PETSC_STATIC_INLINE PETSC_C_STATIC_INLINE
66 #endif
67 
68 #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES) /* For Win32 shared libraries */
69 #  define PETSC_DLLEXPORT __declspec(dllexport)
70 #  define PETSC_DLLIMPORT __declspec(dllimport)
71 #  define PETSC_VISIBILITY_INTERNAL
72 #elif defined(PETSC_USE_VISIBILITY)
73 #  define PETSC_DLLEXPORT __attribute__((visibility ("default")))
74 #  define PETSC_DLLIMPORT __attribute__((visibility ("default")))
75 #  define PETSC_VISIBILITY_INTERNAL __attribute__((visibility ("hidden")))
76 #else
77 #  define PETSC_DLLEXPORT
78 #  define PETSC_DLLIMPORT
79 #  define PETSC_VISIBILITY_INTERNAL
80 #endif
81 
82 #if defined(petsc_EXPORTS)      /* CMake defines this when building the shared library */
83 #  define PETSC_VISIBILITY_PUBLIC PETSC_DLLEXPORT
84 #else  /* Win32 users need this to import symbols from petsc.dll */
85 #  define PETSC_VISIBILITY_PUBLIC PETSC_DLLIMPORT
86 #endif
87 
88 /*
89     Functions tagged with PETSC_EXTERN in the header files are
90   always defined as extern "C" when compiled with C++ so they may be
91   used from C and are always visible in the shared libraries
92 */
93 #if defined(__cplusplus)
94 #define PETSC_EXTERN extern "C" PETSC_VISIBILITY_PUBLIC
95 #define PETSC_EXTERN_TYPEDEF extern "C"
96 #define PETSC_INTERN extern "C" PETSC_VISIBILITY_INTERNAL
97 #else
98 #define PETSC_EXTERN extern PETSC_VISIBILITY_PUBLIC
99 #define PETSC_EXTERN_TYPEDEF
100 #define PETSC_INTERN extern PETSC_VISIBILITY_INTERNAL
101 #endif
102 
103 #include <petscversion.h>
104 #define PETSC_AUTHOR_INFO  "       The PETSc Team\n    petsc-maint@mcs.anl.gov\n http://www.mcs.anl.gov/petsc/\n"
105 
106 /* ========================================================================== */
107 
108 /*
109     Defines the interface to MPI allowing the use of all MPI functions.
110 
111     PETSc does not use the C++ binding of MPI at ALL. The following flag
112     makes sure the C++ bindings are not included. The C++ bindings REQUIRE
113     putting mpi.h before ANY C++ include files, we cannot control this
114     with all PETSc users. Users who want to use the MPI C++ bindings can include
115     mpicxx.h directly in their code
116 */
117 #if !defined(MPICH_SKIP_MPICXX)
118 #  define MPICH_SKIP_MPICXX 1
119 #endif
120 #if !defined(OMPI_SKIP_MPICXX)
121 #  define OMPI_SKIP_MPICXX 1
122 #endif
123 #if !defined(OMPI_WANT_MPI_INTERFACE_WARNING)
124 #  define OMPI_WANT_MPI_INTERFACE_WARNING 0
125 #endif
126 #include <mpi.h>
127 
128 /*
129    Perform various sanity checks that the correct mpi.h is being included at compile time.
130    This usually happens because
131       * either an unexpected mpi.h is in the default compiler path (i.e. in /usr/include) or
132       * an extra include path -I/something (which contains the unexpected mpi.h) is being passed to the compiler
133 */
134 #if defined(PETSC_HAVE_MPIUNI)
135 #  if !defined(__MPIUNI_H)
136 #    error "PETSc was configured with --with-mpi=0 but now appears to be compiling using a different mpi.h"
137 #  endif
138 #elif defined(PETSC_HAVE_MPICH_NUMVERSION)
139 #  if !defined(MPICH_NUMVERSION)
140 #    error "PETSc was configured with MPICH but now appears to be compiling using a non-MPICH mpi.h"
141 #  elif MPICH_NUMVERSION != PETSC_HAVE_MPICH_NUMVERSION
142 #    error "PETSc was configured with one MPICH mpi.h version but now appears to be compiling using a different MPICH mpi.h version"
143 #  endif
144 #elif defined(PETSC_HAVE_OMPI_MAJOR_VERSION)
145 #  if !defined(OMPI_MAJOR_VERSION)
146 #    error "PETSc was configured with OpenMPI but now appears to be compiling using a non-OpenMPI mpi.h"
147 #  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)
148 #    error "PETSc was configured with one OpenMPI mpi.h version but now appears to be compiling using a different OpenMPI mpi.h version"
149 #  endif
150 #endif
151 
152 /*
153     Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler
154     see the top of mpicxx.h in the MPICH2 distribution.
155 */
156 #include <stdio.h>
157 
158 /* MSMPI on 32bit windows requires this yukky hack - that breaks MPI standard compliance */
159 #if !defined(MPIAPI)
160 #define MPIAPI
161 #endif
162 
163 /*
164    Support for Clang (>=3.2) matching type tag arguments with void* buffer types.
165    This allows the compiler to detect cases where the MPI datatype argument passed to a MPI routine
166    does not match the actual type of the argument being passed in
167 */
168 #if defined(__has_attribute) && defined(works_with_const_which_is_not_true)
169 #  if __has_attribute(argument_with_type_tag) && __has_attribute(pointer_with_type_tag) && __has_attribute(type_tag_for_datatype)
170 #    define PetscAttrMPIPointerWithType(bufno,typeno) __attribute__((pointer_with_type_tag(MPI,bufno,typeno)))
171 #    define PetscAttrMPITypeTag(type)                 __attribute__((type_tag_for_datatype(MPI,type)))
172 #    define PetscAttrMPITypeTagLayoutCompatible(type) __attribute__((type_tag_for_datatype(MPI,type,layout_compatible)))
173 #  endif
174 #endif
175 #if !defined(PetscAttrMPIPointerWithType)
176 #  define PetscAttrMPIPointerWithType(bufno,typeno)
177 #  define PetscAttrMPITypeTag(type)
178 #  define PetscAttrMPITypeTagLayoutCompatible(type)
179 #endif
180 
181 /*MC
182     PetscErrorCode - datatype used for return error code from almost all PETSc functions
183 
184     Level: beginner
185 
186 .seealso: CHKERRQ, SETERRQ
187 M*/
188 typedef int PetscErrorCode;
189 
190 /*MC
191 
192     PetscClassId - A unique id used to identify each PETSc class.
193 
194     Notes: Use PetscClassIdRegister() to obtain a new value for a new class being created. Usually
195          XXXInitializePackage() calls it for each class it defines.
196 
197     Developer Notes: Internal integer stored in the _p_PetscObject data structure.
198          These are all computed by an offset from the lowest one, PETSC_SMALLEST_CLASSID.
199 
200     Level: developer
201 
202 .seealso: PetscClassIdRegister(), PetscLogEventRegister(), PetscHeaderCreate()
203 M*/
204 typedef int PetscClassId;
205 
206 
207 /*MC
208     PetscMPIInt - datatype used to represent 'int' parameters to MPI functions.
209 
210     Level: intermediate
211 
212     Notes: usually this is the same as PetscInt, but if PETSc was built with --with-64-bit-indices but
213            standard C/Fortran integers are 32 bit then this is NOT the same as PetscInt it remains 32 bit
214 
215     PetscMPIIntCast(a,&b) checks if the given PetscInt a will fit in a PetscMPIInt, if not it
216       generates a PETSC_ERR_ARG_OUTOFRANGE error.
217 
218 .seealso: PetscBLASInt, PetscInt, PetscMPIIntCast()
219 
220 M*/
221 typedef int PetscMPIInt;
222 
223 /*MC
224     PetscEnum - datatype used to pass enum types within PETSc functions.
225 
226     Level: intermediate
227 
228 .seealso: PetscOptionsGetEnum(), PetscOptionsEnum(), PetscBagRegisterEnum()
229 M*/
230 typedef enum { ENUM_DUMMY } PetscEnum;
231 PETSC_EXTERN MPI_Datatype MPIU_ENUM PetscAttrMPITypeTag(PetscEnum);
232 
233 #if defined(PETSC_HAVE_STDINT_H)
234 #include <stdint.h>
235 #endif
236 
237 /*MC
238     PetscInt - PETSc type that represents integer - used primarily to
239       represent size of arrays and indexing into arrays. Its size can be configured with the option
240       --with-64-bit-indices - to be either 32bit or 64bit [default 32 bit ints]
241 
242    Level: intermediate
243 
244 .seealso: PetscScalar, PetscBLASInt, PetscMPIInt
245 M*/
246 #if defined(PETSC_HAVE_STDINT_H)
247 typedef int64_t Petsc64bitInt;
248 #elif (PETSC_SIZEOF_LONG_LONG == 8)
249 typedef long long Petsc64bitInt;
250 #elif defined(PETSC_HAVE___INT64)
251 typedef __int64 Petsc64bitInt;
252 #else
253 typedef unknown64bit Petsc64bitInt
254 #endif
255 #if defined(PETSC_USE_64BIT_INDICES)
256 typedef Petsc64bitInt PetscInt;
257 #  if defined(PETSC_HAVE_MPI_INT64_T) /* MPI_INT64_T is not guaranteed to be a macro */
258 #    define MPIU_INT MPI_LONG_LONG_INT
259 #  else
260 #    define MPIU_INT MPI_LONG_LONG_INT
261 #  endif
262 #else
263 typedef int PetscInt;
264 #define MPIU_INT MPI_INT
265 #endif
266 #if defined(PETSC_HAVE_MPI_INT64_T)
267 #  define MPIU_INT64 MPI_INT64_T
268 #else
269 #  define MPIU_INT64 MPI_LONG_LONG_INT
270 #endif
271 
272 
273 /*MC
274     PetscBLASInt - datatype used to represent 'int' parameters to BLAS/LAPACK functions.
275 
276     Level: intermediate
277 
278     Notes: usually this is the same as PetscInt, but if PETSc was built with --with-64-bit-indices but
279            standard C/Fortran integers are 32 bit then this is NOT the same as PetscInt it remains 32 bit
280            (except on very rare BLAS/LAPACK implementations that support 64 bit integers see the note below).
281 
282     PetscErrorCode PetscBLASIntCast(a,&b) checks if the given PetscInt a will fit in a PetscBLASInt, if not it
283       generates a PETSC_ERR_ARG_OUTOFRANGE error
284 
285     Installation Notes: The 64bit versions of MATLAB ship with BLAS and LAPACK that use 64 bit integers for sizes etc,
286      if you run ./configure with the option
287      --with-blas-lapack-lib=[/Applications/MATLAB_R2010b.app/bin/maci64/libmwblas.dylib,/Applications/MATLAB_R2010b.app/bin/maci64/libmwlapack.dylib]
288      but you need to also use --known-64-bit-blas-indices.
289 
290         MKL also ships with 64 bit integer versions of the BLAS and LAPACK, if you select those you must also ./configure with --known-64-bit-blas-indices
291 
292      Developer Notes: Eventually ./configure should automatically determine the size of the integers used by BLAS/LAPACK.
293 
294      External packages such as hypre, ML, SuperLU etc do not provide any support for passing 64 bit integers to BLAS/LAPACK so cannot
295      be used with PETSc if you have set PetscBLASInt to long int.
296 
297 .seealso: PetscMPIInt, PetscInt, PetscBLASIntCast()
298 
299 M*/
300 #if defined(PETSC_HAVE_64BIT_BLAS_INDICES)
301 typedef Petsc64bitInt PetscBLASInt;
302 #else
303 typedef int PetscBLASInt;
304 #endif
305 
306 /*EC
307 
308     PetscPrecision - indicates what precision the object is using. This is currently not used.
309 
310     Level: advanced
311 
312 .seealso: PetscObjectSetPrecision()
313 E*/
314 typedef enum { PETSC_PRECISION_SINGLE=4,PETSC_PRECISION_DOUBLE=8 } PetscPrecision;
315 PETSC_EXTERN const char *PetscPrecisions[];
316 
317 /*
318     For the rare cases when one needs to send a size_t object with MPI
319 */
320 #if (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_INT)
321 #define MPIU_SIZE_T MPI_UNSIGNED
322 #elif  (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_LONG)
323 #define MPIU_SIZE_T MPI_UNSIGNED_LONG
324 #elif  (PETSC_SIZEOF_SIZE_T) == (PETSC_SIZEOF_LONG_LONG)
325 #define MPIU_SIZE_T MPI_UNSIGNED_LONG_LONG
326 #else
327 #error "Unknown size for size_t! Send us a bugreport at petsc-maint@mcs.anl.gov"
328 #endif
329 
330 /*
331       You can use PETSC_STDOUT as a replacement of stdout. You can also change
332     the value of PETSC_STDOUT to redirect all standard output elsewhere
333 */
334 PETSC_EXTERN FILE* PETSC_STDOUT;
335 
336 /*
337       You can use PETSC_STDERR as a replacement of stderr. You can also change
338     the value of PETSC_STDERR to redirect all standard error elsewhere
339 */
340 PETSC_EXTERN FILE* PETSC_STDERR;
341 
342 /*MC
343     PetscUnlikely - hints the compiler that the given condition is usually FALSE
344 
345     Synopsis:
346     #include <petscsys.h>
347     PetscBool  PetscUnlikely(PetscBool  cond)
348 
349     Not Collective
350 
351     Input Parameters:
352 .   cond - condition or expression
353 
354     Note: This returns the same truth value, it is only a hint to compilers that the resulting
355     branch is unlikely.
356 
357     Level: advanced
358 
359 .seealso: PetscLikely(), CHKERRQ
360 M*/
361 
362 /*MC
363     PetscLikely - hints the compiler that the given condition is usually TRUE
364 
365     Synopsis:
366     #include <petscsys.h>
367     PetscBool  PetscUnlikely(PetscBool  cond)
368 
369     Not Collective
370 
371     Input Parameters:
372 .   cond - condition or expression
373 
374     Note: This returns the same truth value, it is only a hint to compilers that the resulting
375     branch is likely.
376 
377     Level: advanced
378 
379 .seealso: PetscUnlikely()
380 M*/
381 #if defined(PETSC_HAVE_BUILTIN_EXPECT)
382 #  define PetscUnlikely(cond)   __builtin_expect(!!(cond),0)
383 #  define PetscLikely(cond)     __builtin_expect(!!(cond),1)
384 #else
385 #  define PetscUnlikely(cond)   (cond)
386 #  define PetscLikely(cond)     (cond)
387 #endif
388 
389 /*
390     Declare extern C stuff after including external header files
391 */
392 
393 
394 /*
395        Basic PETSc constants
396 */
397 
398 /*E
399     PetscBool  - Logical variable. Actually an int in C and a logical in Fortran.
400 
401    Level: beginner
402 
403    Developer Note: Why have PetscBool , why not use bool in C? The problem is that K and R C, C99 and C++ all have different mechanisms for
404       boolean values. It is not easy to have a simple macro that that will work properly in all circumstances with all three mechanisms.
405 
406 .seealso: PETSC_TRUE, PETSC_FALSE, PetscNot()
407 E*/
408 typedef enum { PETSC_FALSE,PETSC_TRUE } PetscBool;
409 PETSC_EXTERN const char *const PetscBools[];
410 PETSC_EXTERN MPI_Datatype MPIU_BOOL PetscAttrMPITypeTag(PetscBool);
411 
412 /*
413     Defines some elementary mathematics functions and constants.
414 */
415 #include <petscmath.h>
416 
417 /*E
418     PetscCopyMode  - Determines how an array passed to certain functions is copied or retained
419 
420    Level: beginner
421 
422 $   PETSC_COPY_VALUES - the array values are copied into new space, the user is free to reuse or delete the passed in array
423 $   PETSC_OWN_POINTER - the array values are NOT copied, the object takes ownership of the array and will free it later, the user cannot change or
424 $                       delete the array. The array MUST have been obtained with PetscMalloc(). Hence this mode cannot be used in Fortran.
425 $   PETSC_USE_POINTER - the array values are NOT copied, the object uses the array but does NOT take ownership of the array. The user cannot use
426                         the array but the user must delete the array after the object is destroyed.
427 
428 E*/
429 typedef enum { PETSC_COPY_VALUES, PETSC_OWN_POINTER, PETSC_USE_POINTER} PetscCopyMode;
430 PETSC_EXTERN const char *const PetscCopyModes[];
431 
432 /*MC
433     PETSC_FALSE - False value of PetscBool
434 
435     Level: beginner
436 
437     Note: Zero integer
438 
439 .seealso: PetscBool, PETSC_TRUE
440 M*/
441 
442 /*MC
443     PETSC_TRUE - True value of PetscBool
444 
445     Level: beginner
446 
447     Note: Nonzero integer
448 
449 .seealso: PetscBool, PETSC_FALSE
450 M*/
451 
452 /*MC
453     PETSC_NULL - standard way of passing in a null or array or pointer. This is deprecated in C/C++ simply use NULL
454 
455    Level: beginner
456 
457    Notes: accepted by many PETSc functions to not set a parameter and instead use
458           some default
459 
460           This macro does not exist in Fortran; you must use PETSC_NULL_INTEGER,
461           PETSC_NULL_DOUBLE_PRECISION, PETSC_NULL_FUNCTION, PETSC_NULL_OBJECT etc
462 
463 .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_IGNORE, PETSC_DETERMINE
464 
465 M*/
466 #define PETSC_NULL           NULL
467 
468 /*MC
469     PETSC_IGNORE - same as NULL, means PETSc will ignore this argument
470 
471    Level: beginner
472 
473    Note: accepted by many PETSc functions to not set a parameter and instead use
474           some default
475 
476    Fortran Notes: This macro does not exist in Fortran; you must use PETSC_NULL_INTEGER,
477           PETSC_NULL_DOUBLE_PRECISION etc
478 
479 .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_NULL, PETSC_DETERMINE
480 
481 M*/
482 #define PETSC_IGNORE         NULL
483 
484 /*MC
485     PETSC_DECIDE - standard way of passing in integer or floating point parameter
486        where you wish PETSc to use the default.
487 
488    Level: beginner
489 
490 .seealso: PETSC_NULL, PETSC_DEFAULT, PETSC_IGNORE, PETSC_DETERMINE
491 
492 M*/
493 #define PETSC_DECIDE  -1
494 
495 /*MC
496     PETSC_DETERMINE - standard way of passing in integer or floating point parameter
497        where you wish PETSc to compute the required value.
498 
499    Level: beginner
500 
501    Developer Note: I would like to use const PetscInt PETSC_DETERMINE = PETSC_DECIDE; but for
502      some reason this is not allowed by the standard even though PETSC_DECIDE is a constant value.
503 
504 .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_IGNORE, VecSetSizes()
505 
506 M*/
507 #define PETSC_DETERMINE PETSC_DECIDE
508 
509 /*MC
510     PETSC_DEFAULT - standard way of passing in integer or floating point parameter
511        where you wish PETSc to use the default.
512 
513    Level: beginner
514 
515    Fortran Notes: You need to use PETSC_DEFAULT_INTEGER or PETSC_DEFAULT_REAL.
516 
517 .seealso: PETSC_DECIDE, PETSC_IGNORE, PETSC_DETERMINE
518 
519 M*/
520 #define PETSC_DEFAULT  -2
521 
522 /*MC
523     PETSC_COMM_WORLD - the equivalent of the MPI_COMM_WORLD communicator which represents
524            all the processs that PETSc knows about.
525 
526    Level: beginner
527 
528    Notes: By default PETSC_COMM_WORLD and MPI_COMM_WORLD are identical unless you wish to
529           run PETSc on ONLY a subset of MPI_COMM_WORLD. In that case create your new (smaller)
530           communicator, call it, say comm, and set PETSC_COMM_WORLD = comm BEFORE calling
531           PetscInitialize(), but after MPI_Init() has been called.
532 
533           The value of PETSC_COMM_WORLD should never be USED/accessed before PetscInitialize()
534           is called because it may not have a valid value yet.
535 
536 .seealso: PETSC_COMM_SELF
537 
538 M*/
539 PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;
540 
541 /*MC
542     PETSC_COMM_SELF - This is always MPI_COMM_SELF
543 
544    Level: beginner
545 
546    Notes: Do not USE/access or set this variable before PetscInitialize() has been called.
547 
548 .seealso: PETSC_COMM_WORLD
549 
550 M*/
551 #define PETSC_COMM_SELF MPI_COMM_SELF
552 
553 PETSC_EXTERN PetscBool PetscBeganMPI;
554 PETSC_EXTERN PetscBool PetscInitializeCalled;
555 PETSC_EXTERN PetscBool PetscFinalizeCalled;
556 PETSC_EXTERN PetscBool PetscCUSPSynchronize;
557 PETSC_EXTERN PetscBool PetscViennaCLSynchronize;
558 PETSC_EXTERN PetscBool PetscCUDASynchronize;
559 
560 PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm));
561 PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*);
562 PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm*);
563 
564 /*MC
565    PetscMalloc - Allocates memory, One should use PetscMalloc1() or PetscCalloc1() usually instead of this
566 
567    Synopsis:
568     #include <petscsys.h>
569    PetscErrorCode PetscMalloc(size_t m,void **result)
570 
571    Not Collective
572 
573    Input Parameter:
574 .  m - number of bytes to allocate
575 
576    Output Parameter:
577 .  result - memory allocated
578 
579    Level: beginner
580 
581    Notes:
582    Memory is always allocated at least double aligned
583 
584    It is safe to allocate size 0 and pass the resulting pointer (which may or may not be NULL) to PetscFree().
585 
586 .seealso: PetscFree(), PetscNew()
587 
588   Concepts: memory allocation
589 
590 M*/
591 #define PetscMalloc(a,b)  ((*PetscTrMalloc)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__,(void**)(b)))
592 
593 /*MC
594    PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment
595 
596    Synopsis:
597     #include <petscsys.h>
598    void *PetscAddrAlign(void *addr)
599 
600    Not Collective
601 
602    Input Parameters:
603 .  addr - address to align (any pointer type)
604 
605    Level: developer
606 
607 .seealso: PetscMallocAlign()
608 
609   Concepts: memory allocation
610 M*/
611 #define PetscAddrAlign(a) (void*)((((PETSC_UINTPTR_T)(a))+(PETSC_MEMALIGN-1)) & ~(PETSC_MEMALIGN-1))
612 
613 /*MC
614    PetscMalloc1 - Allocates an array of memory aligned to PETSC_MEMALIGN
615 
616    Synopsis:
617     #include <petscsys.h>
618    PetscErrorCode PetscMalloc1(size_t m1,type **r1)
619 
620    Not Collective
621 
622    Input Parameter:
623 .  m1 - number of elements to allocate in 1st chunk  (may be zero)
624 
625    Output Parameter:
626 .  r1 - memory allocated in first chunk
627 
628    Level: developer
629 
630 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc1(), PetscMalloc2()
631 
632   Concepts: memory allocation
633 
634 M*/
635 #define PetscMalloc1(m1,r1) PetscMalloc((m1)*sizeof(**(r1)),r1)
636 
637 /*MC
638    PetscCalloc1 - Allocates a cleared (zeroed) array of memory aligned to PETSC_MEMALIGN
639 
640    Synopsis:
641     #include <petscsys.h>
642    PetscErrorCode PetscCalloc1(size_t m1,type **r1)
643 
644    Not Collective
645 
646    Input Parameter:
647 .  m1 - number of elements to allocate in 1st chunk  (may be zero)
648 
649    Output Parameter:
650 .  r1 - memory allocated in first chunk
651 
652    Level: developer
653 
654 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc1(), PetscCalloc2()
655 
656   Concepts: memory allocation
657 
658 M*/
659 #define PetscCalloc1(m1,r1) (PetscMalloc1((m1),r1) || PetscMemzero(*(r1),(m1)*sizeof(**(r1))))
660 
661 /*MC
662    PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN
663 
664    Synopsis:
665     #include <petscsys.h>
666    PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)
667 
668    Not Collective
669 
670    Input Parameter:
671 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
672 -  m2 - number of elements to allocate in 2nd chunk  (may be zero)
673 
674    Output Parameter:
675 +  r1 - memory allocated in first chunk
676 -  r2 - memory allocated in second chunk
677 
678    Level: developer
679 
680 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc1(), PetscCalloc2()
681 
682   Concepts: memory allocation
683 
684 M*/
685 #if !defined(PETSC_USE_MALLOC_COALESCED)
686 #define PetscMalloc2(m1,r1,m2,r2) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)))
687 #else
688 #define PetscMalloc2(m1,r1,m2,r2) ((((m1)+(m2)) ? (*(r2) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(PETSC_MEMALIGN-1),r1)) : 0) \
689                                    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),0) \
690                                    || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0))
691 #endif
692 
693 /*MC
694    PetscCalloc2 - Allocates 2 cleared (zeroed) arrays of memory both aligned to PETSC_MEMALIGN
695 
696    Synopsis:
697     #include <petscsys.h>
698    PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)
699 
700    Not Collective
701 
702    Input Parameter:
703 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
704 -  m2 - number of elements to allocate in 2nd chunk  (may be zero)
705 
706    Output Parameter:
707 +  r1 - memory allocated in first chunk
708 -  r2 - memory allocated in second chunk
709 
710    Level: developer
711 
712 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc1(), PetscMalloc2()
713 
714   Concepts: memory allocation
715 M*/
716 #define PetscCalloc2(m1,r1,m2,r2) (PetscMalloc2((m1),(r1),(m2),(r2)) || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))))
717 
718 /*MC
719    PetscMalloc3 - Allocates 3 arrays of memory, all aligned to PETSC_MEMALIGN
720 
721    Synopsis:
722     #include <petscsys.h>
723    PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
724 
725    Not Collective
726 
727    Input Parameter:
728 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
729 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
730 -  m3 - number of elements to allocate in 3rd chunk  (may be zero)
731 
732    Output Parameter:
733 +  r1 - memory allocated in first chunk
734 .  r2 - memory allocated in second chunk
735 -  r3 - memory allocated in third chunk
736 
737    Level: developer
738 
739 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc3(), PetscFree3()
740 
741   Concepts: memory allocation
742 
743 M*/
744 #if !defined(PETSC_USE_MALLOC_COALESCED)
745 #define PetscMalloc3(m1,r1,m2,r2,m3,r3) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)) || PetscMalloc1((m3),(r3)))
746 #else
747 #define PetscMalloc3(m1,r1,m2,r2,m3,r3) ((((m1)+(m2)+(m3)) ? (*(r2) = 0,*(r3) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(m3)*sizeof(**(r3))+2*(PETSC_MEMALIGN-1),r1)) : 0) \
748                                          || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),0) \
749                                          || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0))
750 #endif
751 
752 /*MC
753    PetscCalloc3 - Allocates 3 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN
754 
755    Synopsis:
756     #include <petscsys.h>
757    PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
758 
759    Not Collective
760 
761    Input Parameter:
762 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
763 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
764 -  m3 - number of elements to allocate in 3rd chunk  (may be zero)
765 
766    Output Parameter:
767 +  r1 - memory allocated in first chunk
768 .  r2 - memory allocated in second chunk
769 -  r3 - memory allocated in third chunk
770 
771    Level: developer
772 
773 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc2(), PetscMalloc3(), PetscFree3()
774 
775   Concepts: memory allocation
776 M*/
777 #define PetscCalloc3(m1,r1,m2,r2,m3,r3)                                 \
778   (PetscMalloc3((m1),(r1),(m2),(r2),(m3),(r3))                          \
779    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))))
780 
781 /*MC
782    PetscMalloc4 - Allocates 4 arrays of memory, all aligned to PETSC_MEMALIGN
783 
784    Synopsis:
785     #include <petscsys.h>
786    PetscErrorCode PetscMalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
787 
788    Not Collective
789 
790    Input Parameter:
791 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
792 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
793 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
794 -  m4 - number of elements to allocate in 4th chunk  (may be zero)
795 
796    Output Parameter:
797 +  r1 - memory allocated in first chunk
798 .  r2 - memory allocated in second chunk
799 .  r3 - memory allocated in third chunk
800 -  r4 - memory allocated in fourth chunk
801 
802    Level: developer
803 
804 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc4(), PetscFree4()
805 
806   Concepts: memory allocation
807 
808 M*/
809 #if !defined(PETSC_USE_MALLOC_COALESCED)
810 #define PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)) || PetscMalloc1((m3),(r3)) || PetscMalloc1((m4),(r4)))
811 #else
812 #define PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4)                           \
813   ((((m1)+(m2)+(m3)+(m4)) ? (*(r2) = 0, *(r3) = 0, *(r4) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(m3)*sizeof(**(r3))+(m4)*sizeof(**(r4))+3*(PETSC_MEMALIGN-1),r1)) : 0) \
814    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),0) \
815    || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0) || (!(m4) ? (*(r4) = 0,0) : 0))
816 #endif
817 
818 /*MC
819    PetscCalloc4 - Allocates 4 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN
820 
821    Synopsis:
822     #include <petscsys.h>
823    PetscErrorCode PetscCalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
824 
825    Not Collective
826 
827    Input Parameter:
828 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
829 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
830 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
831 -  m4 - number of elements to allocate in 4th chunk  (may be zero)
832 
833    Output Parameter:
834 +  r1 - memory allocated in first chunk
835 .  r2 - memory allocated in second chunk
836 .  r3 - memory allocated in third chunk
837 -  r4 - memory allocated in fourth chunk
838 
839    Level: developer
840 
841 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc4(), PetscFree4()
842 
843   Concepts: memory allocation
844 
845 M*/
846 #define PetscCalloc4(m1,r1,m2,r2,m3,r3,m4,r4)                           \
847   (PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4)                                \
848    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
849    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))))
850 
851 /*MC
852    PetscMalloc5 - Allocates 5 arrays of memory, all aligned to PETSC_MEMALIGN
853 
854    Synopsis:
855     #include <petscsys.h>
856    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)
857 
858    Not Collective
859 
860    Input Parameter:
861 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
862 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
863 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
864 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
865 -  m5 - number of elements to allocate in 5th chunk  (may be zero)
866 
867    Output Parameter:
868 +  r1 - memory allocated in first chunk
869 .  r2 - memory allocated in second chunk
870 .  r3 - memory allocated in third chunk
871 .  r4 - memory allocated in fourth chunk
872 -  r5 - memory allocated in fifth chunk
873 
874    Level: developer
875 
876 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc5(), PetscFree5()
877 
878   Concepts: memory allocation
879 
880 M*/
881 #if !defined(PETSC_USE_MALLOC_COALESCED)
882 #define PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)) || PetscMalloc1((m3),(r3)) || PetscMalloc1((m4),(r4)) || PetscMalloc1((m5),(r5)))
883 #else
884 #define PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5)      \
885   ((((m1)+(m2)+(m3)+(m4)+(m5)) ? (*(r2) = 0, *(r3) = 0, *(r4) = 0,*(r5) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(m3)*sizeof(**(r3))+(m4)*sizeof(**(r4))+(m5)*sizeof(**(r5))+4*(PETSC_MEMALIGN-1),r1)) : 0) \
886    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),*(void**)(r5) = PetscAddrAlign(*(r4)+(m4)),0) \
887    || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0) || (!(m4) ? (*(r4) = 0,0) : 0) || (!(m5) ? (*(r5) = 0,0) : 0))
888 #endif
889 
890 /*MC
891    PetscCalloc5 - Allocates 5 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN
892 
893    Synopsis:
894     #include <petscsys.h>
895    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)
896 
897    Not Collective
898 
899    Input Parameter:
900 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
901 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
902 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
903 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
904 -  m5 - number of elements to allocate in 5th chunk  (may be zero)
905 
906    Output Parameter:
907 +  r1 - memory allocated in first chunk
908 .  r2 - memory allocated in second chunk
909 .  r3 - memory allocated in third chunk
910 .  r4 - memory allocated in fourth chunk
911 -  r5 - memory allocated in fifth chunk
912 
913    Level: developer
914 
915 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc5(), PetscFree5()
916 
917   Concepts: memory allocation
918 
919 M*/
920 #define PetscCalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5)                     \
921   (PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5)                          \
922    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
923    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))))
924 
925 /*MC
926    PetscMalloc6 - Allocates 6 arrays of memory, all aligned to PETSC_MEMALIGN
927 
928    Synopsis:
929     #include <petscsys.h>
930    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)
931 
932    Not Collective
933 
934    Input Parameter:
935 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
936 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
937 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
938 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
939 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
940 -  m6 - number of elements to allocate in 6th chunk  (may be zero)
941 
942    Output Parameter:
943 +  r1 - memory allocated in first chunk
944 .  r2 - memory allocated in second chunk
945 .  r3 - memory allocated in third chunk
946 .  r4 - memory allocated in fourth chunk
947 .  r5 - memory allocated in fifth chunk
948 -  r6 - memory allocated in sixth chunk
949 
950    Level: developer
951 
952 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc6(), PetscFree3(), PetscFree4(), PetscFree5(), PetscFree6()
953 
954   Concepts: memory allocation
955 
956 M*/
957 #if !defined(PETSC_USE_MALLOC_COALESCED)
958 #define PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)) || PetscMalloc1((m3),(r3)) || PetscMalloc1((m4),(r4)) || PetscMalloc1((m5),(r5)) || PetscMalloc1((m6),(r6)))
959 #else
960 #define PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) \
961   ((((m1)+(m2)+(m3)+(m4)+(m5)+(m6)) ? (*(r2) = 0, *(r3) = 0, *(r4) = 0,*(r5) = 0,*(r6) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(m3)*sizeof(**(r3))+(m4)*sizeof(**(r4))+(m5)*sizeof(**(r5))+(m6)*sizeof(**(r6))+5*(PETSC_MEMALIGN-1),r1)) : 0) \
962    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),*(void**)(r5) = PetscAddrAlign(*(r4)+(m4)),*(void**)(r6) = PetscAddrAlign(*(r5)+(m5)),0) \
963    || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0) || (!(m4) ? (*(r4) = 0,0) : 0) || (!(m5) ? (*(r5) = 0,0) : 0) || (!(m6) ? (*(r6) = 0,0) : 0))
964 #endif
965 
966 /*MC
967    PetscCalloc6 - Allocates 6 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN
968 
969    Synopsis:
970     #include <petscsys.h>
971    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)
972 
973    Not Collective
974 
975    Input Parameter:
976 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
977 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
978 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
979 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
980 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
981 -  m6 - number of elements to allocate in 6th chunk  (may be zero)
982 
983    Output Parameter:
984 +  r1 - memory allocated in first chunk
985 .  r2 - memory allocated in second chunk
986 .  r3 - memory allocated in third chunk
987 .  r4 - memory allocated in fourth chunk
988 .  r5 - memory allocated in fifth chunk
989 -  r6 - memory allocated in sixth chunk
990 
991    Level: developer
992 
993 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscMalloc6(), PetscFree6()
994 
995   Concepts: memory allocation
996 M*/
997 #define PetscCalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6)               \
998   (PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6)                    \
999    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
1000    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))) || PetscMemzero(*(r6),(m6)*sizeof(**(r6))))
1001 
1002 /*MC
1003    PetscMalloc7 - Allocates 7 arrays of memory, all aligned to PETSC_MEMALIGN
1004 
1005    Synopsis:
1006     #include <petscsys.h>
1007    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)
1008 
1009    Not Collective
1010 
1011    Input Parameter:
1012 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
1013 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
1014 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
1015 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
1016 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
1017 .  m6 - number of elements to allocate in 6th chunk  (may be zero)
1018 -  m7 - number of elements to allocate in 7th chunk  (may be zero)
1019 
1020    Output Parameter:
1021 +  r1 - memory allocated in first chunk
1022 .  r2 - memory allocated in second chunk
1023 .  r3 - memory allocated in third chunk
1024 .  r4 - memory allocated in fourth chunk
1025 .  r5 - memory allocated in fifth chunk
1026 .  r6 - memory allocated in sixth chunk
1027 -  r7 - memory allocated in seventh chunk
1028 
1029    Level: developer
1030 
1031 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc7(), PetscFree7()
1032 
1033   Concepts: memory allocation
1034 
1035 M*/
1036 #if !defined(PETSC_USE_MALLOC_COALESCED)
1037 #define PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) (PetscMalloc1((m1),(r1)) || PetscMalloc1((m2),(r2)) || PetscMalloc1((m3),(r3)) || PetscMalloc1((m4),(r4)) || PetscMalloc1((m5),(r5)) || PetscMalloc1((m6),(r6)) || PetscMalloc1((m7),(r7)))
1038 #else
1039 #define PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) \
1040   ((((m1)+(m2)+(m3)+(m4)+(m5)+(m6)+(m7)) ? (*(r2) = 0, *(r3) = 0, *(r4) = 0,*(r5) = 0,*(r6) = 0,*(r7) = 0,PetscMalloc((m1)*sizeof(**(r1))+(m2)*sizeof(**(r2))+(m3)*sizeof(**(r3))+(m4)*sizeof(**(r4))+(m5)*sizeof(**(r5))+(m6)*sizeof(**(r6))+(m7)*sizeof(**(r7))+6*(PETSC_MEMALIGN-1),r1)) : 0) \
1041    || (*(void**)(r2) = PetscAddrAlign(*(r1)+(m1)),*(void**)(r3) = PetscAddrAlign(*(r2)+(m2)),*(void**)(r4) = PetscAddrAlign(*(r3)+(m3)),*(void**)(r5) = PetscAddrAlign(*(r4)+(m4)),*(void**)(r6) = PetscAddrAlign(*(r5)+(m5)),*(void**)(r7) = PetscAddrAlign(*(r6)+(m6)),0) \
1042    || (!(m1) ? (*(r1) = 0,0) : 0) || (!(m2) ? (*(r2) = 0,0) : 0) || (!(m3) ? (*(r3) = 0,0) : 0) || (!(m4) ? (*(r4) = 0,0) : 0) || (!(m5) ? (*(r5) = 0,0) : 0) || (!(m6) ? (*(r6) = 0,0) : 0) || (!(m7) ? (*(r7) = 0,0) : 0))
1043 #endif
1044 
1045 /*MC
1046    PetscCalloc7 - Allocates 7 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN
1047 
1048    Synopsis:
1049     #include <petscsys.h>
1050    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)
1051 
1052    Not Collective
1053 
1054    Input Parameter:
1055 +  m1 - number of elements to allocate in 1st chunk  (may be zero)
1056 .  m2 - number of elements to allocate in 2nd chunk  (may be zero)
1057 .  m3 - number of elements to allocate in 3rd chunk  (may be zero)
1058 .  m4 - number of elements to allocate in 4th chunk  (may be zero)
1059 .  m5 - number of elements to allocate in 5th chunk  (may be zero)
1060 .  m6 - number of elements to allocate in 6th chunk  (may be zero)
1061 -  m7 - number of elements to allocate in 7th chunk  (may be zero)
1062 
1063    Output Parameter:
1064 +  r1 - memory allocated in first chunk
1065 .  r2 - memory allocated in second chunk
1066 .  r3 - memory allocated in third chunk
1067 .  r4 - memory allocated in fourth chunk
1068 .  r5 - memory allocated in fifth chunk
1069 .  r6 - memory allocated in sixth chunk
1070 -  r7 - memory allocated in seventh chunk
1071 
1072    Level: developer
1073 
1074 .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscMalloc7(), PetscFree7()
1075 
1076   Concepts: memory allocation
1077 M*/
1078 #define PetscCalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7)         \
1079   (PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7)              \
1080    || PetscMemzero(*(r1),(m1)*sizeof(**(r1))) || PetscMemzero(*(r2),(m2)*sizeof(**(r2))) || PetscMemzero(*(r3),(m3)*sizeof(**(r3))) \
1081    || PetscMemzero(*(r4),(m4)*sizeof(**(r4))) || PetscMemzero(*(r5),(m5)*sizeof(**(r5))) || PetscMemzero(*(r6),(m6)*sizeof(**(r6))) \
1082    || PetscMemzero(*(r7),(m7)*sizeof(**(r7))))
1083 
1084 /*MC
1085    PetscNew - Allocates memory of a particular type, zeros the memory! Aligned to PETSC_MEMALIGN
1086 
1087    Synopsis:
1088     #include <petscsys.h>
1089    PetscErrorCode PetscNew(type **result)
1090 
1091    Not Collective
1092 
1093    Output Parameter:
1094 .  result - memory allocated, sized to match pointer type
1095 
1096    Level: beginner
1097 
1098 .seealso: PetscFree(), PetscMalloc(), PetscNewLog()
1099 
1100   Concepts: memory allocation
1101 
1102 M*/
1103 #define PetscNew(b)      PetscCalloc1(1,(b))
1104 
1105 /*MC
1106    PetscNewLog - Allocates memory of a type matching pointer, zeros the memory! Aligned to PETSC_MEMALIGN. Associates the memory allocated
1107          with the given object using PetscLogObjectMemory().
1108 
1109    Synopsis:
1110     #include <petscsys.h>
1111    PetscErrorCode PetscNewLog(PetscObject obj,type **result)
1112 
1113    Not Collective
1114 
1115    Input Parameter:
1116 .  obj - object memory is logged to
1117 
1118    Output Parameter:
1119 .  result - memory allocated, sized to match pointer type
1120 
1121    Level: developer
1122 
1123 .seealso: PetscFree(), PetscMalloc(), PetscNew(), PetscLogObjectMemory()
1124 
1125   Concepts: memory allocation
1126 
1127 M*/
1128 #define PetscNewLog(o,b) (PetscNew((b)) || PetscLogObjectMemory((PetscObject)o,sizeof(**(b))))
1129 
1130 /*MC
1131    PetscFree - Frees memory
1132 
1133    Synopsis:
1134     #include <petscsys.h>
1135    PetscErrorCode PetscFree(void *memory)
1136 
1137    Not Collective
1138 
1139    Input Parameter:
1140 .   memory - memory to free (the pointer is ALWAYS set to 0 upon sucess)
1141 
1142    Level: beginner
1143 
1144    Notes:
1145    Memory must have been obtained with PetscNew() or PetscMalloc().
1146    It is safe to call PetscFree() on a NULL pointer.
1147 
1148 .seealso: PetscNew(), PetscMalloc(), PetscFreeVoid()
1149 
1150   Concepts: memory allocation
1151 
1152 M*/
1153 #define PetscFree(a)   ((*PetscTrFree)((void*)(a),__LINE__,PETSC_FUNCTION_NAME,__FILE__) || ((a) = 0,0))
1154 
1155 /*MC
1156    PetscFreeVoid - Frees memory
1157 
1158    Synopsis:
1159     #include <petscsys.h>
1160    void PetscFreeVoid(void *memory)
1161 
1162    Not Collective
1163 
1164    Input Parameter:
1165 .   memory - memory to free
1166 
1167    Level: beginner
1168 
1169    Notes: This is different from PetscFree() in that no error code is returned
1170 
1171 .seealso: PetscFree(), PetscNew(), PetscMalloc()
1172 
1173   Concepts: memory allocation
1174 
1175 M*/
1176 #define PetscFreeVoid(a) ((*PetscTrFree)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__),(a) = 0)
1177 
1178 
1179 /*MC
1180    PetscFree2 - Frees 2 chunks of memory obtained with PetscMalloc2()
1181 
1182    Synopsis:
1183     #include <petscsys.h>
1184    PetscErrorCode PetscFree2(void *memory1,void *memory2)
1185 
1186    Not Collective
1187 
1188    Input Parameter:
1189 +   memory1 - memory to free
1190 -   memory2 - 2nd memory to free
1191 
1192    Level: developer
1193 
1194    Notes: Memory must have been obtained with PetscMalloc2()
1195 
1196 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree()
1197 
1198   Concepts: memory allocation
1199 
1200 M*/
1201 #if !defined(PETSC_USE_MALLOC_COALESCED)
1202 #define PetscFree2(m1,m2)   (PetscFree(m2) || PetscFree(m1))
1203 #else
1204 #define PetscFree2(m1,m2)   ((m1) ? ((m2)=0,PetscFree(m1)) : ((m1)=0,PetscFree(m2)))
1205 #endif
1206 
1207 /*MC
1208    PetscFree3 - Frees 3 chunks of memory obtained with PetscMalloc3()
1209 
1210    Synopsis:
1211     #include <petscsys.h>
1212    PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)
1213 
1214    Not Collective
1215 
1216    Input Parameter:
1217 +   memory1 - memory to free
1218 .   memory2 - 2nd memory to free
1219 -   memory3 - 3rd memory to free
1220 
1221    Level: developer
1222 
1223    Notes: Memory must have been obtained with PetscMalloc3()
1224 
1225 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3()
1226 
1227   Concepts: memory allocation
1228 
1229 M*/
1230 #if !defined(PETSC_USE_MALLOC_COALESCED)
1231 #define PetscFree3(m1,m2,m3)   (PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1232 #else
1233 #define PetscFree3(m1,m2,m3)   ((m1) ? ((m3)=0,(m2)=0,PetscFree(m1)) : ((m2) ? ((m3)=0,(m1)=0,PetscFree(m2)) : ((m2)=0,(m1)=0,PetscFree(m3))))
1234 #endif
1235 
1236 /*MC
1237    PetscFree4 - Frees 4 chunks of memory obtained with PetscMalloc4()
1238 
1239    Synopsis:
1240     #include <petscsys.h>
1241    PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)
1242 
1243    Not Collective
1244 
1245    Input Parameter:
1246 +   m1 - memory to free
1247 .   m2 - 2nd memory to free
1248 .   m3 - 3rd memory to free
1249 -   m4 - 4th memory to free
1250 
1251    Level: developer
1252 
1253    Notes: Memory must have been obtained with PetscMalloc4()
1254 
1255 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4()
1256 
1257   Concepts: memory allocation
1258 
1259 M*/
1260 #if !defined(PETSC_USE_MALLOC_COALESCED)
1261 #define PetscFree4(m1,m2,m3,m4)   (PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1262 #else
1263 #define PetscFree4(m1,m2,m3,m4)   ((m1) ? ((m4)=0,(m3)=0,(m2)=0,PetscFree(m1)) : ((m2) ? ((m4)=0,(m3)=0,(m1)=0,PetscFree(m2)) : ((m3) ? ((m4)=0,(m2)=0,(m1)=0,PetscFree(m3)) : ((m3)=0,(m2)=0,(m1)=0,PetscFree(m4)))))
1264 #endif
1265 
1266 /*MC
1267    PetscFree5 - Frees 5 chunks of memory obtained with PetscMalloc5()
1268 
1269    Synopsis:
1270     #include <petscsys.h>
1271    PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)
1272 
1273    Not Collective
1274 
1275    Input Parameter:
1276 +   m1 - memory to free
1277 .   m2 - 2nd memory to free
1278 .   m3 - 3rd memory to free
1279 .   m4 - 4th memory to free
1280 -   m5 - 5th memory to free
1281 
1282    Level: developer
1283 
1284    Notes: Memory must have been obtained with PetscMalloc5()
1285 
1286 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5()
1287 
1288   Concepts: memory allocation
1289 
1290 M*/
1291 #if !defined(PETSC_USE_MALLOC_COALESCED)
1292 #define PetscFree5(m1,m2,m3,m4,m5)   (PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1293 #else
1294 #define PetscFree5(m1,m2,m3,m4,m5)   ((m1) ? ((m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1)) : ((m2) ? ((m5)=0,(m4)=0,(m3)=0,(m1)=0,PetscFree(m2)) : ((m3) ? ((m5)=0,(m4)=0,(m2)=0,(m1)=0,PetscFree(m3)) : \
1295                                      ((m4) ? ((m5)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m4)) : ((m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m5))))))
1296 #endif
1297 
1298 
1299 /*MC
1300    PetscFree6 - Frees 6 chunks of memory obtained with PetscMalloc6()
1301 
1302    Synopsis:
1303     #include <petscsys.h>
1304    PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)
1305 
1306    Not Collective
1307 
1308    Input Parameter:
1309 +   m1 - memory to free
1310 .   m2 - 2nd memory to free
1311 .   m3 - 3rd memory to free
1312 .   m4 - 4th memory to free
1313 .   m5 - 5th memory to free
1314 -   m6 - 6th memory to free
1315 
1316 
1317    Level: developer
1318 
1319    Notes: Memory must have been obtained with PetscMalloc6()
1320 
1321 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6()
1322 
1323   Concepts: memory allocation
1324 
1325 M*/
1326 #if !defined(PETSC_USE_MALLOC_COALESCED)
1327 #define PetscFree6(m1,m2,m3,m4,m5,m6)   (PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1328 #else
1329 #define PetscFree6(m1,m2,m3,m4,m5,m6)   ((m1) ? ((m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1)) : ((m2) ? ((m6)=0,(m5)=0,(m4)=0,(m3)=0,(m1)=0,PetscFree(m2)) : \
1330                                         ((m3) ? ((m6)=0,(m5)=0,(m4)=0,(m2)=0,(m1)=0,PetscFree(m3)) : ((m4) ? ((m6)=0,(m5)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m4)) : \
1331                                         ((m5) ? ((m6)=0,(m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m5)) : ((m5)=0,(m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m6)))))))
1332 #endif
1333 
1334 /*MC
1335    PetscFree7 - Frees 7 chunks of memory obtained with PetscMalloc7()
1336 
1337    Synopsis:
1338     #include <petscsys.h>
1339    PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)
1340 
1341    Not Collective
1342 
1343    Input Parameter:
1344 +   m1 - memory to free
1345 .   m2 - 2nd memory to free
1346 .   m3 - 3rd memory to free
1347 .   m4 - 4th memory to free
1348 .   m5 - 5th memory to free
1349 .   m6 - 6th memory to free
1350 -   m7 - 7th memory to free
1351 
1352 
1353    Level: developer
1354 
1355    Notes: Memory must have been obtained with PetscMalloc7()
1356 
1357 .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1358           PetscMalloc7()
1359 
1360   Concepts: memory allocation
1361 
1362 M*/
1363 #if !defined(PETSC_USE_MALLOC_COALESCED)
1364 #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   (PetscFree(m7) || PetscFree(m6) || PetscFree(m5) || PetscFree(m4) || PetscFree(m3) || PetscFree(m2) || PetscFree(m1))
1365 #else
1366 #define PetscFree7(m1,m2,m3,m4,m5,m6,m7)   ((m1) ? ((m7)=0,(m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,PetscFree(m1)) : ((m2) ? ((m7)=0,(m6)=0,(m5)=0,(m4)=0,(m3)=0,(m1)=0,PetscFree(m2)) : \
1367                                            ((m3) ? ((m7)=0,(m6)=0,(m5)=0,(m4)=0,(m2)=0,(m1)=0,PetscFree(m3)) : ((m4) ? ((m7)=0,(m6)=0,(m5)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m4)) : \
1368                                            ((m5) ? ((m7)=0,(m6)=0,(m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m5)) : ((m6) ? ((m7)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m6)) : \
1369                                                    ((m6)=0,(m5)=0,(m4)=0,(m3)=0,(m2)=0,(m1)=0,PetscFree(m7))))))))
1370 #endif
1371 
1372 PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,int,const char[],const char[],void**);
1373 PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1374 PETSC_EXTERN PetscErrorCode PetscMallocSet(PetscErrorCode (*)(size_t,int,const char[],const char[],void**),PetscErrorCode (*)(void*,int,const char[],const char[]));
1375 PETSC_EXTERN PetscErrorCode PetscMallocClear(void);
1376 
1377 /*
1378     PetscLogDouble variables are used to contain double precision numbers
1379   that are not used in the numerical computations, but rather in logging,
1380   timing etc.
1381 */
1382 typedef double PetscLogDouble;
1383 #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE
1384 
1385 /*
1386    Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1387 */
1388 PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1389 PETSC_EXTERN PetscErrorCode PetscMallocDumpLog(FILE *);
1390 PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1391 PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1392 PETSC_EXTERN PetscErrorCode PetscMallocDebug(PetscBool);
1393 PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*);
1394 PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1395 PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLog(void);
1396 PETSC_EXTERN PetscErrorCode PetscMallocSetDumpLogThreshold(PetscLogDouble);
1397 PETSC_EXTERN PetscErrorCode PetscMallocGetDumpLog(PetscBool*);
1398 
1399 /*E
1400     PetscDataType - Used for handling different basic data types.
1401 
1402    Level: beginner
1403 
1404    Developer comment: It would be nice if we could always just use MPI Datatypes, why can we not?
1405 
1406 .seealso: PetscBinaryRead(), PetscBinaryWrite(), PetscDataTypeToMPIDataType(),
1407           PetscDataTypeGetSize()
1408 
1409 E*/
1410 typedef enum {PETSC_INT = 0,PETSC_DOUBLE = 1,PETSC_COMPLEX = 2, PETSC_LONG = 3 ,PETSC_SHORT = 4,PETSC_FLOAT = 5,
1411               PETSC_CHAR = 6,PETSC_BIT_LOGICAL = 7,PETSC_ENUM = 8,PETSC_BOOL=9, PETSC___FLOAT128 = 10,PETSC_OBJECT = 11, PETSC_FUNCTION = 12, PETSC_STRING = 12} PetscDataType;
1412 PETSC_EXTERN const char *const PetscDataTypes[];
1413 
1414 #if defined(PETSC_USE_COMPLEX)
1415 #define  PETSC_SCALAR  PETSC_COMPLEX
1416 #else
1417 #if defined(PETSC_USE_REAL_SINGLE)
1418 #define  PETSC_SCALAR  PETSC_FLOAT
1419 #elif defined(PETSC_USE_REAL___FLOAT128)
1420 #define  PETSC_SCALAR  PETSC___FLOAT128
1421 #else
1422 #define  PETSC_SCALAR  PETSC_DOUBLE
1423 #endif
1424 #endif
1425 #if defined(PETSC_USE_REAL_SINGLE)
1426 #define  PETSC_REAL  PETSC_FLOAT
1427 #elif defined(PETSC_USE_REAL___FLOAT128)
1428 #define  PETSC_REAL  PETSC___FLOAT128
1429 #else
1430 #define  PETSC_REAL  PETSC_DOUBLE
1431 #endif
1432 #define  PETSC_FORTRANADDR  PETSC_LONG
1433 
1434 PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1435 PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1436 PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1437 PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);
1438 
1439 /*
1440     Basic memory and string operations. These are usually simple wrappers
1441    around the basic Unix system calls, but a few of them have additional
1442    functionality and/or error checking.
1443 */
1444 PETSC_EXTERN PetscErrorCode PetscBitMemcpy(void*,PetscInt,const void*,PetscInt,PetscInt,PetscDataType);
1445 PETSC_EXTERN PetscErrorCode PetscMemmove(void*,void *,size_t);
1446 PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool  *);
1447 PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1448 PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***);
1449 PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1450 PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool  *);
1451 PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool  *);
1452 PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1453 PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1454 PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1455 PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1456 PETSC_EXTERN PetscErrorCode PetscStrncat(char[],const char[],size_t);
1457 PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1458 PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1459 PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1460 PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1461 PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1462 PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1463 PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1464 PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1465 PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*);
1466 PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1467 PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1468 PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1469 PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1470 PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt,const char *const*,char***);
1471 PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt,char***);
1472 PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);
1473 
1474 PETSC_EXTERN void PetscStrcmpNoError(const char[],const char[],PetscBool  *);
1475 
1476 /*S
1477     PetscToken - 'Token' used for managing tokenizing strings
1478 
1479   Level: intermediate
1480 
1481 .seealso: PetscTokenCreate(), PetscTokenFind(), PetscTokenDestroy()
1482 S*/
1483 typedef struct _p_PetscToken* PetscToken;
1484 
1485 PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1486 PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1487 PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);
1488 
1489 PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1490 PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);
1491 
1492 /*
1493    These are MPI operations for MPI_Allreduce() etc
1494 */
1495 PETSC_EXTERN MPI_Op PetscMaxSum_Op;
1496 #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128)
1497 PETSC_EXTERN MPI_Op MPIU_SUM;
1498 #else
1499 #define MPIU_SUM MPI_SUM
1500 #endif
1501 #if defined(PETSC_USE_REAL___FLOAT128)
1502 PETSC_EXTERN MPI_Op MPIU_MAX;
1503 PETSC_EXTERN MPI_Op MPIU_MIN;
1504 #else
1505 #define MPIU_MAX MPI_MAX
1506 #define MPIU_MIN MPI_MIN
1507 #endif
1508 PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);
1509 
1510 PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1511 PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1512 
1513 /*S
1514      PetscObject - any PETSc object, PetscViewer, Mat, Vec, KSP etc
1515 
1516    Level: beginner
1517 
1518    Note: This is the base class from which all PETSc objects are derived from.
1519 
1520 .seealso:  PetscObjectDestroy(), PetscObjectView(), PetscObjectGetName(), PetscObjectSetName(), PetscObjectReference(), PetscObjectDereference()
1521 S*/
1522 typedef struct _p_PetscObject* PetscObject;
1523 
1524 /*MC
1525     PetscObjectId - unique integer Id for a PetscObject
1526 
1527     Level: developer
1528 
1529     Notes: Unlike pointer values, object ids are never reused.
1530 
1531 .seealso: PetscObjectState, PetscObjectGetId()
1532 M*/
1533 typedef Petsc64bitInt PetscObjectId;
1534 
1535 /*MC
1536     PetscObjectState - integer state for a PetscObject
1537 
1538     Level: developer
1539 
1540     Notes:
1541     Object state is always-increasing and (for objects that track state) can be used to determine if an object has
1542     changed since the last time you interacted with it.  It is 64-bit so that it will not overflow for a very long time.
1543 
1544 .seealso: PetscObjectId, PetscObjectStateGet(), PetscObjectStateIncrease(), PetscObjectStateSet()
1545 M*/
1546 typedef Petsc64bitInt PetscObjectState;
1547 
1548 /*S
1549      PetscFunctionList - Linked list of functions, possibly stored in dynamic libraries, accessed
1550       by string name
1551 
1552    Level: advanced
1553 
1554 .seealso:  PetscFunctionListAdd(), PetscFunctionListDestroy(), PetscOpFlist
1555 S*/
1556 typedef struct _n_PetscFunctionList *PetscFunctionList;
1557 
1558 /*E
1559   PetscFileMode - Access mode for a file.
1560 
1561   Level: beginner
1562 
1563 $  FILE_MODE_READ - open a file at its beginning for reading
1564 $  FILE_MODE_WRITE - open a file at its beginning for writing (will create if the file does not exist)
1565 $  FILE_MODE_APPEND - open a file at end for writing
1566 $  FILE_MODE_UPDATE - open a file for updating, meaning for reading and writing
1567 $  FILE_MODE_APPEND_UPDATE - open a file for updating, meaning for reading and writing, at the end
1568 
1569 .seealso: PetscViewerFileSetMode()
1570 E*/
1571 typedef enum {FILE_MODE_READ, FILE_MODE_WRITE, FILE_MODE_APPEND, FILE_MODE_UPDATE, FILE_MODE_APPEND_UPDATE} PetscFileMode;
1572 extern const char *const PetscFileModes[];
1573 
1574 /*
1575     Defines PETSc error handling.
1576 */
1577 #include <petscerror.h>
1578 
1579 #define PETSC_SMALLEST_CLASSID  1211211
1580 PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1581 PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1582 PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);
1583 
1584 /*
1585    Routines that get memory usage information from the OS
1586 */
1587 PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1588 PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1589 PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1590 PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);
1591 
1592 PETSC_EXTERN PetscErrorCode PetscInfoAllow(PetscBool ,const char []);
1593 PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);
1594 
1595 /*
1596    Initialization of PETSc
1597 */
1598 PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1599 PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1600 PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1601 PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1602 PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1603 PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1604 PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1605 PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1606 PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1607 PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);
1608 
1609 PETSC_EXTERN PetscErrorCode PetscEnd(void);
1610 PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);
1611 
1612 PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1613 PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1614 PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1615 PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);
1616 
1617 /*
1618      These are so that in extern C code we can caste function pointers to non-extern C
1619    function pointers. Since the regular C++ code expects its function pointers to be C++
1620 */
1621 PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1622 PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1623 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);
1624 
1625 /*
1626     Functions that can act on any PETSc object.
1627 */
1628 PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1629 PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1630 PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1631 PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1632 PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1633 PETSC_EXTERN PetscErrorCode PetscObjectSetPrecision(PetscObject,PetscPrecision);
1634 PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1635 PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1636 PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1637 PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1638 PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1639 PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1640 PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1641 PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1642 PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1643 PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt *);
1644 PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1645 PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1646 PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject *);
1647 PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void));
1648 #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d))
1649 PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1650 PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1651 PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt *);
1652 
1653 #include <petscviewertypes.h>
1654 #include <petscoptions.h>
1655 
1656 PETSC_EXTERN PetscErrorCode PetscObjectsListGetGlobalNumbering(MPI_Comm,PetscInt,PetscObject*,PetscInt*,PetscInt*);
1657 
1658 PETSC_EXTERN PetscErrorCode PetscMemoryShowUsage(PetscViewer,const char[]);
1659 PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer,const char[]);
1660 PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1661 PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1662 #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1663 PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1664 PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1665 PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1666 PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1667 PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1668 PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1669 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1670 PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1671 PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,PetscObject,const char[]);
1672 PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1673 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1674 PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1675 PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1676 PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);
1677 
1678 #if defined(PETSC_HAVE_SAWS)
1679 PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1680 PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1681 PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1682 PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1683 PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1684 PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1685 PETSC_EXTERN void           PetscStackSAWsGrantAccess(void);
1686 PETSC_EXTERN void           PetscStackSAWsTakeAccess(void);
1687 PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1688 PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);
1689 
1690 #else
1691 #define PetscSAWsBlock()                        0
1692 #define PetscObjectSAWsViewOff(obj)             0
1693 #define PetscObjectSAWsSetBlock(obj,flg)        0
1694 #define PetscObjectSAWsBlock(obj)               0
1695 #define PetscObjectSAWsGrantAccess(obj)         0
1696 #define PetscObjectSAWsTakeAccess(obj)          0
1697 #define PetscStackViewSAWs()                    0
1698 #define PetscStackSAWsViewOff()                 0
1699 #define PetscStackSAWsTakeAccess()
1700 #define PetscStackSAWsGrantAccess()
1701 
1702 #endif
1703 
1704 typedef void* PetscDLHandle;
1705 typedef enum {PETSC_DL_DECIDE=0,PETSC_DL_NOW=1,PETSC_DL_LOCAL=2} PetscDLMode;
1706 PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1707 PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1708 PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);
1709 
1710 #if defined(PETSC_USE_DEBUG)
1711 PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1712 #endif
1713 PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);
1714 
1715 /*S
1716      PetscObjectList - Linked list of PETSc objects, each accessable by string name
1717 
1718    Level: developer
1719 
1720    Notes: Used by PetscObjectCompose() and PetscObjectQuery()
1721 
1722 .seealso:  PetscObjectListAdd(), PetscObjectListDestroy(), PetscObjectListFind(), PetscObjectCompose(), PetscObjectQuery(), PetscFunctionList
1723 S*/
1724 typedef struct _n_PetscObjectList *PetscObjectList;
1725 
1726 PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1727 PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1728 PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1729 PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1730 PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1731 PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);
1732 
1733 /*
1734     Dynamic library lists. Lists of names of routines in objects or in dynamic
1735   link libraries that will be loaded as needed.
1736 */
1737 
1738 #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1739 PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1740 PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1741 #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1742 PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1743 PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[]);
1744 PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1745 PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1746 PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);
1747 
1748 /*S
1749      PetscDLLibrary - Linked list of dynamics libraries to search for functions
1750 
1751    Level: advanced
1752 
1753 .seealso:  PetscDLLibraryOpen()
1754 S*/
1755 typedef struct _n_PetscDLLibrary *PetscDLLibrary;
1756 PETSC_EXTERN PetscDLLibrary  PetscDLLibrariesLoaded;
1757 PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1758 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1759 PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1760 PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1761 PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool  *);
1762 PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1763 PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);
1764 
1765 /*
1766      Useful utility routines
1767 */
1768 PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1769 PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1770 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1771 PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1772 PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1773 PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);
1774 
1775 /*
1776     PetscNot - negates a logical type value and returns result as a PetscBool
1777 
1778     Notes: This is useful in cases like
1779 $     int        *a;
1780 $     PetscBool  flag = PetscNot(a)
1781      where !a would not return a PetscBool because we cannot provide a cast from int to PetscBool in C.
1782 */
1783 #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)
1784 
1785 /*MC
1786     PetscHelpPrintf - Prints help messages.
1787 
1788    Synopsis:
1789     #include <petscsys.h>
1790      PetscErrorCode (*PetscHelpPrintf)(const char format[],...);
1791 
1792     Not Collective
1793 
1794     Input Parameters:
1795 .   format - the usual printf() format string
1796 
1797    Level: developer
1798 
1799     Fortran Note:
1800     This routine is not supported in Fortran.
1801 
1802     Concepts: help messages^printing
1803     Concepts: printing^help messages
1804 
1805 .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1806 M*/
1807 PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);
1808 
1809 /*
1810      Defines PETSc profiling.
1811 */
1812 #include <petsclog.h>
1813 
1814 /*
1815       Simple PETSc parallel IO for ASCII printing
1816 */
1817 PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1818 PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1819 PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1820 PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1821 PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1822 PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1823 PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);
1824 
1825 /* These are used internally by PETSc ASCII IO routines*/
1826 #include <stdarg.h>
1827 PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list);
1828 PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE*,const char[],va_list);
1829 PETSC_EXTERN PetscErrorCode PetscVFPrintfDefault(FILE*,const char[],va_list);
1830 
1831 #if defined(PETSC_HAVE_MATLAB_ENGINE)
1832 PETSC_EXTERN PetscErrorCode PetscVFPrintf_Matlab(FILE*,const char[],va_list);
1833 #endif
1834 
1835 #if defined(PETSC_HAVE_CLOSURES)
1836 PETSC_EXTERN PetscErrorCode PetscVFPrintfSetClosure(int (^)(const char*));
1837 #endif
1838 
1839 PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1840 PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1841 PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);
1842 
1843 #if defined(PETSC_HAVE_POPEN)
1844 PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1845 PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*,int*);
1846 PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1847 #endif
1848 
1849 PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1850 PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1851 PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*);
1852 PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1853 PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1854 PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1855 PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);
1856 
1857 PETSC_EXTERN PetscErrorCode PetscPopUpSelect(MPI_Comm,const char*,const char*,int,const char**,int*);
1858 
1859 /*S
1860      PetscContainer - Simple PETSc object that contains a pointer to any required data
1861 
1862    Level: advanced
1863 
1864 .seealso:  PetscObject, PetscContainerCreate()
1865 S*/
1866 PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1867 typedef struct _p_PetscContainer*  PetscContainer;
1868 PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void **);
1869 PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void *);
1870 PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1871 PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer *);
1872 PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer, PetscErrorCode (*)(void*));
1873 
1874 /*
1875    For use in debuggers
1876 */
1877 PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1878 PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1879 PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1880 PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1881 PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);
1882 
1883 #include <stddef.h>
1884 #include <string.h>             /* for memcpy, memset */
1885 #if defined(PETSC_HAVE_STDLIB_H)
1886 #include <stdlib.h>
1887 #endif
1888 
1889 #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1890 #include <xmmintrin.h>
1891 #endif
1892 
1893 #undef __FUNCT__
1894 #define __FUNCT__ "PetscMemcpy"
1895 /*@C
1896    PetscMemcpy - Copies n bytes, beginning at location b, to the space
1897    beginning at location a. The two memory regions CANNOT overlap, use
1898    PetscMemmove() in that case.
1899 
1900    Not Collective
1901 
1902    Input Parameters:
1903 +  b - pointer to initial memory space
1904 -  n - length (in bytes) of space to copy
1905 
1906    Output Parameter:
1907 .  a - pointer to copy space
1908 
1909    Level: intermediate
1910 
1911    Compile Option:
1912     PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1913                                   for memory copies on double precision values.
1914     PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1915                                   for memory copies on double precision values.
1916     PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1917                                   for memory copies on double precision values.
1918 
1919    Note:
1920    This routine is analogous to memcpy().
1921 
1922    Developer Note: this is inlined for fastest performance
1923 
1924   Concepts: memory^copying
1925   Concepts: copying^memory
1926 
1927 .seealso: PetscMemmove()
1928 
1929 @*/
1930 PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n)
1931 {
1932 #if defined(PETSC_USE_DEBUG)
1933   size_t al = (size_t) a,bl = (size_t) b;
1934   size_t nl = (size_t) n;
1935   PetscFunctionBegin;
1936   if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1937   if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1938 #else
1939   PetscFunctionBegin;
1940 #endif
1941   if (a != b && n > 0) {
1942 #if defined(PETSC_USE_DEBUG)
1943     if ((al > bl && (al - bl) < nl) || (bl - al) < nl)  SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1944               or make sure your copy regions and lengths are correct. \n\
1945               Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1946 #endif
1947 #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1948    if (!(a % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1949       size_t len = n/sizeof(PetscScalar);
1950 #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1951       PetscBLASInt   one = 1,blen;
1952       PetscErrorCode ierr;
1953       ierr = PetscBLASIntCast(len,&blen);CHKERRQ(ierr);
1954       PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one));
1955 #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1956       fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1957 #else
1958       size_t      i;
1959       PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1960       for (i=0; i<len; i++) y[i] = x[i];
1961 #endif
1962     } else {
1963       memcpy((char*)(a),(char*)(b),n);
1964     }
1965 #else
1966     memcpy((char*)(a),(char*)(b),n);
1967 #endif
1968   }
1969   PetscFunctionReturn(0);
1970 }
1971 
1972 /*@C
1973    PetscMemzero - Zeros the specified memory.
1974 
1975    Not Collective
1976 
1977    Input Parameters:
1978 +  a - pointer to beginning memory location
1979 -  n - length (in bytes) of memory to initialize
1980 
1981    Level: intermediate
1982 
1983    Compile Option:
1984    PETSC_PREFER_BZERO - on certain machines (the IBM RS6000) the bzero() routine happens
1985   to be faster than the memset() routine. This flag causes the bzero() routine to be used.
1986 
1987    Developer Note: this is inlined for fastest performance
1988 
1989    Concepts: memory^zeroing
1990    Concepts: zeroing^memory
1991 
1992 .seealso: PetscMemcpy()
1993 @*/
1994 PETSC_STATIC_INLINE PetscErrorCode  PetscMemzero(void *a,size_t n)
1995 {
1996   if (n > 0) {
1997 #if defined(PETSC_USE_DEBUG)
1998     if (!a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer");
1999 #endif
2000 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
2001     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
2002       size_t      i,len = n/sizeof(PetscScalar);
2003       PetscScalar *x = (PetscScalar*)a;
2004       for (i=0; i<len; i++) x[i] = 0.0;
2005     } else {
2006 #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
2007     if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
2008       PetscInt len = n/sizeof(PetscScalar);
2009       fortranzero_(&len,(PetscScalar*)a);
2010     } else {
2011 #endif
2012 #if defined(PETSC_PREFER_BZERO)
2013       bzero((char *)a,n);
2014 #else
2015       memset((char*)a,0,n);
2016 #endif
2017 #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
2018     }
2019 #endif
2020   }
2021   return 0;
2022 }
2023 
2024 /*MC
2025    PetscPrefetchBlock - Prefetches a block of memory
2026 
2027    Synopsis:
2028     #include <petscsys.h>
2029     void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)
2030 
2031    Not Collective
2032 
2033    Input Parameters:
2034 +  a - pointer to first element to fetch (any type but usually PetscInt or PetscScalar)
2035 .  n - number of elements to fetch
2036 .  rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors)
2037 -  t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note
2038 
2039    Level: developer
2040 
2041    Notes:
2042    The last two arguments (rw and t) must be compile-time constants.
2043 
2044    Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality.  Not all architectures offer
2045    equivalent locality hints, but the following macros are always defined to their closest analogue.
2046 +  PETSC_PREFETCH_HINT_NTA - Non-temporal.  Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched).
2047 .  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.
2048 .  PETSC_PREFETCH_HINT_T1 - Fetch to level 2 and higher (not L1).
2049 -  PETSC_PREFETCH_HINT_T2 - Fetch to high-level cache only.  (On many systems, T0 and T1 are equivalent.)
2050 
2051    This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid
2052    address).
2053 
2054    Concepts: memory
2055 M*/
2056 #define PetscPrefetchBlock(a,n,rw,t) do {                               \
2057     const char *_p = (const char*)(a),*_end = (const char*)((a)+(n));   \
2058     for ( ; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
2059   } while (0)
2060 
2061 /*
2062       Determine if some of the kernel computation routines use
2063    Fortran (rather than C) for the numerical calculations. On some machines
2064    and compilers (like complex numbers) the Fortran version of the routines
2065    is faster than the C/C++ versions. The flag --with-fortran-kernels
2066    should be used with ./configure to turn these on.
2067 */
2068 #if defined(PETSC_USE_FORTRAN_KERNELS)
2069 
2070 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
2071 #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
2072 #endif
2073 
2074 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
2075 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
2076 #endif
2077 
2078 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
2079 #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
2080 #endif
2081 
2082 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
2083 #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
2084 #endif
2085 
2086 #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
2087 #define PETSC_USE_FORTRAN_KERNEL_NORM
2088 #endif
2089 
2090 #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
2091 #define PETSC_USE_FORTRAN_KERNEL_MAXPY
2092 #endif
2093 
2094 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
2095 #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
2096 #endif
2097 
2098 #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
2099 #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
2100 #endif
2101 
2102 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
2103 #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
2104 #endif
2105 
2106 #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
2107 #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
2108 #endif
2109 
2110 #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
2111 #define PETSC_USE_FORTRAN_KERNEL_MDOT
2112 #endif
2113 
2114 #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
2115 #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
2116 #endif
2117 
2118 #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
2119 #define PETSC_USE_FORTRAN_KERNEL_AYPX
2120 #endif
2121 
2122 #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
2123 #define PETSC_USE_FORTRAN_KERNEL_WAXPY
2124 #endif
2125 
2126 #endif
2127 
2128 /*
2129     Macros for indicating code that should be compiled with a C interface,
2130    rather than a C++ interface. Any routines that are dynamically loaded
2131    (such as the PCCreate_XXX() routines) must be wrapped so that the name
2132    mangler does not change the functions symbol name. This just hides the
2133    ugly extern "C" {} wrappers.
2134 */
2135 #if defined(__cplusplus)
2136 #define EXTERN_C_BEGIN extern "C" {
2137 #define EXTERN_C_END }
2138 #else
2139 #define EXTERN_C_BEGIN
2140 #define EXTERN_C_END
2141 #endif
2142 
2143 /* --------------------------------------------------------------------*/
2144 
2145 /*MC
2146     MPI_Comm - the basic object used by MPI to determine which processes are involved in a
2147         communication
2148 
2149    Level: beginner
2150 
2151    Note: This manual page is a place-holder because MPICH does not have a manual page for MPI_Comm
2152 
2153 .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
2154 M*/
2155 
2156 /*MC
2157     PetscScalar - PETSc type that represents either a double precision real number, a double precision
2158        complex number, a single precision real number, a long double or an int - if the code is configured
2159        with --with-scalar-type=real,complex --with-precision=single,double,__float128
2160 
2161    Level: beginner
2162 
2163 .seealso: PetscReal, MPIU_SCALAR, PetscInt, MPIU_REAL
2164 M*/
2165 
2166 /*MC
2167     PetscComplex - PETSc type that represents a complex number with precision matching that of PetscReal.
2168 
2169    Synopsis:
2170    #include <petscsys.h>
2171    PetscComplex number = 1. + 2.*PETSC_i;
2172 
2173    Level: beginner
2174 
2175    Note:
2176    Complex numbers are automatically available if PETSc was able to find a working complex implementation
2177 
2178 .seealso: PetscReal, PetscComplex, MPIU_COMPLEX, PetscInt, PETSC_i
2179 M*/
2180 
2181 /*MC
2182     PetscReal - PETSc type that represents a real number version of PetscScalar
2183 
2184    Level: beginner
2185 
2186 .seealso: PetscScalar
2187 M*/
2188 
2189 /*MC
2190     MPIU_SCALAR - MPI datatype corresponding to PetscScalar
2191 
2192    Level: beginner
2193 
2194     Note: In MPI calls that require an MPI datatype that matches a PetscScalar or array of PetscScalars
2195           pass this value
2196 
2197 .seealso: PetscReal, PetscScalar, MPIU_INT
2198 M*/
2199 
2200 #if defined(PETSC_HAVE_MPIIO)
2201 #if !defined(PETSC_WORDS_BIGENDIAN)
2202 PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2203 PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2204 #else
2205 #define MPIU_File_write_all(a,b,c,d,e) MPI_File_write_all(a,b,c,d,e)
2206 #define MPIU_File_read_all(a,b,c,d,e) MPI_File_read_all(a,b,c,d,e)
2207 #endif
2208 #endif
2209 
2210 /* the following petsc_static_inline require petscerror.h */
2211 
2212 /* Limit MPI to 32-bits */
2213 #define PETSC_MPI_INT_MAX  2147483647
2214 #define PETSC_MPI_INT_MIN -2147483647
2215 /* Limit BLAS to 32-bits */
2216 #define PETSC_BLAS_INT_MAX  2147483647
2217 #define PETSC_BLAS_INT_MIN -2147483647
2218 
2219 #undef __FUNCT__
2220 #define __FUNCT__ "PetscBLASIntCast"
2221 /*@C
2222     PetscBLASIntCast - casts a PetscInt (which may be 64 bits in size) to a PetscBLASInt (which may be 32 bits in size), generates an
2223          error if the PetscBLASInt is not large enough to hold the number.
2224 
2225    Not Collective
2226 
2227    Input Parameter:
2228 .     a - the PetscInt value
2229 
2230    Output Parameter:
2231 .     b - the resulting PetscBLASInt value
2232 
2233    Level: advanced
2234 
2235 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast()
2236 @*/
2237 PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
2238 {
2239   PetscFunctionBegin;
2240   *b =  (PetscBLASInt)(a);
2241 #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
2242   if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK");
2243 #endif
2244   PetscFunctionReturn(0);
2245 }
2246 
2247 #undef __FUNCT__
2248 #define __FUNCT__ "PetscMPIIntCast"
2249 /*@C
2250     PetscMPIIntCast - casts a PetscInt (which may be 64 bits in size) to a PetscMPIInt (which may be 32 bits in size), generates an
2251          error if the PetscMPIInt is not large enough to hold the number.
2252 
2253    Not Collective
2254 
2255    Input Parameter:
2256 .     a - the PetscInt value
2257 
2258    Output Parameter:
2259 .     b - the resulting PetscMPIInt value
2260 
2261    Level: advanced
2262 
2263 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast()
2264 @*/
2265 PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
2266 {
2267   PetscFunctionBegin;
2268   *b =  (PetscMPIInt)(a);
2269 #if defined(PETSC_USE_64BIT_INDICES)
2270   if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for MPI");
2271 #endif
2272   PetscFunctionReturn(0);
2273 }
2274 
2275 #define PetscIntMult64bit(a,b)   ((Petsc64bitInt)(a))*((Petsc64bitInt)(b))
2276 
2277 #undef __FUNCT__
2278 #define __FUNCT__ "PetscRealIntMultTruncate"
2279 /*@C
2280 
2281    PetscRealIntMultTruncate - Computes the product of a positive PetscReal and a positive PetscInt and truncates the value to slightly less than the maximal possible value
2282 
2283    Not Collective
2284 
2285    Input Parameter:
2286 +     a - the PetscReal value
2287 -     b - the second value
2288 
2289    Output Parameter:
2290 .     c - the result as a PetscInt value
2291 
2292    Use PetscIntMult64bit() to compute the product of two PetscInt as a Petsc64bitInt
2293    Use PetscIntMultTruncate() to compute the product of two positive PetscInt and truncate to fit a PetscInt
2294    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
2295 
2296    Developers Note: We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks.
2297 
2298    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2299 
2300    Level: advanced
2301 
2302 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64()
2303 @*/
2304 PETSC_STATIC_INLINE PetscInt PetscRealIntMultTruncate(PetscReal a,PetscInt b)
2305 {
2306   Petsc64bitInt r;
2307 
2308   r  =  (Petsc64bitInt) (a*b);
2309   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2310   return (PetscInt) r;
2311 }
2312 
2313 #undef __FUNCT__
2314 #define __FUNCT__ "PetscIntMultTruncate"
2315 /*@C
2316 
2317    PetscIntMultTruncate - Computes the product of two positive PetscInt and truncates the value to slightly less than the maximal possible value
2318 
2319    Not Collective
2320 
2321    Input Parameter:
2322 +     a - the PetscInt value
2323 -     b - the second value
2324 
2325    Output Parameter:
2326 .     c - the result as a PetscInt value
2327 
2328    Use PetscIntMult64bit() to compute the product of two PetscInt as a Petsc64bitInt
2329    Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2330    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
2331 
2332    Developers Note: We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks.
2333 
2334    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2335 
2336    Level: advanced
2337 
2338 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64()
2339 @*/
2340 PETSC_STATIC_INLINE PetscInt PetscIntMultTruncate(PetscInt a,PetscInt b)
2341 {
2342   Petsc64bitInt r;
2343 
2344   r  =  PetscIntMult64bit(a,b);
2345   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2346   return (PetscInt) r;
2347 }
2348 
2349 #undef __FUNCT__
2350 #define __FUNCT__ "PetscIntSumTruncate"
2351 /*@C
2352 
2353    PetscIntSumTruncate - Computes the sum of two positive PetscInt and truncates the value to slightly less than the maximal possible value
2354 
2355    Not Collective
2356 
2357    Input Parameter:
2358 +     a - the PetscInt value
2359 -     b - the second value
2360 
2361    Output Parameter:
2362 .     c - the result as a PetscInt value
2363 
2364    Use PetscIntMult64bit() to compute the product of two PetscInt as a Petsc64bitInt
2365    Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2366    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
2367 
2368    This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2369 
2370    Level: advanced
2371 
2372 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64()
2373 @*/
2374 PETSC_STATIC_INLINE PetscInt PetscIntSumTruncate(PetscInt a,PetscInt b)
2375 {
2376   Petsc64bitInt r;
2377 
2378   r  =  ((Petsc64bitInt)a) + ((Petsc64bitInt)b);
2379   if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2380   return (PetscInt) r;
2381 }
2382 
2383 #undef __FUNCT__
2384 #define __FUNCT__ "PetscIntMultError"
2385 /*@C
2386 
2387    PetscIntMultError - Computes the product of two positive PetscInt and generates an error with overflow.
2388 
2389    Not Collective
2390 
2391    Input Parameter:
2392 +     a - the PetscInt value
2393 -     b - the second value
2394 
2395    Output Parameter:ma
2396 .     c - the result as a PetscInt value
2397 
2398    Use PetscIntMult64bit() to compute the product of two 32 bit PetscInt and store in a Petsc64bitInt
2399    Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt
2400 
2401    Developers Note: We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks.
2402 
2403    Level: advanced
2404 
2405 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64()
2406 @*/
2407 PETSC_STATIC_INLINE PetscErrorCode PetscIntMultError(PetscInt a,PetscInt b,PetscInt *result)
2408 {
2409   Petsc64bitInt r;
2410 
2411   PetscFunctionBegin;
2412   r  =  PetscIntMult64bit(a,b);
2413 #if !defined(PETSC_USE_64BIT_INDICES)
2414   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);
2415 #endif
2416   if (result) *result = (PetscInt) r;
2417   PetscFunctionReturn(0);
2418 }
2419 
2420  #undef __FUNCT__
2421 #define __FUNCT__ "PetscIntSumError"
2422 /*@C
2423 
2424    PetscIntSumError - Computes the product of two positive PetscInt and generates an error with overflow.
2425 
2426    Not Collective
2427 
2428    Input Parameter:
2429 +     a - the PetscInt value
2430 -     b - the second value
2431 
2432    Output Parameter:ma
2433 .     c - the result as a PetscInt value
2434 
2435    Use PetscIntMult64bit() to compute the product of two 32 bit PetscInt and store in a Petsc64bitInt
2436    Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt
2437 
2438    Level: advanced
2439 
2440 .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64()
2441 @*/
2442 PETSC_STATIC_INLINE PetscErrorCode PetscIntSumError(PetscInt a,PetscInt b,PetscInt *result)
2443 {
2444   Petsc64bitInt r;
2445 
2446   PetscFunctionBegin;
2447   r  =  ((Petsc64bitInt)a) + ((Petsc64bitInt)b);
2448 #if !defined(PETSC_USE_64BIT_INDICES)
2449   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);
2450 #endif
2451   if (result) *result = (PetscInt) r;
2452   PetscFunctionReturn(0);
2453 }
2454 
2455 /*
2456      The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2457 */
2458 #if defined(hz)
2459 #undef hz
2460 #endif
2461 
2462 /*  For arrays that contain filenames or paths */
2463 
2464 
2465 #if defined(PETSC_HAVE_LIMITS_H)
2466 #include <limits.h>
2467 #endif
2468 #if defined(PETSC_HAVE_SYS_PARAM_H)
2469 #include <sys/param.h>
2470 #endif
2471 #if defined(PETSC_HAVE_SYS_TYPES_H)
2472 #include <sys/types.h>
2473 #endif
2474 #if defined(MAXPATHLEN)
2475 #  define PETSC_MAX_PATH_LEN     MAXPATHLEN
2476 #elif defined(MAX_PATH)
2477 #  define PETSC_MAX_PATH_LEN     MAX_PATH
2478 #elif defined(_MAX_PATH)
2479 #  define PETSC_MAX_PATH_LEN     _MAX_PATH
2480 #else
2481 #  define PETSC_MAX_PATH_LEN     4096
2482 #endif
2483 
2484 /*MC
2485 
2486     UsingFortran - Fortran can be used with PETSc in four distinct approaches
2487 
2488 $    1) classic Fortran 77 style
2489 $#include "petsc/finclude/petscXXX.h" to work with material from the XXX component of PETSc
2490 $       XXX variablename
2491 $      You cannot use this approach if you wish to use the Fortran 90 specific PETSc routines
2492 $      which end in F90; such as VecGetArrayF90()
2493 $
2494 $    2) classic Fortran 90 style
2495 $#include "petsc/finclude/petscXXX.h"
2496 $#include "petsc/finclude/petscXXX.h90" to work with material from the XXX component of PETSc
2497 $       XXX variablename
2498 $
2499 $    3) Using Fortran modules
2500 $#include "petsc/finclude/petscXXXdef.h"
2501 $         use petscXXXX
2502 $       XXX variablename
2503 $
2504 $    4) Use Fortran modules and Fortran data types for PETSc types
2505 $#include "petsc/finclude/petscXXXdef.h"
2506 $         use petscXXXX
2507 $       type(XXX) variablename
2508 $      To use this approach you must ./configure PETSc with the additional
2509 $      option --with-fortran-datatypes You cannot use the type(XXX) declaration approach without using Fortran modules
2510 
2511     Finally if you absolutely do not want to use any #include you can use either
2512 
2513 $    3a) skip the #include BUT you cannot use any PETSc data type names like Vec, Mat, PetscInt, PetscErrorCode etc
2514 $        and you must declare the variables as integer, for example
2515 $        integer variablename
2516 $
2517 $    4a) skip the #include, you use the object types like type(Vec) type(Mat) but cannot use the data type
2518 $        names like PetscErrorCode, PetscInt etc. again for those you must use integer
2519 
2520    We recommend either 2 or 3. Approaches 2 and 3 provide type checking for most PETSc function calls; 4 has type checking
2521 for only a few PETSc functions.
2522 
2523    Fortran type checking with interfaces is strick, this means you cannot pass a scalar value when an array value
2524 is expected (even though it is legal Fortran). For example when setting a single value in a matrix with MatSetValues()
2525 you cannot have something like
2526 $      PetscInt row,col
2527 $      PetscScalar val
2528 $        ...
2529 $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)
2530 You must instead have
2531 $      PetscInt row(1),col(1)
2532 $      PetscScalar val(1)
2533 $        ...
2534 $      call MatSetValues(mat,1,row,1,col,val,INSERT_VALUES,ierr)
2535 
2536 
2537     See the example src/vec/vec/examples/tutorials/ex20f90.F90 for an example that can use all four approaches
2538 
2539     Developer Notes: The petsc/finclude/petscXXXdef.h contain all the #defines (would be typedefs in C code) these
2540      automatically include their predecessors; for example petsc/finclude/petscvecdef.h includes petsc/finclude/petscisdef.h
2541 
2542      The petsc/finclude/petscXXXX.h contain all the parameter statements for that package. These automatically include
2543      their petsc/finclude/petscXXXdef.h file but DO NOT automatically include their predecessors;  for example
2544      petsc/finclude/petscvec.h does NOT automatically include petsc/finclude/petscis.h
2545 
2546      The petsc/finclude/ftn-custom/petscXXXdef.h90 are not intended to be used directly in code, they define the
2547      Fortran data type type(XXX) (for example type(Vec)) when PETSc is ./configure with the --with-fortran-datatypes option.
2548 
2549      The petsc/finclude/ftn-custom/petscXXX.h90 (not included directly by code) contain interface definitions for
2550      the PETSc Fortran stubs that have different bindings then their C version (for example VecGetArrayF90).
2551 
2552      The petsc/finclude/ftn-auto/petscXXX.h90 (not included directly by code) contain interface definitions generated
2553      automatically by "make allfortranstubs".
2554 
2555      The petsc/finclude/petscXXX.h90 includes the custom petsc/finclude/ftn-custom/petscXXX.h90 and if ./configure
2556      was run with --with-fortran-interfaces it also includes the petsc/finclude/ftn-auto/petscXXX.h90 These DO NOT automatically
2557      include their predecessors
2558 
2559     Level: beginner
2560 
2561 M*/
2562 
2563 PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2564 PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2565 PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2566 PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2567 PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2568 PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2569 PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2570 
2571 PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2572 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2573 PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2574 PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2575 PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2576 PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2577 PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt*,PetscInt*,PetscInt*);
2578 PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2579 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2580 PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2581 PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2582 PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt,PetscInt[],void*,size_t,void*);
2583 PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2584 PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2585 PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2586 PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2587 PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2588 PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2589 PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt*,const PetscInt*,PetscInt,const PetscInt*,const PetscInt*,PetscInt*,PetscInt**,PetscInt**);
2590 PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt,const PetscInt*,PetscInt,const PetscInt*,PetscInt*,PetscInt**);
2591 PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt,const PetscMPIInt[],PetscInt,const PetscMPIInt[],PetscInt*,PetscMPIInt**);
2592 
2593 PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2594 PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);
2595 
2596 /*J
2597     PetscRandomType - String with the name of a PETSc randomizer
2598 
2599    Level: beginner
2600 
2601    Notes: to use the SPRNG you must have ./configure PETSc
2602    with the option --download-sprng
2603 
2604 .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2605 J*/
2606 typedef const char* PetscRandomType;
2607 #define PETSCRAND       "rand"
2608 #define PETSCRAND48     "rand48"
2609 #define PETSCSPRNG      "sprng"
2610 #define PETSCRANDER48   "rander48"
2611 
2612 /* Logging support */
2613 PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;
2614 
2615 PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);
2616 
2617 /*S
2618      PetscRandom - Abstract PETSc object that manages generating random numbers
2619 
2620    Level: intermediate
2621 
2622   Concepts: random numbers
2623 
2624 .seealso:  PetscRandomCreate(), PetscRandomGetValue(), PetscRandomType
2625 S*/
2626 typedef struct _p_PetscRandom*   PetscRandom;
2627 
2628 /* Dynamic creation and loading functions */
2629 PETSC_EXTERN PetscFunctionList PetscRandomList;
2630 
2631 PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2632 PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2633 PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2634 PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType*);
2635 PETSC_STATIC_INLINE PetscErrorCode PetscRandomViewFromOptions(PetscRandom A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
2636 PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);
2637 
2638 PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2639 PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2640 PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2641 PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2642 PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2643 PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2644 PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2645 PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2646 PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);
2647 
2648 PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2649 PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2650 PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2651 PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2652 PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2653 PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2654 PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);
2655 PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2656 PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);
2657 
2658 PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscDataType);
2659 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscDataType);
2660 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,void*,PetscInt,PetscDataType,PetscBool );
2661 PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,void*,PetscInt,PetscDataType,PetscBool );
2662 PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2663 PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2664 PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool  *);
2665 PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool  *);
2666 PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2667 PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2668 PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2669 PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);
2670 
2671 /*
2672    In binary files variables are stored using the following lengths,
2673   regardless of how they are stored in memory on any one particular
2674   machine. Use these rather then sizeof() in computing sizes for
2675   PetscBinarySeek().
2676 */
2677 #define PETSC_BINARY_INT_SIZE   (32/8)
2678 #define PETSC_BINARY_FLOAT_SIZE  (32/8)
2679 #define PETSC_BINARY_CHAR_SIZE  (8/8)
2680 #define PETSC_BINARY_SHORT_SIZE  (16/8)
2681 #define PETSC_BINARY_DOUBLE_SIZE  (64/8)
2682 #define PETSC_BINARY_SCALAR_SIZE  sizeof(PetscScalar)
2683 
2684 /*E
2685   PetscBinarySeekType - argument to PetscBinarySeek()
2686 
2687   Level: advanced
2688 
2689 .seealso: PetscBinarySeek(), PetscBinarySynchronizedSeek()
2690 E*/
2691 typedef enum {PETSC_BINARY_SEEK_SET = 0,PETSC_BINARY_SEEK_CUR = 1,PETSC_BINARY_SEEK_END = 2} PetscBinarySeekType;
2692 PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2693 PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2694 PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);
2695 
2696 PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2697 PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool );
2698 PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2699 PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2700 PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2701 PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2702 
2703 PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2704 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2705 PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2706 PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2707 PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2708 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt*,const void*,PetscMPIInt*,PetscMPIInt**,void*)
2709   PetscAttrMPIPointerWithType(6,3);
2710 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2711                                                     PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2712                                                     PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2713   PetscAttrMPIPointerWithType(6,3);
2714 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2715                                                        MPI_Request**,MPI_Request**,
2716                                                        PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2717                                                        PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2718   PetscAttrMPIPointerWithType(6,3);
2719 
2720 /*E
2721     PetscBuildTwoSidedType - algorithm for setting up two-sided communication
2722 
2723 $  PETSC_BUILDTWOSIDED_ALLREDUCE - classical algorithm using an MPI_Allreduce with
2724 $      a buffer of length equal to the communicator size. Not memory-scalable due to
2725 $      the large reduction size. Requires only MPI-1.
2726 $  PETSC_BUILDTWOSIDED_IBARRIER - nonblocking algorithm based on MPI_Issend and MPI_Ibarrier.
2727 $      Proved communication-optimal in Hoefler, Siebert, and Lumsdaine (2010). Requires MPI-3.
2728 $  PETSC_BUILDTWOSIDED_REDSCATTER - similar to above, but use more optimized function
2729 $      that only communicates the part of the reduction that is necessary.  Requires MPI-2.
2730 
2731    Level: developer
2732 
2733 .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedSetType(), PetscCommBuildTwoSidedGetType()
2734 E*/
2735 typedef enum {
2736   PETSC_BUILDTWOSIDED_NOTSET = -1,
2737   PETSC_BUILDTWOSIDED_ALLREDUCE = 0,
2738   PETSC_BUILDTWOSIDED_IBARRIER = 1,
2739   PETSC_BUILDTWOSIDED_REDSCATTER = 2
2740   /* Updates here must be accompanied by updates in finclude/petscsys.h and the string array in mpits.c */
2741 } PetscBuildTwoSidedType;
2742 PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2743 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2744 PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);
2745 
2746 PETSC_EXTERN PetscErrorCode PetscSSEIsEnabled(MPI_Comm,PetscBool  *,PetscBool  *);
2747 
2748 /*E
2749   InsertMode - Whether entries are inserted or added into vectors or matrices
2750 
2751   Level: beginner
2752 
2753 .seealso: VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2754           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(),
2755           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd()
2756 E*/
2757  typedef enum {NOT_SET_VALUES, INSERT_VALUES, ADD_VALUES, MAX_VALUES, INSERT_ALL_VALUES, ADD_ALL_VALUES, INSERT_BC_VALUES, ADD_BC_VALUES} InsertMode;
2758 
2759 /*MC
2760     INSERT_VALUES - Put a value into a vector or matrix, overwrites any previous value
2761 
2762     Level: beginner
2763 
2764 .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2765           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), ADD_VALUES,
2766           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES
2767 
2768 M*/
2769 
2770 /*MC
2771     ADD_VALUES - Adds a value into a vector or matrix, if there previously was no value, just puts the
2772                 value into that location
2773 
2774     Level: beginner
2775 
2776 .seealso: InsertMode, VecSetValues(), MatSetValues(), VecSetValue(), VecSetValuesBlocked(),
2777           VecSetValuesLocal(), VecSetValuesBlockedLocal(), MatSetValuesBlocked(), INSERT_VALUES,
2778           MatSetValuesBlockedLocal(), MatSetValuesLocal(), VecScatterBegin(), VecScatterEnd(), MAX_VALUES
2779 
2780 M*/
2781 
2782 /*MC
2783     MAX_VALUES - Puts the maximum of the scattered/gathered value and the current value into each location
2784 
2785     Level: beginner
2786 
2787 .seealso: InsertMode, VecScatterBegin(), VecScatterEnd(), ADD_VALUES, INSERT_VALUES
2788 
2789 M*/
2790 
2791 PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);
2792 
2793 typedef enum {PETSC_SUBCOMM_GENERAL=0,PETSC_SUBCOMM_CONTIGUOUS=1,PETSC_SUBCOMM_INTERLACED=2} PetscSubcommType;
2794 PETSC_EXTERN const char *const PetscSubcommTypes[];
2795 
2796 /*S
2797    PetscSubcomm - A decomposition of an MPI communicator into subcommunicators
2798 
2799    Notes: After a call to PetscSubcommSetType(), PetscSubcommSetTypeGeneral(), or PetscSubcommSetFromOptions() one may call
2800 $     PetscSubcommChild() returns the associated subcommunicator on this process
2801 $     PetscSubcommContiguousParent() returns a parent communitor but with all child of the same subcommunicator having contiquous rank
2802 
2803    Sample Usage:
2804        PetscSubcommCreate()
2805        PetscSubcommSetNumber()
2806        PetscSubcommSetType(PETSC_SUBCOMM_INTERLACED);
2807        ccomm = PetscSubcommChild()
2808        PetscSubcommDestroy()
2809 
2810    Level: advanced
2811 
2812    Concepts: communicator, create
2813 
2814    Notes:
2815 $   PETSC_SUBCOMM_GENERAL - similar to MPI_Comm_split() each process sets the new communicator (color) they will belong to and the order within that communicator
2816 $   PETSC_SUBCOMM_CONTIGUOUS - each new communicator contains a set of process with contiquous ranks in the original MPI communicator
2817 $   PETSC_SUBCOMM_INTERLACED - each new communictor contains a set of processes equally far apart in rank from the others in that new communicator
2818 
2819    Examaple: Consider a communicator with six processes split into 3 subcommunicators.
2820 $     PETSC_SUBCOMM_CONTIGUOUS - the first communicator contains rank 0,1  the second rank 2,3 and the third rank 4,5 in the original ordering of the original communicator
2821 $     PETSC_SUBCOMM_INTERLACED - the first communicator contains rank 0,3, the second 1,4 and the third 2,5
2822 
2823    Developer Notes: This is used in objects such as PCREDUNDANT() to manage the subcommunicators on which the redundant computations
2824       are performed.
2825 
2826 
2827 .seealso: PetscSubcommCreate(), PetscSubcommSetNumber(), PetscSubcommSetType(), PetscSubcommView(), PetscSubcommSetFromOptions()
2828 
2829 S*/
2830 typedef struct _n_PetscSubcomm* PetscSubcomm;
2831 
2832 struct _n_PetscSubcomm {
2833   MPI_Comm         parent;           /* parent communicator */
2834   MPI_Comm         dupparent;        /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2835   MPI_Comm         child;            /* the sub-communicator */
2836   PetscMPIInt      n;                /* num of subcommunicators under the parent communicator */
2837   PetscMPIInt      color;            /* color of processors belong to this communicator */
2838   PetscMPIInt      *subsize;         /* size of subcommunicator[color] */
2839   PetscSubcommType type;
2840 };
2841 
2842 PETSC_STATIC_INLINE MPI_Comm PetscSubcommParent(PetscSubcomm scomm) {return scomm->parent;}
2843 PETSC_STATIC_INLINE MPI_Comm PetscSubcommChild(PetscSubcomm scomm) {return scomm->child;}
2844 PETSC_STATIC_INLINE MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm) {return scomm->dupparent;}
2845 PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2846 PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2847 PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2848 PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2849 PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2850 PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2851 PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2852 
2853 /*S
2854    PetscSegBuffer - a segmented extendable buffer
2855 
2856    Level: developer
2857 
2858 .seealso: PetscSegBufferCreate(), PetscSegBufferGet(), PetscSegBufferExtract(), PetscSegBufferDestroy()
2859 S*/
2860 typedef struct _n_PetscSegBuffer *PetscSegBuffer;
2861 PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2862 PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2863 PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2864 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2865 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2866 PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2867 PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2868 PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);
2869 
2870 /* Type-safe wrapper to encourage use of PETSC_RESTRICT. Does not use PetscFunctionBegin because the error handling
2871  * prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as
2872  * possible. */
2873 PETSC_STATIC_INLINE PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg,PetscInt count,PetscInt *PETSC_RESTRICT *slot) {return PetscSegBufferGet(seg,(size_t)count,(void**)slot);}
2874 
2875 PETSC_EXTERN PetscSegBuffer PetscCitationsList;
2876 #undef __FUNCT__
2877 #define __FUNCT__ "PetscCitationsRegister"
2878 /*@C
2879       PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code.
2880 
2881      Not Collective - only what is registered on rank 0 of PETSC_COMM_WORLD will be printed
2882 
2883      Input Parameters:
2884 +      cite - the bibtex item, formated to displayed on multiple lines nicely
2885 -      set - a boolean variable initially set to PETSC_FALSE; this is used to insure only a single registration of the citation
2886 
2887    Level: intermediate
2888 
2889      Options Database:
2890 .     -citations [filenmae]   - print out the bibtex entries for the given computation
2891 @*/
2892 PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2893 {
2894   size_t         len;
2895   char           *vstring;
2896   PetscErrorCode ierr;
2897 
2898   PetscFunctionBegin;
2899   if (set && *set) PetscFunctionReturn(0);
2900   ierr = PetscStrlen(cit,&len);CHKERRQ(ierr);
2901   ierr = PetscSegBufferGet(PetscCitationsList,len,&vstring);CHKERRQ(ierr);
2902   ierr = PetscMemcpy(vstring,cit,len);CHKERRQ(ierr);
2903   if (set) *set = PETSC_TRUE;
2904   PetscFunctionReturn(0);
2905 }
2906 
2907 PETSC_EXTERN PetscErrorCode PetscURLShorten(const char[],char[],size_t);
2908 PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm,char[],char[],size_t);
2909 PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm,const char[],char[],size_t);
2910 PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm,const char[],const char []);
2911 
2912 PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm,char[],char[],size_t);
2913 PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm,const char[],char[],char[],size_t);
2914 
2915 PETSC_EXTERN PetscErrorCode PetscTextBelt(MPI_Comm,const char[],const char[],PetscBool*);
2916 
2917 PETSC_EXTERN PetscErrorCode PetscPullJSONValue(const char[],const char[],char[],size_t,PetscBool*);
2918 PETSC_EXTERN PetscErrorCode PetscPushJSONValue(char[],const char[],const char[],size_t);
2919 
2920 
2921 #if defined(PETSC_USE_DEBUG)
2922 /*
2923    Verify that all processes in the communicator have called this from the same line of code
2924  */
2925 PETSC_EXTERN PetscErrorCode PetscAllreduceBarrierCheck(MPI_Comm,PetscMPIInt,int,const char*,const char *);
2926 #define MPIU_Allreduce(a,b,c,d,e,fcomm) (PetscAllreduceBarrierCheck(fcomm,c,__LINE__,__FUNCT__,__FILE__) || MPI_Allreduce(a,b,c,d,e,fcomm))
2927 #else
2928 #define MPIU_Allreduce(a,b,c,d,e,fcomm) MPI_Allreduce(a,b,c,d,e,fcomm)
2929 #endif
2930 
2931 /* Reset __FUNCT__ in case the user does not define it themselves */
2932 #undef __FUNCT__
2933 #define __FUNCT__ "User provided function"
2934 
2935 #endif
2936