Actual source code: petscsys.h
petsc-master 2021-01-19
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>
17: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP) || defined(PETSC_HAVE_KOKKOS)
18: #define PETSC_HAVE_DEVICE
19: #endif
21: #if defined(PETSC_DESIRE_FEATURE_TEST_MACROS)
22: /*
23: Feature test macros must be included before headers defined by IEEE Std 1003.1-2001
24: We only turn these in PETSc source files that require them by setting PETSC_DESIRE_FEATURE_TEST_MACROS
25: */
26: # if defined(PETSC__POSIX_C_SOURCE_200112L) && !defined(_POSIX_C_SOURCE)
27: # define _POSIX_C_SOURCE 200112L
28: # endif
29: # if defined(PETSC__BSD_SOURCE) && !defined(_BSD_SOURCE)
30: # define _BSD_SOURCE
31: # endif
32: # if defined(PETSC__DEFAULT_SOURCE) && !defined(_DEFAULT_SOURCE)
33: # define _DEFAULT_SOURCE
34: # endif
35: # if defined(PETSC__GNU_SOURCE) && !defined(_GNU_SOURCE)
36: # define _GNU_SOURCE
37: # endif
38: #endif
40: #include <petscsystypes.h>
42: /* ========================================================================== */
43: /*
44: This facilitates using the C version of PETSc from C++ and the C++ version from C.
45: */
46: #if defined(__cplusplus)
47: # define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_CXX
48: #else
49: # define PETSC_FUNCTION_NAME PETSC_FUNCTION_NAME_C
50: #endif
52: /* ========================================================================== */
53: /*
54: Since PETSc manages its own extern "C" handling users should never include PETSc include
55: files within extern "C". This will generate a compiler error if a user does put the include
56: file within an extern "C".
57: */
58: #if defined(__cplusplus)
59: void assert_never_put_petsc_headers_inside_an_extern_c(int); void assert_never_put_petsc_headers_inside_an_extern_c(double);
60: #endif
62: #if defined(__cplusplus)
63: # define PETSC_RESTRICT PETSC_CXX_RESTRICT
64: #else
65: # define PETSC_RESTRICT PETSC_C_RESTRICT
66: #endif
68: #if defined(__cplusplus)
69: # define PETSC_INLINE PETSC_CXX_INLINE
70: #else
71: # define PETSC_INLINE PETSC_C_INLINE
72: #endif
74: #define PETSC_STATIC_INLINE static PETSC_INLINE
76: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES) /* For Win32 shared libraries */
77: # define __declspec(dllexport)
78: # define PETSC_DLLIMPORT __declspec(dllimport)
79: # define PETSC_VISIBILITY_INTERNAL
80: #elif defined(PETSC_USE_VISIBILITY_CXX) && defined(__cplusplus)
81: # define __attribute__((visibility ("default")))
82: # define PETSC_DLLIMPORT __attribute__((visibility ("default")))
83: # define PETSC_VISIBILITY_INTERNAL __attribute__((visibility ("hidden")))
84: #elif defined(PETSC_USE_VISIBILITY_C) && !defined(__cplusplus)
85: # define __attribute__((visibility ("default")))
86: # define PETSC_DLLIMPORT __attribute__((visibility ("default")))
87: # define PETSC_VISIBILITY_INTERNAL __attribute__((visibility ("hidden")))
88: #else
89: # define
90: # define PETSC_DLLIMPORT
91: # define PETSC_VISIBILITY_INTERNAL
92: #endif
94: #if defined(petsc_EXPORTS) /* CMake defines this when building the shared library */
95: # define PETSC_VISIBILITY_PUBLIC
96: #else /* Win32 users need this to import symbols from petsc.dll */
97: # define PETSC_VISIBILITY_PUBLIC PETSC_DLLIMPORT
98: #endif
100: /*
101: Functions tagged with PETSC_EXTERN in the header files are
102: always defined as extern "C" when compiled with C++ so they may be
103: used from C and are always visible in the shared libraries
104: */
105: #if defined(__cplusplus)
106: # define PETSC_EXTERN extern "C" PETSC_VISIBILITY_PUBLIC
107: # define PETSC_EXTERN_TYPEDEF extern "C"
108: # define PETSC_INTERN extern "C" PETSC_VISIBILITY_INTERNAL
109: #else
110: # define PETSC_EXTERN extern PETSC_VISIBILITY_PUBLIC
111: # define PETSC_EXTERN_TYPEDEF
112: # define PETSC_INTERN extern PETSC_VISIBILITY_INTERNAL
113: #endif
115: #include <petscversion.h>
116: #define PETSC_AUTHOR_INFO " The PETSc Team\n petsc-maint@mcs.anl.gov\n https://www.mcs.anl.gov/petsc/\n"
118: /* ========================================================================== */
120: /*
121: Defines the interface to MPI allowing the use of all MPI functions.
123: PETSc does not use the C++ binding of MPI at ALL. The following flag
124: makes sure the C++ bindings are not included. The C++ bindings REQUIRE
125: putting mpi.h before ANY C++ include files, we cannot control this
126: with all PETSc users. Users who want to use the MPI C++ bindings can include
127: mpicxx.h directly in their code
128: */
129: #if !defined(MPICH_SKIP_MPICXX)
130: # define MPICH_SKIP_MPICXX 1
131: #endif
132: #if !defined(OMPI_SKIP_MPICXX)
133: # define OMPI_SKIP_MPICXX 1
134: #endif
135: #if defined(PETSC_HAVE_MPIUNI)
136: # include <petsc/mpiuni/mpi.h>
137: #else
138: # include <mpi.h>
139: #endif
141: /*MC
142: PetscDefined - determine whether a boolean macro is defined
144: Notes:
145: The prefix "PETSC_" is added to the argument.
147: Typical usage is within normal code,
149: $ if (PetscDefined(USE_DEBUG)) { ... }
151: but can also be used in the preprocessor,
153: $ #if PetscDefined(USE_DEBUG)
154: $ ...
155: $ #else
157: Either way, it evaluates true if PETSC_USE_DEBUG is defined (merely defined or defined to 1), and false if PETSC_USE_DEBUG is undefined. This macro
158: should not be used if its argument may be defined to a non-empty value other than 1.
160: To avoid prepending "PETSC_", say to add custom checks in user code, one can use e.g.
162: $ #define FooDefined(d) PetscDefined_(FOO_ ## d)
164: Developer Notes:
165: Getting something that works in C and CPP for an arg that may or may not be defined is tricky. Here, if we have
166: "#define PETSC_HAVE_BOOGER 1" we match on the placeholder define, insert the "0," for arg1 and generate the triplet
167: (0, 1, 0). Then the last step cherry picks the 2nd arg (a one). When PETSC_HAVE_BOOGER is not defined, we generate
168: a (... 1, 0) pair, and when the last step cherry picks the 2nd arg, we get a zero.
170: Our extra expansion via PetscDefined__take_second_expand() is needed with MSVC, which has a nonconforming
171: implementation of variadic macros.
173: Level: developer
174: M*/
175: #if !defined(PETSC_SKIP_VARIADIC_MACROS)
176: # define PetscDefined_arg_1 shift,
177: # define PetscDefined_arg_ shift,
178: # define PetscDefined__take_second_expanded(ignored, val, ...) val
179: # define PetscDefined__take_second_expand(args) PetscDefined__take_second_expanded args
180: # define PetscDefined__take_second(...) PetscDefined__take_second_expand((__VA_ARGS__))
181: # define PetscDefined___(arg1_or_junk) PetscDefined__take_second(arg1_or_junk 1, 0, at_)
182: # define PetscDefined__(value) PetscDefined___(PetscDefined_arg_ ## value)
183: # define PetscDefined_(d) PetscDefined__(d)
184: # define PetscDefined(d) PetscDefined_(PETSC_ ## d)
185: #endif
187: /*
188: Perform various sanity checks that the correct mpi.h is being included at compile time.
189: This usually happens because
190: * either an unexpected mpi.h is in the default compiler path (i.e. in /usr/include) or
191: * an extra include path -I/something (which contains the unexpected mpi.h) is being passed to the compiler
192: */
193: #if defined(PETSC_HAVE_MPIUNI)
194: # if !defined(MPIUNI_H)
195: # error "PETSc was configured with --with-mpi=0 but now appears to be compiling using a different mpi.h"
196: # endif
197: #elif defined(PETSC_HAVE_I_MPI_NUMVERSION)
198: # if !defined(I_MPI_NUMVERSION)
199: # error "PETSc was configured with I_MPI but now appears to be compiling using a non-I_MPI mpi.h"
200: # elif I_MPI_NUMVERSION != PETSC_HAVE_I_MPI_NUMVERSION
201: # error "PETSc was configured with one I_MPI mpi.h version but now appears to be compiling using a different I_MPI mpi.h version"
202: # endif
203: #elif defined(PETSC_HAVE_MVAPICH2_NUMVERSION)
204: # if !defined(MVAPICH2_NUMVERSION)
205: # error "PETSc was configured with MVAPICH2 but now appears to be compiling using a non-MVAPICH2 mpi.h"
206: # elif MVAPICH2_NUMVERSION != PETSC_HAVE_MVAPICH2_NUMVERSION
207: # error "PETSc was configured with one MVAPICH2 mpi.h version but now appears to be compiling using a different MVAPICH2 mpi.h version"
208: # endif
209: #elif defined(PETSC_HAVE_MPICH_NUMVERSION)
210: # if !defined(MPICH_NUMVERSION) || defined(MVAPICH2_NUMVERSION) || defined(I_MPI_NUMVERSION)
211: # error "PETSc was configured with MPICH but now appears to be compiling using a non-MPICH mpi.h"
212: # elif (MPICH_NUMVERSION/100000 != PETSC_HAVE_MPICH_NUMVERSION/100000) || (MPICH_NUMVERSION%100000/1000 < PETSC_HAVE_MPICH_NUMVERSION%100000/1000)
213: # error "PETSc was configured with one MPICH mpi.h version but now appears to be compiling using a different MPICH mpi.h version"
214: # endif
215: #elif defined(PETSC_HAVE_OMPI_MAJOR_VERSION)
216: # if !defined(OMPI_MAJOR_VERSION)
217: # error "PETSc was configured with OpenMPI but now appears to be compiling using a non-OpenMPI mpi.h"
218: # 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)
219: # error "PETSc was configured with one OpenMPI mpi.h version but now appears to be compiling using a different OpenMPI mpi.h version"
220: # endif
221: #elif defined(PETSC_HAVE_MSMPI_VERSION)
222: # if !defined(MSMPI_VER)
223: # error "PETSc was configured with MSMPI but now appears to be compiling using a non-MSMPI mpi.h"
224: # elif (MSMPI_VER != PETSC_HAVE_MSMPI_VERSION)
225: # error "PETSc was configured with one MSMPI mpi.h version but now appears to be compiling using a different MSMPI mpi.h version"
226: # endif
227: #elif defined(OMPI_MAJOR_VERSION) || defined(MPICH_NUMVERSION) || defined(MSMPI_VER)
228: # error "PETSc was configured with undetermined MPI - but now appears to be compiling using any of OpenMPI, MS-MPI or a MPICH variant"
229: #endif
231: /*
232: Need to put stdio.h AFTER mpi.h for MPICH2 with C++ compiler
233: see the top of mpicxx.h in the MPICH2 distribution.
234: */
235: #include <stdio.h>
237: /* MSMPI on 32bit windows requires this yukky hack - that breaks MPI standard compliance */
238: #if !defined(MPIAPI)
239: #define MPIAPI
240: #endif
242: /*
243: Support for Clang (>=3.2) matching type tag arguments with void* buffer types.
244: This allows the compiler to detect cases where the MPI datatype argument passed to a MPI routine
245: does not match the actual type of the argument being passed in
246: */
247: #if defined(__has_attribute) && defined(works_with_const_which_is_not_true)
248: # if __has_attribute(argument_with_type_tag) && __has_attribute(pointer_with_type_tag) && __has_attribute(type_tag_for_datatype)
249: # define PetscAttrMPIPointerWithType(bufno,typeno) __attribute__((pointer_with_type_tag(MPI,bufno,typeno)))
250: # define PetscAttrMPITypeTag(type) __attribute__((type_tag_for_datatype(MPI,type)))
251: # define PetscAttrMPITypeTagLayoutCompatible(type) __attribute__((type_tag_for_datatype(MPI,type,layout_compatible)))
252: # endif
253: #endif
254: #if !defined(PetscAttrMPIPointerWithType)
255: # define PetscAttrMPIPointerWithType(bufno,typeno)
256: # define PetscAttrMPITypeTag(type)
257: # define PetscAttrMPITypeTagLayoutCompatible(type)
258: #endif
260: PETSC_EXTERN MPI_Datatype MPIU_ENUM PetscAttrMPITypeTag(PetscEnum);
261: PETSC_EXTERN MPI_Datatype MPIU_BOOL PetscAttrMPITypeTag(PetscBool);
263: /*MC
264: MPIU_INT - MPI datatype corresponding to PetscInt
266: Notes:
267: In MPI calls that require an MPI datatype that matches a PetscInt or array of PetscInt values, pass this value.
269: Level: beginner
271: .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_SCALAR, MPIU_COMPLEX
272: M*/
274: #if defined(PETSC_HAVE_STDINT_H) && defined(PETSC_HAVE_INTTYPES_H) && defined(PETSC_HAVE_MPI_INT64_T) /* MPI_INT64_T is not guaranteed to be a macro */
275: # define MPIU_INT64 MPI_INT64_T
276: # define PetscInt64_FMT PRId64
277: #elif (PETSC_SIZEOF_LONG_LONG == 8)
278: # define MPIU_INT64 MPI_LONG_LONG_INT
279: # define PetscInt64_FMT "lld"
280: #elif defined(PETSC_HAVE___INT64)
281: # define MPIU_INT64 MPI_INT64_T
282: # define PetscInt64_FMT "ld"
283: #else
284: # error "cannot determine PetscInt64 type"
285: #endif
287: PETSC_EXTERN MPI_Datatype MPIU_FORTRANADDR;
289: #if defined(PETSC_USE_64BIT_INDICES)
290: # define MPIU_INT MPIU_INT64
291: # define PetscInt_FMT PetscInt64_FMT
292: #else
293: # define MPIU_INT MPI_INT
294: # define PetscInt_FMT "d"
295: #endif
297: /*
298: For the rare cases when one needs to send a size_t object with MPI
299: */
300: PETSC_EXTERN MPI_Datatype MPIU_SIZE_T;
302: /*
303: You can use PETSC_STDOUT as a replacement of stdout. You can also change
304: the value of PETSC_STDOUT to redirect all standard output elsewhere
305: */
306: PETSC_EXTERN FILE* PETSC_STDOUT;
308: /*
309: You can use PETSC_STDERR as a replacement of stderr. You can also change
310: the value of PETSC_STDERR to redirect all standard error elsewhere
311: */
312: PETSC_EXTERN FILE* PETSC_STDERR;
314: /*MC
315: PetscUnlikely - hints the compiler that the given condition is usually FALSE
317: Synopsis:
318: #include <petscsys.h>
319: PetscBool PetscUnlikely(PetscBool cond)
321: Not Collective
323: Input Parameters:
324: . cond - condition or expression
326: Notes:
327: This returns the same truth value, it is only a hint to compilers that the resulting
328: branch is unlikely.
330: Level: advanced
332: .seealso: PetscUnlikelyDebug(), PetscLikely(), CHKERRQ
333: M*/
335: /*MC
336: PetscLikely - hints the compiler that the given condition is usually TRUE
338: Synopsis:
339: #include <petscsys.h>
340: PetscBool PetscLikely(PetscBool cond)
342: Not Collective
344: Input Parameters:
345: . cond - condition or expression
347: Notes:
348: This returns the same truth value, it is only a hint to compilers that the resulting
349: branch is likely.
351: Level: advanced
353: .seealso: PetscUnlikely()
354: M*/
355: #if defined(PETSC_HAVE_BUILTIN_EXPECT)
356: # define PetscUnlikely(cond) __builtin_expect(!!(cond),0)
357: # define PetscLikely(cond) __builtin_expect(!!(cond),1)
358: #else
359: # define PetscUnlikely(cond) (cond)
360: # define PetscLikely(cond) (cond)
361: #endif
363: /*MC
364: PetscUnlikelyDebug - hints the compiler that the given condition is usually FALSE, eliding the check in optimized mode
366: Synopsis:
367: #include <petscsys.h>
368: PetscBool PetscUnlikelyDebug(PetscBool cond)
370: Not Collective
372: Input Parameters:
373: . cond - condition or expression
375: Notes:
376: This returns the same truth value, it is only a hint to compilers that the resulting
377: branch is unlikely. When compiled in optimized mode, it always returns false.
379: Level: advanced
381: .seealso: PetscUnlikely(), CHKERRQ, SETERRQ
382: M*/
383: #define PetscUnlikelyDebug(cond) (PetscDefined(USE_DEBUG) && PetscUnlikely(cond))
385: /* PetscPragmaSIMD - from CeedPragmaSIMD */
387: #if defined(__INTEL_COMPILER) && !defined(_WIN32)
388: # define PetscPragmaSIMD _Pragma("vector")
389: #elif defined(__GNUC__) && __GNUC__ >= 5 && !defined(__PGI)
390: # define PetscPragmaSIMD _Pragma("GCC ivdep")
391: #elif defined(_OPENMP) && _OPENMP >= 201307 && !defined(_WIN32)
392: # define PetscPragmaSIMD _Pragma("omp simd")
393: #elif defined(_OPENMP) && _OPENMP >= 201307 && defined(_WIN32)
394: # define PetscPragmaSIMD __pragma(omp simd)
395: #elif defined(PETSC_HAVE_CRAY_VECTOR)
396: # define PetscPragmaSIMD _Pragma("_CRI ivdep")
397: #else
398: # define PetscPragmaSIMD
399: #endif
401: /*
402: Declare extern C stuff after including external header files
403: */
405: PETSC_EXTERN const char *const PetscBools[];
407: /*
408: Defines elementary mathematics functions and constants.
409: */
410: #include <petscmath.h>
412: PETSC_EXTERN const char *const PetscCopyModes[];
414: /*MC
415: PETSC_IGNORE - same as NULL, means PETSc will ignore this argument
417: Level: beginner
419: Note:
420: Accepted by many PETSc functions to not set a parameter and instead use
421: some default
423: Fortran Notes:
424: This macro does not exist in Fortran; you must use PETSC_NULL_INTEGER,
425: PETSC_NULL_DOUBLE_PRECISION etc
427: .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_DETERMINE
429: M*/
430: #define PETSC_IGNORE NULL
432: /* This is deprecated */
433: #define PETSC_NULL NULL
435: /*MC
436: PETSC_DECIDE - standard way of passing in integer or floating point parameter
437: where you wish PETSc to use the default.
439: Level: beginner
441: .seealso: PETSC_DEFAULT, PETSC_IGNORE, PETSC_DETERMINE
443: M*/
444: #define PETSC_DECIDE -1
446: /*MC
447: PETSC_DETERMINE - standard way of passing in integer or floating point parameter
448: where you wish PETSc to compute the required value.
450: Level: beginner
452: Developer Note:
453: I would like to use const PetscInt PETSC_DETERMINE = PETSC_DECIDE; but for
454: some reason this is not allowed by the standard even though PETSC_DECIDE is a constant value.
456: .seealso: PETSC_DECIDE, PETSC_DEFAULT, PETSC_IGNORE, VecSetSizes()
458: M*/
459: #define PETSC_DETERMINE PETSC_DECIDE
461: /*MC
462: PETSC_DEFAULT - standard way of passing in integer or floating point parameter
463: where you wish PETSc to use the default.
465: Level: beginner
467: Fortran Notes:
468: You need to use PETSC_DEFAULT_INTEGER or PETSC_DEFAULT_REAL.
470: .seealso: PETSC_DECIDE, PETSC_IGNORE, PETSC_DETERMINE
472: M*/
473: #define PETSC_DEFAULT -2
475: /*MC
476: PETSC_COMM_WORLD - the equivalent of the MPI_COMM_WORLD communicator which represents
477: all the processes that PETSc knows about.
479: Level: beginner
481: Notes:
482: By default PETSC_COMM_WORLD and MPI_COMM_WORLD are identical unless you wish to
483: run PETSc on ONLY a subset of MPI_COMM_WORLD. In that case create your new (smaller)
484: communicator, call it, say comm, and set PETSC_COMM_WORLD = comm BEFORE calling
485: PetscInitialize(), but after MPI_Init() has been called.
487: The value of PETSC_COMM_WORLD should never be USED/accessed before PetscInitialize()
488: is called because it may not have a valid value yet.
490: .seealso: PETSC_COMM_SELF
492: M*/
493: PETSC_EXTERN MPI_Comm PETSC_COMM_WORLD;
495: /*MC
496: PETSC_COMM_SELF - This is always MPI_COMM_SELF
498: Level: beginner
500: Notes:
501: Do not USE/access or set this variable before PetscInitialize() has been called.
503: .seealso: PETSC_COMM_WORLD
505: M*/
506: #define PETSC_COMM_SELF MPI_COMM_SELF
508: /*MC
509: PETSC_MPI_THREAD_REQUIRED - the required threading support used if PETSc initializes
510: MPI with MPI_Init_thread.
512: Level: beginner
514: Notes:
515: By default PETSC_MPI_THREAD_REQUIRED equals MPI_THREAD_FUNNELED.
517: .seealso: PetscInitialize()
519: M*/
520: PETSC_EXTERN PetscMPIInt PETSC_MPI_THREAD_REQUIRED;
522: PETSC_EXTERN PetscBool PetscBeganMPI;
523: PETSC_EXTERN PetscBool PetscErrorHandlingInitialized;
524: PETSC_EXTERN PetscBool PetscInitializeCalled;
525: PETSC_EXTERN PetscBool PetscFinalizeCalled;
526: PETSC_EXTERN PetscBool PetscViennaCLSynchronize;
528: PETSC_EXTERN PetscErrorCode PetscSetHelpVersionFunctions(PetscErrorCode (*)(MPI_Comm),PetscErrorCode (*)(MPI_Comm));
529: PETSC_EXTERN PetscErrorCode PetscCommDuplicate(MPI_Comm,MPI_Comm*,int*);
530: PETSC_EXTERN PetscErrorCode PetscCommDestroy(MPI_Comm*);
532: #if defined(PETSC_HAVE_CUDA)
533: PETSC_EXTERN PetscBool PetscCUDASynchronize;
534: PETSC_EXTERN PetscErrorCode PetscCUDAInitialize(MPI_Comm,PetscInt);
535: PETSC_EXTERN PetscErrorCode PetscCUDAInitializeCheck(void);
536: #endif
538: #if defined(PETSC_HAVE_HIP)
539: PETSC_EXTERN PetscBool PetscHIPSynchronize;
540: PETSC_EXTERN PetscErrorCode PetscHIPInitialize(MPI_Comm,PetscInt);
541: PETSC_EXTERN PetscErrorCode PetscHIPInitializeCheck(void);
542: #endif
544: #if defined(PETSC_HAVE_KOKKOS)
545: PETSC_EXTERN PetscErrorCode PetscKokkosInitializeCheck(void); /* Initialize Kokkos if not yet. */
546: #endif
548: #if defined(PETSC_HAVE_ELEMENTAL)
549: PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
550: PETSC_EXTERN PetscErrorCode PetscElementalInitialized(PetscBool*);
551: PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
552: #endif
554: /*MC
555: PetscMalloc - Allocates memory, One should use PetscNew(), PetscMalloc1() or PetscCalloc1() usually instead of this
557: Synopsis:
558: #include <petscsys.h>
559: PetscErrorCode PetscMalloc(size_t m,void **result)
561: Not Collective
563: Input Parameter:
564: . m - number of bytes to allocate
566: Output Parameter:
567: . result - memory allocated
569: Level: beginner
571: Notes:
572: Memory is always allocated at least double aligned
574: It is safe to allocate size 0 and pass the resulting pointer (which may or may not be NULL) to PetscFree().
576: .seealso: PetscFree(), PetscNew()
578: M*/
579: #define PetscMalloc(a,b) ((*PetscTrMalloc)((a),PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(void**)(b)))
581: /*MC
582: PetscRealloc - Rellocates memory
584: Synopsis:
585: #include <petscsys.h>
586: PetscErrorCode PetscRealloc(size_t m,void **result)
588: Not Collective
590: Input Parameters:
591: + m - number of bytes to allocate
592: - result - previous memory
594: Output Parameter:
595: . result - new memory allocated
597: Level: developer
599: Notes:
600: Memory is always allocated at least double aligned
602: .seealso: PetscMalloc(), PetscFree(), PetscNew()
604: M*/
605: #define PetscRealloc(a,b) ((*PetscTrRealloc)((a),__LINE__,PETSC_FUNCTION_NAME,__FILE__,(void**)(b)))
607: /*MC
608: PetscAddrAlign - Rounds up an address to PETSC_MEMALIGN alignment
610: Synopsis:
611: #include <petscsys.h>
612: void *PetscAddrAlign(void *addr)
614: Not Collective
616: Input Parameters:
617: . addr - address to align (any pointer type)
619: Level: developer
621: .seealso: PetscMallocAlign()
623: M*/
624: #define PetscAddrAlign(a) (void*)((((PETSC_UINTPTR_T)(a))+(PETSC_MEMALIGN-1)) & ~(PETSC_MEMALIGN-1))
626: /*MC
627: PetscMalloc1 - Allocates an array of memory aligned to PETSC_MEMALIGN
629: Synopsis:
630: #include <petscsys.h>
631: PetscErrorCode PetscMalloc1(size_t m1,type **r1)
633: Not Collective
635: Input Parameter:
636: . m1 - number of elements to allocate (may be zero)
638: Output Parameter:
639: . r1 - memory allocated
641: Note:
642: This uses the sizeof() of the memory type requested to determine the total memory to be allocated, therefore you should not
643: multiply the number of elements requested by the sizeof() the type. For example use
644: $ PetscInt *id;
645: $ PetscMalloc1(10,&id);
646: not
647: $ PetscInt *id;
648: $ PetscMalloc1(10*sizeof(PetscInt),&id);
650: Does not zero the memory allocated, use PetscCalloc1() to obtain memory that has been zeroed.
652: Level: beginner
654: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc1(), PetscMalloc2()
656: M*/
657: #define PetscMalloc1(m1,r1) PetscMallocA(1,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1))
659: /*MC
660: PetscCalloc1 - Allocates a cleared (zeroed) array of memory aligned to PETSC_MEMALIGN
662: Synopsis:
663: #include <petscsys.h>
664: PetscErrorCode PetscCalloc1(size_t m1,type **r1)
666: Not Collective
668: Input Parameter:
669: . m1 - number of elements to allocate in 1st chunk (may be zero)
671: Output Parameter:
672: . r1 - memory allocated
674: Notes:
675: See PetsMalloc1() for more details on usage.
677: Level: beginner
679: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc1(), PetscCalloc2()
681: M*/
682: #define PetscCalloc1(m1,r1) PetscMallocA(1,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1))
684: /*MC
685: PetscMalloc2 - Allocates 2 arrays of memory both aligned to PETSC_MEMALIGN
687: Synopsis:
688: #include <petscsys.h>
689: PetscErrorCode PetscMalloc2(size_t m1,type **r1,size_t m2,type **r2)
691: Not Collective
693: Input Parameter:
694: + m1 - number of elements to allocate in 1st chunk (may be zero)
695: - m2 - number of elements to allocate in 2nd chunk (may be zero)
697: Output Parameter:
698: + r1 - memory allocated in first chunk
699: - r2 - memory allocated in second chunk
701: Level: developer
703: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc1(), PetscCalloc2()
705: M*/
706: #define PetscMalloc2(m1,r1,m2,r2) PetscMallocA(2,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2))
708: /*MC
709: PetscCalloc2 - Allocates 2 cleared (zeroed) arrays of memory both aligned to PETSC_MEMALIGN
711: Synopsis:
712: #include <petscsys.h>
713: PetscErrorCode PetscCalloc2(size_t m1,type **r1,size_t m2,type **r2)
715: Not Collective
717: Input Parameter:
718: + m1 - number of elements to allocate in 1st chunk (may be zero)
719: - m2 - number of elements to allocate in 2nd chunk (may be zero)
721: Output Parameter:
722: + r1 - memory allocated in first chunk
723: - r2 - memory allocated in second chunk
725: Level: developer
727: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc1(), PetscMalloc2()
729: M*/
730: #define PetscCalloc2(m1,r1,m2,r2) PetscMallocA(2,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2))
732: /*MC
733: PetscMalloc3 - Allocates 3 arrays of memory, all aligned to PETSC_MEMALIGN
735: Synopsis:
736: #include <petscsys.h>
737: PetscErrorCode PetscMalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
739: Not Collective
741: Input Parameter:
742: + m1 - number of elements to allocate in 1st chunk (may be zero)
743: . m2 - number of elements to allocate in 2nd chunk (may be zero)
744: - m3 - number of elements to allocate in 3rd chunk (may be zero)
746: Output Parameter:
747: + r1 - memory allocated in first chunk
748: . r2 - memory allocated in second chunk
749: - r3 - memory allocated in third chunk
751: Level: developer
753: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc3(), PetscFree3()
755: M*/
756: #define PetscMalloc3(m1,r1,m2,r2,m3,r3) PetscMallocA(3,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3))
758: /*MC
759: PetscCalloc3 - Allocates 3 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN
761: Synopsis:
762: #include <petscsys.h>
763: PetscErrorCode PetscCalloc3(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3)
765: Not Collective
767: Input Parameter:
768: + m1 - number of elements to allocate in 1st chunk (may be zero)
769: . m2 - number of elements to allocate in 2nd chunk (may be zero)
770: - m3 - number of elements to allocate in 3rd chunk (may be zero)
772: Output Parameter:
773: + r1 - memory allocated in first chunk
774: . r2 - memory allocated in second chunk
775: - r3 - memory allocated in third chunk
777: Level: developer
779: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscCalloc2(), PetscMalloc3(), PetscFree3()
781: M*/
782: #define PetscCalloc3(m1,r1,m2,r2,m3,r3) PetscMallocA(3,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3))
784: /*MC
785: PetscMalloc4 - Allocates 4 arrays of memory, all aligned to PETSC_MEMALIGN
787: Synopsis:
788: #include <petscsys.h>
789: PetscErrorCode PetscMalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
791: Not Collective
793: Input Parameter:
794: + m1 - number of elements to allocate in 1st chunk (may be zero)
795: . m2 - number of elements to allocate in 2nd chunk (may be zero)
796: . m3 - number of elements to allocate in 3rd chunk (may be zero)
797: - m4 - number of elements to allocate in 4th chunk (may be zero)
799: Output Parameter:
800: + r1 - memory allocated in first chunk
801: . r2 - memory allocated in second chunk
802: . r3 - memory allocated in third chunk
803: - r4 - memory allocated in fourth chunk
805: Level: developer
807: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc4(), PetscFree4()
809: M*/
810: #define PetscMalloc4(m1,r1,m2,r2,m3,r3,m4,r4) PetscMallocA(4,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4))
812: /*MC
813: PetscCalloc4 - Allocates 4 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN
815: Synopsis:
816: #include <petscsys.h>
817: PetscErrorCode PetscCalloc4(size_t m1,type **r1,size_t m2,type **r2,size_t m3,type **r3,size_t m4,type **r4)
819: Not Collective
821: Input Parameters:
822: + m1 - number of elements to allocate in 1st chunk (may be zero)
823: . m2 - number of elements to allocate in 2nd chunk (may be zero)
824: . m3 - number of elements to allocate in 3rd chunk (may be zero)
825: - m4 - number of elements to allocate in 4th chunk (may be zero)
827: Output Parameters:
828: + r1 - memory allocated in first chunk
829: . r2 - memory allocated in second chunk
830: . r3 - memory allocated in third chunk
831: - r4 - memory allocated in fourth chunk
833: Level: developer
835: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc4(), PetscFree4()
837: M*/
838: #define PetscCalloc4(m1,r1,m2,r2,m3,r3,m4,r4) PetscMallocA(4,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4))
840: /*MC
841: PetscMalloc5 - Allocates 5 arrays of memory, all aligned to PETSC_MEMALIGN
843: Synopsis:
844: #include <petscsys.h>
845: 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)
847: Not Collective
849: Input Parameters:
850: + m1 - number of elements to allocate in 1st chunk (may be zero)
851: . m2 - number of elements to allocate in 2nd chunk (may be zero)
852: . m3 - number of elements to allocate in 3rd chunk (may be zero)
853: . m4 - number of elements to allocate in 4th chunk (may be zero)
854: - m5 - number of elements to allocate in 5th chunk (may be zero)
856: Output Parameters:
857: + r1 - memory allocated in first chunk
858: . r2 - memory allocated in second chunk
859: . r3 - memory allocated in third chunk
860: . r4 - memory allocated in fourth chunk
861: - r5 - memory allocated in fifth chunk
863: Level: developer
865: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc5(), PetscFree5()
867: M*/
868: #define PetscMalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5) PetscMallocA(5,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4),(size_t)(m5)*sizeof(**(r5)),(r5))
870: /*MC
871: PetscCalloc5 - Allocates 5 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN
873: Synopsis:
874: #include <petscsys.h>
875: 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)
877: Not Collective
879: Input Parameters:
880: + m1 - number of elements to allocate in 1st chunk (may be zero)
881: . m2 - number of elements to allocate in 2nd chunk (may be zero)
882: . m3 - number of elements to allocate in 3rd chunk (may be zero)
883: . m4 - number of elements to allocate in 4th chunk (may be zero)
884: - m5 - number of elements to allocate in 5th chunk (may be zero)
886: Output Parameters:
887: + r1 - memory allocated in first chunk
888: . r2 - memory allocated in second chunk
889: . r3 - memory allocated in third chunk
890: . r4 - memory allocated in fourth chunk
891: - r5 - memory allocated in fifth chunk
893: Level: developer
895: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc5(), PetscFree5()
897: M*/
898: #define PetscCalloc5(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5) PetscMallocA(5,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4),(size_t)(m5)*sizeof(**(r5)),(r5))
900: /*MC
901: PetscMalloc6 - Allocates 6 arrays of memory, all aligned to PETSC_MEMALIGN
903: Synopsis:
904: #include <petscsys.h>
905: 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)
907: Not Collective
909: Input Parameters:
910: + m1 - number of elements to allocate in 1st chunk (may be zero)
911: . m2 - number of elements to allocate in 2nd chunk (may be zero)
912: . m3 - number of elements to allocate in 3rd chunk (may be zero)
913: . m4 - number of elements to allocate in 4th chunk (may be zero)
914: . m5 - number of elements to allocate in 5th chunk (may be zero)
915: - m6 - number of elements to allocate in 6th chunk (may be zero)
917: Output Parameteasr:
918: + r1 - memory allocated in first chunk
919: . r2 - memory allocated in second chunk
920: . r3 - memory allocated in third chunk
921: . r4 - memory allocated in fourth chunk
922: . r5 - memory allocated in fifth chunk
923: - r6 - memory allocated in sixth chunk
925: Level: developer
927: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc6(), PetscFree3(), PetscFree4(), PetscFree5(), PetscFree6()
929: M*/
930: #define PetscMalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) PetscMallocA(6,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4),(size_t)(m5)*sizeof(**(r5)),(r5),(size_t)(m6)*sizeof(**(r6)),(r6))
932: /*MC
933: PetscCalloc6 - Allocates 6 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN
935: Synopsis:
936: #include <petscsys.h>
937: 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)
939: Not Collective
941: Input Parameters:
942: + m1 - number of elements to allocate in 1st chunk (may be zero)
943: . m2 - number of elements to allocate in 2nd chunk (may be zero)
944: . m3 - number of elements to allocate in 3rd chunk (may be zero)
945: . m4 - number of elements to allocate in 4th chunk (may be zero)
946: . m5 - number of elements to allocate in 5th chunk (may be zero)
947: - m6 - number of elements to allocate in 6th chunk (may be zero)
949: Output Parameters:
950: + r1 - memory allocated in first chunk
951: . r2 - memory allocated in second chunk
952: . r3 - memory allocated in third chunk
953: . r4 - memory allocated in fourth chunk
954: . r5 - memory allocated in fifth chunk
955: - r6 - memory allocated in sixth chunk
957: Level: developer
959: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscMalloc6(), PetscFree6()
961: M*/
962: #define PetscCalloc6(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6) PetscMallocA(6,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4),(size_t)(m5)*sizeof(**(r5)),(r5),(size_t)(m6)*sizeof(**(r6)),(r6))
964: /*MC
965: PetscMalloc7 - Allocates 7 arrays of memory, all aligned to PETSC_MEMALIGN
967: Synopsis:
968: #include <petscsys.h>
969: 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)
971: Not Collective
973: Input Parameters:
974: + m1 - number of elements to allocate in 1st chunk (may be zero)
975: . m2 - number of elements to allocate in 2nd chunk (may be zero)
976: . m3 - number of elements to allocate in 3rd chunk (may be zero)
977: . m4 - number of elements to allocate in 4th chunk (may be zero)
978: . m5 - number of elements to allocate in 5th chunk (may be zero)
979: . m6 - number of elements to allocate in 6th chunk (may be zero)
980: - m7 - number of elements to allocate in 7th chunk (may be zero)
982: Output Parameters:
983: + r1 - memory allocated in first chunk
984: . r2 - memory allocated in second chunk
985: . r3 - memory allocated in third chunk
986: . r4 - memory allocated in fourth chunk
987: . r5 - memory allocated in fifth chunk
988: . r6 - memory allocated in sixth chunk
989: - r7 - memory allocated in seventh chunk
991: Level: developer
993: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscCalloc7(), PetscFree7()
995: M*/
996: #define PetscMalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) PetscMallocA(7,PETSC_FALSE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4),(size_t)(m5)*sizeof(**(r5)),(r5),(size_t)(m6)*sizeof(**(r6)),(r6),(size_t)(m7)*sizeof(**(r7)),(r7))
998: /*MC
999: PetscCalloc7 - Allocates 7 cleared (zeroed) arrays of memory, all aligned to PETSC_MEMALIGN
1001: Synopsis:
1002: #include <petscsys.h>
1003: 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)
1005: Not Collective
1007: Input Parameters:
1008: + m1 - number of elements to allocate in 1st chunk (may be zero)
1009: . m2 - number of elements to allocate in 2nd chunk (may be zero)
1010: . m3 - number of elements to allocate in 3rd chunk (may be zero)
1011: . m4 - number of elements to allocate in 4th chunk (may be zero)
1012: . m5 - number of elements to allocate in 5th chunk (may be zero)
1013: . m6 - number of elements to allocate in 6th chunk (may be zero)
1014: - m7 - number of elements to allocate in 7th chunk (may be zero)
1016: Output Parameters:
1017: + r1 - memory allocated in first chunk
1018: . r2 - memory allocated in second chunk
1019: . r3 - memory allocated in third chunk
1020: . r4 - memory allocated in fourth chunk
1021: . r5 - memory allocated in fifth chunk
1022: . r6 - memory allocated in sixth chunk
1023: - r7 - memory allocated in seventh chunk
1025: Level: developer
1027: .seealso: PetscFree(), PetscNew(), PetscMalloc(), PetscMalloc2(), PetscMalloc7(), PetscFree7()
1029: M*/
1030: #define PetscCalloc7(m1,r1,m2,r2,m3,r3,m4,r4,m5,r5,m6,r6,m7,r7) PetscMallocA(7,PETSC_TRUE,__LINE__,PETSC_FUNCTION_NAME,__FILE__,(size_t)(m1)*sizeof(**(r1)),(r1),(size_t)(m2)*sizeof(**(r2)),(r2),(size_t)(m3)*sizeof(**(r3)),(r3),(size_t)(m4)*sizeof(**(r4)),(r4),(size_t)(m5)*sizeof(**(r5)),(r5),(size_t)(m6)*sizeof(**(r6)),(r6),(size_t)(m7)*sizeof(**(r7)),(r7))
1032: /*MC
1033: PetscNew - Allocates memory of a particular type, zeros the memory! Aligned to PETSC_MEMALIGN
1035: Synopsis:
1036: #include <petscsys.h>
1037: PetscErrorCode PetscNew(type **result)
1039: Not Collective
1041: Output Parameter:
1042: . result - memory allocated, sized to match pointer type
1044: Level: beginner
1046: .seealso: PetscFree(), PetscMalloc(), PetscNewLog(), PetscCalloc1(), PetscMalloc1()
1048: M*/
1049: #define PetscNew(b) PetscCalloc1(1,(b))
1051: /*MC
1052: PetscNewLog - Allocates memory of a type matching pointer, zeros the memory! Aligned to PETSC_MEMALIGN. Associates the memory allocated
1053: with the given object using PetscLogObjectMemory().
1055: Synopsis:
1056: #include <petscsys.h>
1057: PetscErrorCode PetscNewLog(PetscObject obj,type **result)
1059: Not Collective
1061: Input Parameter:
1062: . obj - object memory is logged to
1064: Output Parameter:
1065: . result - memory allocated, sized to match pointer type
1067: Level: developer
1069: .seealso: PetscFree(), PetscMalloc(), PetscNew(), PetscLogObjectMemory(), PetscCalloc1(), PetscMalloc1()
1071: M*/
1072: #define PetscNewLog(o,b) (PetscNew((b)) || PetscLogObjectMemory((PetscObject)o,sizeof(**(b))))
1074: /*MC
1075: PetscFree - Frees memory
1077: Synopsis:
1078: #include <petscsys.h>
1079: PetscErrorCode PetscFree(void *memory)
1081: Not Collective
1083: Input Parameter:
1084: . memory - memory to free (the pointer is ALWAYS set to NULL upon success)
1086: Level: beginner
1088: Notes:
1089: Do not free memory obtained with PetscMalloc2(), PetscCalloc2() etc, they must be freed with PetscFree2() etc.
1091: It is safe to call PetscFree() on a NULL pointer.
1093: .seealso: PetscNew(), PetscMalloc(), PetscNewLog(), PetscMalloc1(), PetscCalloc1()
1095: M*/
1096: #define PetscFree(a) ((*PetscTrFree)((void*)(a),__LINE__,PETSC_FUNCTION_NAME,__FILE__) || ((a) = NULL,0))
1098: /*MC
1099: PetscFree2 - Frees 2 chunks of memory obtained with PetscMalloc2()
1101: Synopsis:
1102: #include <petscsys.h>
1103: PetscErrorCode PetscFree2(void *memory1,void *memory2)
1105: Not Collective
1107: Input Parameters:
1108: + memory1 - memory to free
1109: - memory2 - 2nd memory to free
1111: Level: developer
1113: Notes:
1114: Memory must have been obtained with PetscMalloc2()
1116: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree()
1118: M*/
1119: #define PetscFree2(m1,m2) PetscFreeA(2,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2))
1121: /*MC
1122: PetscFree3 - Frees 3 chunks of memory obtained with PetscMalloc3()
1124: Synopsis:
1125: #include <petscsys.h>
1126: PetscErrorCode PetscFree3(void *memory1,void *memory2,void *memory3)
1128: Not Collective
1130: Input Parameters:
1131: + memory1 - memory to free
1132: . memory2 - 2nd memory to free
1133: - memory3 - 3rd memory to free
1135: Level: developer
1137: Notes:
1138: Memory must have been obtained with PetscMalloc3()
1140: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3()
1142: M*/
1143: #define PetscFree3(m1,m2,m3) PetscFreeA(3,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3))
1145: /*MC
1146: PetscFree4 - Frees 4 chunks of memory obtained with PetscMalloc4()
1148: Synopsis:
1149: #include <petscsys.h>
1150: PetscErrorCode PetscFree4(void *m1,void *m2,void *m3,void *m4)
1152: Not Collective
1154: Input Parameters:
1155: + m1 - memory to free
1156: . m2 - 2nd memory to free
1157: . m3 - 3rd memory to free
1158: - m4 - 4th memory to free
1160: Level: developer
1162: Notes:
1163: Memory must have been obtained with PetscMalloc4()
1165: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4()
1167: M*/
1168: #define PetscFree4(m1,m2,m3,m4) PetscFreeA(4,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4))
1170: /*MC
1171: PetscFree5 - Frees 5 chunks of memory obtained with PetscMalloc5()
1173: Synopsis:
1174: #include <petscsys.h>
1175: PetscErrorCode PetscFree5(void *m1,void *m2,void *m3,void *m4,void *m5)
1177: Not Collective
1179: Input Parameters:
1180: + m1 - memory to free
1181: . m2 - 2nd memory to free
1182: . m3 - 3rd memory to free
1183: . m4 - 4th memory to free
1184: - m5 - 5th memory to free
1186: Level: developer
1188: Notes:
1189: Memory must have been obtained with PetscMalloc5()
1191: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5()
1193: M*/
1194: #define PetscFree5(m1,m2,m3,m4,m5) PetscFreeA(5,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5))
1196: /*MC
1197: PetscFree6 - Frees 6 chunks of memory obtained with PetscMalloc6()
1199: Synopsis:
1200: #include <petscsys.h>
1201: PetscErrorCode PetscFree6(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6)
1203: Not Collective
1205: Input Parameters:
1206: + m1 - memory to free
1207: . m2 - 2nd memory to free
1208: . m3 - 3rd memory to free
1209: . m4 - 4th memory to free
1210: . m5 - 5th memory to free
1211: - m6 - 6th memory to free
1214: Level: developer
1216: Notes:
1217: Memory must have been obtained with PetscMalloc6()
1219: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6()
1221: M*/
1222: #define PetscFree6(m1,m2,m3,m4,m5,m6) PetscFreeA(6,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5),&(m6))
1224: /*MC
1225: PetscFree7 - Frees 7 chunks of memory obtained with PetscMalloc7()
1227: Synopsis:
1228: #include <petscsys.h>
1229: PetscErrorCode PetscFree7(void *m1,void *m2,void *m3,void *m4,void *m5,void *m6,void *m7)
1231: Not Collective
1233: Input Parameters:
1234: + m1 - memory to free
1235: . m2 - 2nd memory to free
1236: . m3 - 3rd memory to free
1237: . m4 - 4th memory to free
1238: . m5 - 5th memory to free
1239: . m6 - 6th memory to free
1240: - m7 - 7th memory to free
1243: Level: developer
1245: Notes:
1246: Memory must have been obtained with PetscMalloc7()
1248: .seealso: PetscNew(), PetscMalloc(), PetscMalloc2(), PetscFree(), PetscMalloc3(), PetscMalloc4(), PetscMalloc5(), PetscMalloc6(),
1249: PetscMalloc7()
1251: M*/
1252: #define PetscFree7(m1,m2,m3,m4,m5,m6,m7) PetscFreeA(7,__LINE__,PETSC_FUNCTION_NAME,__FILE__,&(m1),&(m2),&(m3),&(m4),&(m5),&(m6),&(m7))
1254: PETSC_EXTERN PetscErrorCode PetscMallocA(int,PetscBool,int,const char *,const char *,size_t,void *,...);
1255: PETSC_EXTERN PetscErrorCode PetscFreeA(int,int,const char *,const char *,void *,...);
1256: PETSC_EXTERN PetscErrorCode (*PetscTrMalloc)(size_t,PetscBool,int,const char[],const char[],void**);
1257: PETSC_EXTERN PetscErrorCode (*PetscTrFree)(void*,int,const char[],const char[]);
1258: PETSC_EXTERN PetscErrorCode (*PetscTrRealloc)(size_t,int,const char[],const char[],void**);
1259: PETSC_EXTERN PetscErrorCode PetscMallocSetCoalesce(PetscBool);
1260: PETSC_EXTERN PetscErrorCode PetscMallocSet(PetscErrorCode (*)(size_t,PetscBool,int,const char[],const char[],void**),PetscErrorCode (*)(void*,int,const char[],const char[]),PetscErrorCode (*)(size_t,int,const char[],const char[], void **));
1261: PETSC_EXTERN PetscErrorCode PetscMallocClear(void);
1263: /*
1264: Unlike PetscMallocSet and PetscMallocClear which overwrite the existing settings, these two functions save the previous choice of allocator, and should be used in pair.
1265: */
1266: PETSC_EXTERN PetscErrorCode PetscMallocSetDRAM(void);
1267: PETSC_EXTERN PetscErrorCode PetscMallocResetDRAM(void);
1268: #if defined(PETSC_HAVE_CUDA)
1269: PETSC_EXTERN PetscErrorCode PetscMallocSetCUDAHost(void);
1270: PETSC_EXTERN PetscErrorCode PetscMallocResetCUDAHost(void);
1271: #endif
1272: #if defined(PETSC_HAVE_HIP)
1273: PETSC_EXTERN PetscErrorCode PetscMallocSetHIPHost(void);
1274: PETSC_EXTERN PetscErrorCode PetscMallocResetHIPHost(void);
1275: #endif
1278: #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE
1279: #define MPIU_2PETSCLOGDOUBLE MPI_2DOUBLE_PRECISION
1281: /*
1282: Routines for tracing memory corruption/bleeding with default PETSc memory allocation
1283: */
1284: PETSC_EXTERN PetscErrorCode PetscMallocDump(FILE *);
1285: PETSC_EXTERN PetscErrorCode PetscMallocView(FILE *);
1286: PETSC_EXTERN PetscErrorCode PetscMallocGetCurrentUsage(PetscLogDouble *);
1287: PETSC_EXTERN PetscErrorCode PetscMallocGetMaximumUsage(PetscLogDouble *);
1288: PETSC_EXTERN PetscErrorCode PetscMallocPushMaximumUsage(int);
1289: PETSC_EXTERN PetscErrorCode PetscMallocPopMaximumUsage(int,PetscLogDouble*);
1290: PETSC_EXTERN PetscErrorCode PetscMallocSetDebug(PetscBool,PetscBool);
1291: PETSC_EXTERN PetscErrorCode PetscMallocGetDebug(PetscBool*,PetscBool*,PetscBool*);
1292: PETSC_EXTERN PetscErrorCode PetscMallocValidate(int,const char[],const char[]);
1293: PETSC_EXTERN PetscErrorCode PetscMallocViewSet(PetscLogDouble);
1294: PETSC_EXTERN PetscErrorCode PetscMallocViewGet(PetscBool*);
1295: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeSet(PetscBool);
1296: PETSC_EXTERN PetscErrorCode PetscMallocLogRequestedSizeGet(PetscBool*);
1298: PETSC_EXTERN const char *const PetscDataTypes[];
1299: PETSC_EXTERN PetscErrorCode PetscDataTypeToMPIDataType(PetscDataType,MPI_Datatype*);
1300: PETSC_EXTERN PetscErrorCode PetscMPIDataTypeToPetscDataType(MPI_Datatype,PetscDataType*);
1301: PETSC_EXTERN PetscErrorCode PetscDataTypeGetSize(PetscDataType,size_t*);
1302: PETSC_EXTERN PetscErrorCode PetscDataTypeFromString(const char*,PetscDataType*,PetscBool*);
1304: /*
1305: Basic memory and string operations. These are usually simple wrappers
1306: around the basic Unix system calls, but a few of them have additional
1307: functionality and/or error checking.
1308: */
1309: PETSC_EXTERN PetscErrorCode PetscMemcmp(const void*,const void*,size_t,PetscBool *);
1310: PETSC_EXTERN PetscErrorCode PetscStrlen(const char[],size_t*);
1311: PETSC_EXTERN PetscErrorCode PetscStrToArray(const char[],char,int*,char ***);
1312: PETSC_EXTERN PetscErrorCode PetscStrToArrayDestroy(int,char **);
1313: PETSC_EXTERN PetscErrorCode PetscStrcmp(const char[],const char[],PetscBool *);
1314: PETSC_EXTERN PetscErrorCode PetscStrgrt(const char[],const char[],PetscBool *);
1315: PETSC_EXTERN PetscErrorCode PetscStrcasecmp(const char[],const char[],PetscBool *);
1316: PETSC_EXTERN PetscErrorCode PetscStrncmp(const char[],const char[],size_t,PetscBool *);
1317: PETSC_EXTERN PetscErrorCode PetscStrcpy(char[],const char[]);
1318: PETSC_EXTERN PetscErrorCode PetscStrcat(char[],const char[]);
1319: PETSC_EXTERN PetscErrorCode PetscStrlcat(char[],const char[],size_t);
1320: PETSC_EXTERN PetscErrorCode PetscStrncpy(char[],const char[],size_t);
1321: PETSC_EXTERN PetscErrorCode PetscStrchr(const char[],char,char *[]);
1322: PETSC_EXTERN PetscErrorCode PetscStrtolower(char[]);
1323: PETSC_EXTERN PetscErrorCode PetscStrtoupper(char[]);
1324: PETSC_EXTERN PetscErrorCode PetscStrrchr(const char[],char,char *[]);
1325: PETSC_EXTERN PetscErrorCode PetscStrstr(const char[],const char[],char *[]);
1326: PETSC_EXTERN PetscErrorCode PetscStrrstr(const char[],const char[],char *[]);
1327: PETSC_EXTERN PetscErrorCode PetscStrendswith(const char[],const char[],PetscBool*);
1328: PETSC_EXTERN PetscErrorCode PetscStrbeginswith(const char[],const char[],PetscBool*);
1329: PETSC_EXTERN PetscErrorCode PetscStrendswithwhich(const char[],const char *const*,PetscInt*);
1330: PETSC_EXTERN PetscErrorCode PetscStrallocpy(const char[],char *[]);
1331: PETSC_EXTERN PetscErrorCode PetscStrArrayallocpy(const char *const*,char***);
1332: PETSC_EXTERN PetscErrorCode PetscStrArrayDestroy(char***);
1333: PETSC_EXTERN PetscErrorCode PetscStrNArrayallocpy(PetscInt,const char *const*,char***);
1334: PETSC_EXTERN PetscErrorCode PetscStrNArrayDestroy(PetscInt,char***);
1335: PETSC_EXTERN PetscErrorCode PetscStrreplace(MPI_Comm,const char[],char[],size_t);
1337: PETSC_EXTERN void PetscStrcmpNoError(const char[],const char[],PetscBool *);
1339: PETSC_EXTERN PetscErrorCode PetscTokenCreate(const char[],const char,PetscToken*);
1340: PETSC_EXTERN PetscErrorCode PetscTokenFind(PetscToken,char *[]);
1341: PETSC_EXTERN PetscErrorCode PetscTokenDestroy(PetscToken*);
1343: PETSC_EXTERN PetscErrorCode PetscStrInList(const char[],const char[],char,PetscBool*);
1344: PETSC_EXTERN PetscErrorCode PetscEListFind(PetscInt,const char *const*,const char*,PetscInt*,PetscBool*);
1345: PETSC_EXTERN PetscErrorCode PetscEnumFind(const char *const*,const char*,PetscEnum*,PetscBool*);
1347: /*
1348: These are MPI operations for MPI_Allreduce() etc
1349: */
1350: PETSC_EXTERN MPI_Op MPIU_MAXSUM_OP;
1351: #if (defined(PETSC_HAVE_COMPLEX) && !defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)) || defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1352: PETSC_EXTERN MPI_Op MPIU_SUM;
1353: #else
1354: #define MPIU_SUM MPI_SUM
1355: #endif
1356: #if defined(PETSC_USE_REAL___FLOAT128) || defined(PETSC_USE_REAL___FP16)
1357: PETSC_EXTERN MPI_Op MPIU_MAX;
1358: PETSC_EXTERN MPI_Op MPIU_MIN;
1359: #else
1360: #define MPIU_MAX MPI_MAX
1361: #define MPIU_MIN MPI_MIN
1362: #endif
1363: PETSC_EXTERN PetscErrorCode PetscMaxSum(MPI_Comm,const PetscInt[],PetscInt*,PetscInt*);
1365: PETSC_EXTERN PetscErrorCode MPIULong_Send(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1366: PETSC_EXTERN PetscErrorCode MPIULong_Recv(void*,PetscInt,MPI_Datatype,PetscMPIInt,PetscMPIInt,MPI_Comm);
1368: PETSC_EXTERN const char *const PetscFileModes[];
1370: /*
1371: Defines PETSc error handling.
1372: */
1373: #include <petscerror.h>
1375: #define PETSC_SMALLEST_CLASSID 1211211
1376: PETSC_EXTERN PetscClassId PETSC_LARGEST_CLASSID;
1377: PETSC_EXTERN PetscClassId PETSC_OBJECT_CLASSID;
1378: PETSC_EXTERN PetscErrorCode PetscClassIdRegister(const char[],PetscClassId *);
1379: PETSC_EXTERN PetscErrorCode PetscObjectGetId(PetscObject,PetscObjectId*);
1380: PETSC_EXTERN PetscErrorCode PetscObjectCompareId(PetscObject,PetscObjectId,PetscBool*);
1383: /*
1384: Routines that get memory usage information from the OS
1385: */
1386: PETSC_EXTERN PetscErrorCode PetscMemoryGetCurrentUsage(PetscLogDouble *);
1387: PETSC_EXTERN PetscErrorCode PetscMemoryGetMaximumUsage(PetscLogDouble *);
1388: PETSC_EXTERN PetscErrorCode PetscMemorySetGetMaximumUsage(void);
1389: PETSC_EXTERN PetscErrorCode PetscMemoryTrace(const char[]);
1391: PETSC_EXTERN PetscErrorCode PetscSleep(PetscReal);
1393: /*
1394: Initialization of PETSc
1395: */
1396: PETSC_EXTERN PetscErrorCode PetscInitialize(int*,char***,const char[],const char[]);
1397: PETSC_EXTERN PetscErrorCode PetscInitializeNoPointers(int,char**,const char[],const char[]);
1398: PETSC_EXTERN PetscErrorCode PetscInitializeNoArguments(void);
1399: PETSC_EXTERN PetscErrorCode PetscInitialized(PetscBool *);
1400: PETSC_EXTERN PetscErrorCode PetscFinalized(PetscBool *);
1401: PETSC_EXTERN PetscErrorCode PetscFinalize(void);
1402: PETSC_EXTERN PetscErrorCode PetscInitializeFortran(void);
1403: PETSC_EXTERN PetscErrorCode PetscGetArgs(int*,char ***);
1404: PETSC_EXTERN PetscErrorCode PetscGetArguments(char ***);
1405: PETSC_EXTERN PetscErrorCode PetscFreeArguments(char **);
1407: PETSC_EXTERN PetscErrorCode PetscEnd(void);
1408: PETSC_EXTERN PetscErrorCode PetscSysInitializePackage(void);
1410: PETSC_EXTERN PetscErrorCode PetscPythonInitialize(const char[],const char[]);
1411: PETSC_EXTERN PetscErrorCode PetscPythonFinalize(void);
1412: PETSC_EXTERN PetscErrorCode PetscPythonPrintError(void);
1413: PETSC_EXTERN PetscErrorCode PetscPythonMonitorSet(PetscObject,const char[]);
1415: PETSC_EXTERN PetscErrorCode PetscMonitorCompare(PetscErrorCode (*)(void),void *,PetscErrorCode (*)(void**),PetscErrorCode (*)(void),void *,PetscErrorCode (*)(void**),PetscBool *);
1418: /*
1419: These are so that in extern C code we can caste function pointers to non-extern C
1420: function pointers. Since the regular C++ code expects its function pointers to be C++
1421: */
1422: PETSC_EXTERN_TYPEDEF typedef void (**PetscVoidStarFunction)(void);
1423: PETSC_EXTERN_TYPEDEF typedef void (*PetscVoidFunction)(void);
1424: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*PetscErrorCodeFunction)(void);
1426: /*
1427: Functions that can act on any PETSc object.
1428: */
1429: PETSC_EXTERN PetscErrorCode PetscObjectDestroy(PetscObject*);
1430: PETSC_EXTERN PetscErrorCode PetscObjectGetComm(PetscObject,MPI_Comm *);
1431: PETSC_EXTERN PetscErrorCode PetscObjectGetClassId(PetscObject,PetscClassId *);
1432: PETSC_EXTERN PetscErrorCode PetscObjectGetClassName(PetscObject,const char *[]);
1433: PETSC_EXTERN PetscErrorCode PetscObjectSetType(PetscObject,const char []);
1434: PETSC_EXTERN PetscErrorCode PetscObjectGetType(PetscObject,const char *[]);
1435: PETSC_EXTERN PetscErrorCode PetscObjectSetName(PetscObject,const char[]);
1436: PETSC_EXTERN PetscErrorCode PetscObjectGetName(PetscObject,const char*[]);
1437: PETSC_EXTERN PetscErrorCode PetscObjectSetTabLevel(PetscObject,PetscInt);
1438: PETSC_EXTERN PetscErrorCode PetscObjectGetTabLevel(PetscObject,PetscInt*);
1439: PETSC_EXTERN PetscErrorCode PetscObjectIncrementTabLevel(PetscObject,PetscObject,PetscInt);
1440: PETSC_EXTERN PetscErrorCode PetscObjectReference(PetscObject);
1441: PETSC_EXTERN PetscErrorCode PetscObjectGetReference(PetscObject,PetscInt*);
1442: PETSC_EXTERN PetscErrorCode PetscObjectDereference(PetscObject);
1443: PETSC_EXTERN PetscErrorCode PetscObjectGetNewTag(PetscObject,PetscMPIInt*);
1444: PETSC_EXTERN PetscErrorCode PetscObjectCompose(PetscObject,const char[],PetscObject);
1445: PETSC_EXTERN PetscErrorCode PetscObjectRemoveReference(PetscObject,const char[]);
1446: PETSC_EXTERN PetscErrorCode PetscObjectQuery(PetscObject,const char[],PetscObject*);
1447: PETSC_EXTERN PetscErrorCode PetscObjectComposeFunction_Private(PetscObject,const char[],void (*)(void));
1448: #define PetscObjectComposeFunction(a,b,d) PetscObjectComposeFunction_Private(a,b,(PetscVoidFunction)(d))
1449: PETSC_EXTERN PetscErrorCode PetscObjectSetFromOptions(PetscObject);
1450: PETSC_EXTERN PetscErrorCode PetscObjectSetUp(PetscObject);
1451: PETSC_EXTERN PetscErrorCode PetscObjectSetPrintedOptions(PetscObject);
1452: PETSC_EXTERN PetscErrorCode PetscObjectInheritPrintedOptions(PetscObject,PetscObject);
1453: PETSC_EXTERN PetscErrorCode PetscCommGetNewTag(MPI_Comm,PetscMPIInt*);
1455: #include <petscviewertypes.h>
1456: #include <petscoptions.h>
1458: PETSC_EXTERN PetscErrorCode PetscMallocTraceSet(PetscViewer,PetscBool,PetscLogDouble);
1459: PETSC_EXTERN PetscErrorCode PetscMallocTraceGet(PetscBool*);
1461: PETSC_EXTERN PetscErrorCode PetscObjectsListGetGlobalNumbering(MPI_Comm,PetscInt,PetscObject*,PetscInt*,PetscInt*);
1463: PETSC_EXTERN PetscErrorCode PetscMemoryView(PetscViewer,const char[]);
1464: PETSC_EXTERN PetscErrorCode PetscObjectPrintClassNamePrefixType(PetscObject,PetscViewer);
1465: PETSC_EXTERN PetscErrorCode PetscObjectView(PetscObject,PetscViewer);
1466: #define PetscObjectQueryFunction(obj,name,fptr) PetscObjectQueryFunction_Private((obj),(name),(PetscVoidFunction*)(fptr))
1467: PETSC_EXTERN PetscErrorCode PetscObjectQueryFunction_Private(PetscObject,const char[],void (**)(void));
1468: PETSC_EXTERN PetscErrorCode PetscObjectSetOptionsPrefix(PetscObject,const char[]);
1469: PETSC_EXTERN PetscErrorCode PetscObjectAppendOptionsPrefix(PetscObject,const char[]);
1470: PETSC_EXTERN PetscErrorCode PetscObjectPrependOptionsPrefix(PetscObject,const char[]);
1471: PETSC_EXTERN PetscErrorCode PetscObjectGetOptionsPrefix(PetscObject,const char*[]);
1472: PETSC_EXTERN PetscErrorCode PetscObjectChangeTypeName(PetscObject,const char[]);
1473: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroy(PetscObject);
1474: PETSC_EXTERN PetscErrorCode PetscObjectRegisterDestroyAll(void);
1475: PETSC_EXTERN PetscErrorCode PetscObjectViewFromOptions(PetscObject,PetscObject,const char[]);
1476: PETSC_EXTERN PetscErrorCode PetscObjectName(PetscObject);
1477: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompare(PetscObject,const char[],PetscBool *);
1478: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompare(PetscObject,const char[],PetscBool *);
1479: PETSC_EXTERN PetscErrorCode PetscObjectTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1480: PETSC_EXTERN PetscErrorCode PetscObjectBaseTypeCompareAny(PetscObject,PetscBool*,const char[],...);
1481: PETSC_EXTERN PetscErrorCode PetscRegisterFinalize(PetscErrorCode (*)(void));
1482: PETSC_EXTERN PetscErrorCode PetscRegisterFinalizeAll(void);
1484: #if defined(PETSC_HAVE_SAWS)
1485: PETSC_EXTERN PetscErrorCode PetscSAWsBlock(void);
1486: PETSC_EXTERN PetscErrorCode PetscObjectSAWsViewOff(PetscObject);
1487: PETSC_EXTERN PetscErrorCode PetscObjectSAWsSetBlock(PetscObject,PetscBool);
1488: PETSC_EXTERN PetscErrorCode PetscObjectSAWsBlock(PetscObject);
1489: PETSC_EXTERN PetscErrorCode PetscObjectSAWsGrantAccess(PetscObject);
1490: PETSC_EXTERN PetscErrorCode PetscObjectSAWsTakeAccess(PetscObject);
1491: PETSC_EXTERN void PetscStackSAWsGrantAccess(void);
1492: PETSC_EXTERN void PetscStackSAWsTakeAccess(void);
1493: PETSC_EXTERN PetscErrorCode PetscStackViewSAWs(void);
1494: PETSC_EXTERN PetscErrorCode PetscStackSAWsViewOff(void);
1496: #else
1497: #define PetscSAWsBlock() 0
1498: #define PetscObjectSAWsViewOff(obj) 0
1499: #define PetscObjectSAWsSetBlock(obj,flg) 0
1500: #define PetscObjectSAWsBlock(obj) 0
1501: #define PetscObjectSAWsGrantAccess(obj) 0
1502: #define PetscObjectSAWsTakeAccess(obj) 0
1503: #define PetscStackViewSAWs() 0
1504: #define PetscStackSAWsViewOff() 0
1505: #define PetscStackSAWsTakeAccess()
1506: #define PetscStackSAWsGrantAccess()
1508: #endif
1510: PETSC_EXTERN PetscErrorCode PetscDLOpen(const char[],PetscDLMode,PetscDLHandle *);
1511: PETSC_EXTERN PetscErrorCode PetscDLClose(PetscDLHandle *);
1512: PETSC_EXTERN PetscErrorCode PetscDLSym(PetscDLHandle,const char[],void **);
1514: #if defined(PETSC_USE_DEBUG)
1515: PETSC_EXTERN PetscErrorCode PetscMallocGetStack(void*,PetscStack**);
1516: #endif
1517: PETSC_EXTERN PetscErrorCode PetscObjectsDump(FILE*,PetscBool);
1519: PETSC_EXTERN PetscErrorCode PetscObjectListDestroy(PetscObjectList*);
1520: PETSC_EXTERN PetscErrorCode PetscObjectListFind(PetscObjectList,const char[],PetscObject*);
1521: PETSC_EXTERN PetscErrorCode PetscObjectListReverseFind(PetscObjectList,PetscObject,char**,PetscBool*);
1522: PETSC_EXTERN PetscErrorCode PetscObjectListAdd(PetscObjectList *,const char[],PetscObject);
1523: PETSC_EXTERN PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *,const char[]);
1524: PETSC_EXTERN PetscErrorCode PetscObjectListDuplicate(PetscObjectList,PetscObjectList *);
1526: /*
1527: Dynamic library lists. Lists of names of routines in objects or in dynamic
1528: link libraries that will be loaded as needed.
1529: */
1531: #define PetscFunctionListAdd(list,name,fptr) PetscFunctionListAdd_Private((list),(name),(PetscVoidFunction)(fptr))
1532: PETSC_EXTERN PetscErrorCode PetscFunctionListAdd_Private(PetscFunctionList*,const char[],void (*)(void));
1533: PETSC_EXTERN PetscErrorCode PetscFunctionListDestroy(PetscFunctionList*);
1534: #define PetscFunctionListFind(list,name,fptr) PetscFunctionListFind_Private((list),(name),(PetscVoidFunction*)(fptr))
1535: PETSC_EXTERN PetscErrorCode PetscFunctionListFind_Private(PetscFunctionList,const char[],void (**)(void));
1536: PETSC_EXTERN PetscErrorCode PetscFunctionListPrintTypes(MPI_Comm,FILE*,const char[],const char[],const char[],const char[],PetscFunctionList,const char[],const char[]);
1537: PETSC_EXTERN PetscErrorCode PetscFunctionListDuplicate(PetscFunctionList,PetscFunctionList *);
1538: PETSC_EXTERN PetscErrorCode PetscFunctionListView(PetscFunctionList,PetscViewer);
1539: PETSC_EXTERN PetscErrorCode PetscFunctionListGet(PetscFunctionList,const char ***,int*);
1541: PETSC_EXTERN PetscDLLibrary PetscDLLibrariesLoaded;
1542: PETSC_EXTERN PetscErrorCode PetscDLLibraryAppend(MPI_Comm,PetscDLLibrary *,const char[]);
1543: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrepend(MPI_Comm,PetscDLLibrary *,const char[]);
1544: PETSC_EXTERN PetscErrorCode PetscDLLibrarySym(MPI_Comm,PetscDLLibrary *,const char[],const char[],void **);
1545: PETSC_EXTERN PetscErrorCode PetscDLLibraryPrintPath(PetscDLLibrary);
1546: PETSC_EXTERN PetscErrorCode PetscDLLibraryRetrieve(MPI_Comm,const char[],char *,size_t,PetscBool *);
1547: PETSC_EXTERN PetscErrorCode PetscDLLibraryOpen(MPI_Comm,const char[],PetscDLLibrary *);
1548: PETSC_EXTERN PetscErrorCode PetscDLLibraryClose(PetscDLLibrary);
1550: /*
1551: Useful utility routines
1552: */
1553: PETSC_EXTERN PetscErrorCode PetscSplitOwnership(MPI_Comm,PetscInt*,PetscInt*);
1554: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm,PetscInt,PetscInt*,PetscInt*);
1555: PETSC_EXTERN PetscErrorCode PetscSplitOwnershipEqual(MPI_Comm,PetscInt*,PetscInt*);
1556: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseBegin(MPI_Comm,PetscMPIInt);
1557: PETSC_EXTERN PetscErrorCode PetscSequentialPhaseEnd(MPI_Comm,PetscMPIInt);
1558: PETSC_EXTERN PetscErrorCode PetscBarrier(PetscObject);
1559: PETSC_EXTERN PetscErrorCode PetscMPIDump(FILE*);
1560: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxInt(MPI_Comm,PetscInt[],PetscInt[]);
1561: PETSC_EXTERN PetscErrorCode PetscGlobalMinMaxReal(MPI_Comm,PetscReal[],PetscReal[]);
1563: /*MC
1564: PetscNot - negates a logical type value and returns result as a PetscBool
1566: Notes:
1567: This is useful in cases like
1568: $ int *a;
1569: $ PetscBool flag = PetscNot(a)
1570: where !a would not return a PetscBool because we cannot provide a cast from int to PetscBool in C.
1572: Level: beginner
1574: .seealso : PetscBool, PETSC_TRUE, PETSC_FALSE
1575: M*/
1576: #define PetscNot(a) ((a) ? PETSC_FALSE : PETSC_TRUE)
1578: /*MC
1579: PetscHelpPrintf - Prints help messages.
1581: Synopsis:
1582: #include <petscsys.h>
1583: PetscErrorCode (*PetscHelpPrintf)(MPI_Comm comm, const char format[],args);
1585: Collective on comm
1587: Input Parameters:
1588: + comm - the MPI communicator over which the help message is printed
1589: . format - the usual printf() format string
1590: - args - arguments to be printed
1592: Level: developer
1594: Fortran Note:
1595: This routine is not supported in Fortran.
1597: Note:
1598: You can change how help messages are printed by replacing the function pointer with a function that does not simply write to stdout.
1600: To use, write your own function, for example,
1601: $PetscErrorCode mypetschelpprintf(MPI_Comm comm,const char format[],....)
1602: ${
1603: $ return(0);
1604: $}
1605: then do the assigment
1606: $ PetscHelpPrintf = mypetschelpprintf;
1607: You can do the assignment before PetscInitialize().
1609: The default routine used is called PetscHelpPrintfDefault().
1611: .seealso: PetscFPrintf(), PetscSynchronizedPrintf(), PetscErrorPrintf()
1612: M*/
1613: PETSC_EXTERN PetscErrorCode (*PetscHelpPrintf)(MPI_Comm,const char[],...);
1615: /*
1616: Defines PETSc profiling.
1617: */
1618: #include <petsclog.h>
1620: /*
1621: Simple PETSc parallel IO for ASCII printing
1622: */
1623: PETSC_EXTERN PetscErrorCode PetscFixFilename(const char[],char[]);
1624: PETSC_EXTERN PetscErrorCode PetscFOpen(MPI_Comm,const char[],const char[],FILE**);
1625: PETSC_EXTERN PetscErrorCode PetscFClose(MPI_Comm,FILE*);
1626: PETSC_EXTERN PetscErrorCode PetscFPrintf(MPI_Comm,FILE*,const char[],...);
1627: PETSC_EXTERN PetscErrorCode PetscPrintf(MPI_Comm,const char[],...);
1628: PETSC_EXTERN PetscErrorCode PetscSNPrintf(char*,size_t,const char [],...);
1629: PETSC_EXTERN PetscErrorCode PetscSNPrintfCount(char*,size_t,const char [],size_t*,...);
1630: PETSC_EXTERN PetscErrorCode PetscFormatRealArray(char[],size_t,const char*,PetscInt,const PetscReal[]);
1632: PETSC_EXTERN PetscErrorCode PetscErrorPrintfDefault(const char [],...);
1633: PETSC_EXTERN PetscErrorCode PetscErrorPrintfNone(const char [],...);
1634: PETSC_EXTERN PetscErrorCode PetscHelpPrintfDefault(MPI_Comm,const char [],...);
1636: PETSC_EXTERN PetscErrorCode PetscFormatConvertGetSize(const char*,size_t*);
1637: PETSC_EXTERN PetscErrorCode PetscFormatConvert(const char*,char *);
1639: #if defined(PETSC_HAVE_POPEN)
1640: PETSC_EXTERN PetscErrorCode PetscPOpen(MPI_Comm,const char[],const char[],const char[],FILE **);
1641: PETSC_EXTERN PetscErrorCode PetscPClose(MPI_Comm,FILE*);
1642: PETSC_EXTERN PetscErrorCode PetscPOpenSetMachine(const char[]);
1643: #endif
1645: PETSC_EXTERN PetscErrorCode PetscSynchronizedPrintf(MPI_Comm,const char[],...);
1646: PETSC_EXTERN PetscErrorCode PetscSynchronizedFPrintf(MPI_Comm,FILE*,const char[],...);
1647: PETSC_EXTERN PetscErrorCode PetscSynchronizedFlush(MPI_Comm,FILE*);
1648: PETSC_EXTERN PetscErrorCode PetscSynchronizedFGets(MPI_Comm,FILE*,size_t,char[]);
1649: PETSC_EXTERN PetscErrorCode PetscStartMatlab(MPI_Comm,const char[],const char[],FILE**);
1650: PETSC_EXTERN PetscErrorCode PetscStartJava(MPI_Comm,const char[],const char[],FILE**);
1651: PETSC_EXTERN PetscErrorCode PetscGetPetscDir(const char*[]);
1653: PETSC_EXTERN PetscClassId PETSC_CONTAINER_CLASSID;
1654: PETSC_EXTERN PetscErrorCode PetscContainerGetPointer(PetscContainer,void**);
1655: PETSC_EXTERN PetscErrorCode PetscContainerSetPointer(PetscContainer,void*);
1656: PETSC_EXTERN PetscErrorCode PetscContainerDestroy(PetscContainer*);
1657: PETSC_EXTERN PetscErrorCode PetscContainerCreate(MPI_Comm,PetscContainer*);
1658: PETSC_EXTERN PetscErrorCode PetscContainerSetUserDestroy(PetscContainer,PetscErrorCode (*)(void*));
1659: PETSC_EXTERN PetscErrorCode PetscContainerUserDestroyDefault(void*);
1661: /*
1662: For use in debuggers
1663: */
1664: PETSC_EXTERN PetscMPIInt PetscGlobalRank;
1665: PETSC_EXTERN PetscMPIInt PetscGlobalSize;
1666: PETSC_EXTERN PetscErrorCode PetscIntView(PetscInt,const PetscInt[],PetscViewer);
1667: PETSC_EXTERN PetscErrorCode PetscRealView(PetscInt,const PetscReal[],PetscViewer);
1668: PETSC_EXTERN PetscErrorCode PetscScalarView(PetscInt,const PetscScalar[],PetscViewer);
1670: #include <stddef.h>
1671: #include <string.h> /* for memcpy, memset */
1672: #include <stdlib.h>
1674: #if defined(PETSC_HAVE_XMMINTRIN_H) && !defined(__CUDACC__)
1675: #include <xmmintrin.h>
1676: #endif
1678: /*@C
1679: PetscMemmove - Copies n bytes, beginning at location b, to the space
1680: beginning at location a. Copying between regions that overlap will
1681: take place correctly. Use PetscMemcpy() if the locations do not overlap
1683: Not Collective
1685: Input Parameters:
1686: + b - pointer to initial memory space
1687: . a - pointer to copy space
1688: - n - length (in bytes) of space to copy
1690: Level: intermediate
1692: Note:
1693: PetscArraymove() is preferred
1694: This routine is analogous to memmove().
1696: Developers Note: This is inlined for performance
1698: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscStrallocpy(),
1699: PetscArraymove()
1700: @*/
1701: PETSC_STATIC_INLINE PetscErrorCode PetscMemmove(void *a,const void *b,size_t n)
1702: {
1704: if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to null pointer");
1705: if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1706: #if !defined(PETSC_HAVE_MEMMOVE)
1707: if (a < b) {
1708: if (a <= b - n) memcpy(a,b,n);
1709: else {
1710: memcpy(a,b,(int)(b - a));
1711: PetscMemmove(b,b + (int)(b - a),n - (int)(b - a));
1712: }
1713: } else {
1714: if (b <= a - n) memcpy(a,b,n);
1715: else {
1716: memcpy(b + n,b + (n - (int)(a - b)),(int)(a - b));
1717: PetscMemmove(a,b,n - (int)(a - b));
1718: }
1719: }
1720: #else
1721: memmove((char*)(a),(char*)(b),n);
1722: #endif
1723: return(0);
1724: }
1726: /*@C
1727: PetscMemcpy - Copies n bytes, beginning at location b, to the space
1728: beginning at location a. The two memory regions CANNOT overlap, use
1729: PetscMemmove() in that case.
1731: Not Collective
1733: Input Parameters:
1734: + b - pointer to initial memory space
1735: - n - length (in bytes) of space to copy
1737: Output Parameter:
1738: . a - pointer to copy space
1740: Level: intermediate
1742: Compile Option:
1743: PETSC_PREFER_DCOPY_FOR_MEMCPY will cause the BLAS dcopy() routine to be used
1744: for memory copies on double precision values.
1745: PETSC_PREFER_COPY_FOR_MEMCPY will cause C code to be used
1746: for memory copies on double precision values.
1747: PETSC_PREFER_FORTRAN_FORMEMCPY will cause Fortran code to be used
1748: for memory copies on double precision values.
1750: Note:
1751: Prefer PetscArraycpy()
1752: This routine is analogous to memcpy().
1753: Not available from Fortran
1755: Developer Note: this is inlined for fastest performance
1757: .seealso: PetscMemzero(), PetscMemcmp(), PetscArrayzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy()
1759: @*/
1760: PETSC_STATIC_INLINE PetscErrorCode PetscMemcpy(void *a,const void *b,size_t n)
1761: {
1762: #if defined(PETSC_USE_DEBUG)
1763: size_t al = (size_t) a,bl = (size_t) b;
1764: size_t nl = (size_t) n;
1766: if (n > 0 && !b) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy from a null pointer");
1767: if (n > 0 && !a) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to copy to a null pointer");
1768: #else
1770: #endif
1771: if (a != b && n > 0) {
1772: #if defined(PETSC_USE_DEBUG)
1773: if ((al > bl && (al - bl) < nl) || (bl - al) < nl) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Memory regions overlap: either use PetscMemmov()\n\
1774: or make sure your copy regions and lengths are correct. \n\
1775: Length (bytes) %ld first address %ld second address %ld",nl,al,bl);
1776: #endif
1777: #if (defined(PETSC_PREFER_DCOPY_FOR_MEMCPY) || defined(PETSC_PREFER_COPY_FOR_MEMCPY) || defined(PETSC_PREFER_FORTRAN_FORMEMCPY))
1778: if (!(a % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1779: size_t len = n/sizeof(PetscScalar);
1780: #if defined(PETSC_PREFER_DCOPY_FOR_MEMCPY)
1781: PetscBLASInt one = 1,blen;
1783: PetscBLASIntCast(len,&blen);
1784: PetscStackCallBLAS("BLAScopy",BLAScopy_(&blen,(PetscScalar *)b,&one,(PetscScalar *)a,&one));
1785: #elif defined(PETSC_PREFER_FORTRAN_FORMEMCPY)
1786: fortrancopy_(&len,(PetscScalar*)b,(PetscScalar*)a);
1787: #else
1788: size_t i;
1789: PetscScalar *x = (PetscScalar*)b, *y = (PetscScalar*)a;
1790: for (i=0; i<len; i++) y[i] = x[i];
1791: #endif
1792: } else {
1793: memcpy((char*)(a),(char*)(b),n);
1794: }
1795: #else
1796: memcpy((char*)(a),(char*)(b),n);
1797: #endif
1798: }
1799: return(0);
1800: }
1802: /*@C
1803: PetscMemzero - Zeros the specified memory.
1805: Not Collective
1807: Input Parameters:
1808: + a - pointer to beginning memory location
1809: - n - length (in bytes) of memory to initialize
1811: Level: intermediate
1813: Compile Option:
1814: PETSC_PREFER_BZERO - on certain machines (the IBM RS6000) the bzero() routine happens
1815: to be faster than the memset() routine. This flag causes the bzero() routine to be used.
1817: Not available from Fortran
1818: Prefer PetscArrayzero()
1820: Developer Note: this is inlined for fastest performance
1822: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy()
1823: @*/
1824: PETSC_STATIC_INLINE PetscErrorCode PetscMemzero(void *a,size_t n)
1825: {
1826: if (n > 0) {
1827: #if defined(PETSC_USE_DEBUG)
1828: if (!a) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Trying to zero at a null pointer with %zu bytes",n);
1829: #endif
1830: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO)
1831: if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1832: size_t i,len = n/sizeof(PetscScalar);
1833: PetscScalar *x = (PetscScalar*)a;
1834: for (i=0; i<len; i++) x[i] = 0.0;
1835: } else {
1836: #elif defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1837: if (!(((long) a) % sizeof(PetscScalar)) && !(n % sizeof(PetscScalar))) {
1838: PetscInt len = n/sizeof(PetscScalar);
1839: fortranzero_(&len,(PetscScalar*)a);
1840: } else {
1841: #endif
1842: #if defined(PETSC_PREFER_BZERO)
1843: bzero((char *)a,n);
1844: #else
1845: memset((char*)a,0,n);
1846: #endif
1847: #if defined(PETSC_PREFER_ZERO_FOR_MEMZERO) || defined(PETSC_PREFER_FORTRAN_FOR_MEMZERO)
1848: }
1849: #endif
1850: }
1851: return 0;
1852: }
1854: /*MC
1855: PetscArraycmp - Compares two arrays in memory.
1857: Synopsis:
1858: #include <petscsys.h>
1859: PetscErrorCode PetscArraycmp(const anytype *str1,const anytype *str2,size_t cnt,PetscBool *e)
1861: Not Collective
1863: Input Parameters:
1864: + str1 - First array
1865: . str2 - Second array
1866: - cnt - Count of the array, not in bytes, but number of entries in the arrays
1868: Output Parameters:
1869: . e - PETSC_TRUE if equal else PETSC_FALSE.
1871: Level: intermediate
1873: Note:
1874: This routine is a preferred replacement to PetscMemcmp()
1875: The arrays must be of the same type
1877: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy(),
1878: PetscArraymove()
1879: M*/
1880: #define PetscArraycmp(str1,str2,cnt,e) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcmp(str1,str2,(size_t)(cnt)*sizeof(*(str1)),e))
1882: /*MC
1883: PetscArraymove - Copies from one array in memory to another, the arrays may overlap. Use PetscArraycpy() when the arrays
1884: do not overlap
1886: Synopsis:
1887: #include <petscsys.h>
1888: PetscErrorCode PetscArraymove(anytype *str1,const anytype *str2,size_t cnt)
1890: Not Collective
1892: Input Parameters:
1893: + str1 - First array
1894: . str2 - Second array
1895: - cnt - Count of the array, not in bytes, but number of entries in the arrays
1897: Level: intermediate
1899: Note:
1900: This routine is a preferred replacement to PetscMemmove()
1901: The arrays must be of the same type
1903: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycpy(), PetscMemmove(), PetscArraycmp(), PetscStrallocpy()
1904: M*/
1905: #define PetscArraymove(str1,str2,cnt) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemmove(str1,str2,(size_t)(cnt)*sizeof(*(str1))))
1907: /*MC
1908: PetscArraycpy - Copies from one array in memory to another
1910: Synopsis:
1911: #include <petscsys.h>
1912: PetscErrorCode PetscArraycpy(anytype *str1,const anytype *str2,size_t cnt)
1914: Not Collective
1916: Input Parameters:
1917: + str1 - First array (destination)
1918: . str2 - Second array (source)
1919: - cnt - Count of the array, not in bytes, but number of entries in the arrays
1921: Level: intermediate
1923: Note:
1924: This routine is a preferred replacement to PetscMemcpy()
1925: The arrays must be of the same type
1927: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraymove(), PetscMemmove(), PetscArraycmp(), PetscStrallocpy()
1928: M*/
1929: #define PetscArraycpy(str1,str2,cnt) ((sizeof(*(str1)) != sizeof(*(str2))) || PetscMemcpy(str1,str2,(size_t)(cnt)*sizeof(*(str1))))
1931: /*MC
1932: PetscArrayzero - Zeros an array in memory.
1934: Synopsis:
1935: #include <petscsys.h>
1936: PetscErrorCode PetscArrayzero(anytype *str1,size_t cnt)
1938: Not Collective
1940: Input Parameters:
1941: + str1 - array
1942: - cnt - Count of the array, not in bytes, but number of entries in the array
1944: Level: intermediate
1946: Note:
1947: This routine is a preferred replacement to PetscMemzero()
1949: .seealso: PetscMemcpy(), PetscMemcmp(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscMemmove(), PetscStrallocpy(), PetscArraymove()
1950: M*/
1951: #define PetscArrayzero(str1,cnt) PetscMemzero(str1,(size_t)(cnt)*sizeof(*(str1)))
1953: /*MC
1954: PetscPrefetchBlock - Prefetches a block of memory
1956: Synopsis:
1957: #include <petscsys.h>
1958: void PetscPrefetchBlock(const anytype *a,size_t n,int rw,int t)
1960: Not Collective
1962: Input Parameters:
1963: + a - pointer to first element to fetch (any type but usually PetscInt or PetscScalar)
1964: . n - number of elements to fetch
1965: . rw - 1 if the memory will be written to, otherwise 0 (ignored by many processors)
1966: - t - temporal locality (PETSC_PREFETCH_HINT_{NTA,T0,T1,T2}), see note
1968: Level: developer
1970: Notes:
1971: The last two arguments (rw and t) must be compile-time constants.
1973: Adopting Intel's x86/x86-64 conventions, there are four levels of temporal locality. Not all architectures offer
1974: equivalent locality hints, but the following macros are always defined to their closest analogue.
1975: + PETSC_PREFETCH_HINT_NTA - Non-temporal. Prefetches directly to L1, evicts to memory (skips higher level cache unless it was already there when prefetched).
1976: . 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.
1977: . PETSC_PREFETCH_HINT_T1 - Fetch to level 2 and higher (not L1).
1978: - PETSC_PREFETCH_HINT_T2 - Fetch to high-level cache only. (On many systems, T0 and T1 are equivalent.)
1980: This function does nothing on architectures that do not support prefetch and never errors (even if passed an invalid
1981: address).
1983: M*/
1984: #define PetscPrefetchBlock(a,n,rw,t) do { \
1985: const char *_p = (const char*)(a),*_end = (const char*)((a)+(n)); \
1986: for (; _p < _end; _p += PETSC_LEVEL1_DCACHE_LINESIZE) PETSC_Prefetch(_p,(rw),(t)); \
1987: } while (0)
1989: /*
1990: Determine if some of the kernel computation routines use
1991: Fortran (rather than C) for the numerical calculations. On some machines
1992: and compilers (like complex numbers) the Fortran version of the routines
1993: is faster than the C/C++ versions. The flag --with-fortran-kernels
1994: should be used with ./configure to turn these on.
1995: */
1996: #if defined(PETSC_USE_FORTRAN_KERNELS)
1998: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTCRL)
1999: #define PETSC_USE_FORTRAN_KERNEL_MULTCRL
2000: #endif
2002: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM)
2003: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJPERM
2004: #endif
2006: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
2007: #define PETSC_USE_FORTRAN_KERNEL_MULTAIJ
2008: #endif
2010: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
2011: #define PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ
2012: #endif
2014: #if !defined(PETSC_USE_FORTRAN_KERNEL_NORM)
2015: #define PETSC_USE_FORTRAN_KERNEL_NORM
2016: #endif
2018: #if !defined(PETSC_USE_FORTRAN_KERNEL_MAXPY)
2019: #define PETSC_USE_FORTRAN_KERNEL_MAXPY
2020: #endif
2022: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
2023: #define PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ
2024: #endif
2026: #if !defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
2027: #define PETSC_USE_FORTRAN_KERNEL_RELAXAIJ
2028: #endif
2030: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
2031: #define PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ
2032: #endif
2034: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
2035: #define PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ
2036: #endif
2038: #if !defined(PETSC_USE_FORTRAN_KERNEL_MDOT)
2039: #define PETSC_USE_FORTRAN_KERNEL_MDOT
2040: #endif
2042: #if !defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
2043: #define PETSC_USE_FORTRAN_KERNEL_XTIMESY
2044: #endif
2046: #if !defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
2047: #define PETSC_USE_FORTRAN_KERNEL_AYPX
2048: #endif
2050: #if !defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
2051: #define PETSC_USE_FORTRAN_KERNEL_WAXPY
2052: #endif
2054: #endif
2056: /*
2057: Macros for indicating code that should be compiled with a C interface,
2058: rather than a C++ interface. Any routines that are dynamically loaded
2059: (such as the PCCreate_XXX() routines) must be wrapped so that the name
2060: mangler does not change the functions symbol name. This just hides the
2061: ugly extern "C" {} wrappers.
2062: */
2063: #if defined(__cplusplus)
2064: # define EXTERN_C_BEGIN extern "C" {
2065: # define EXTERN_C_END }
2066: #else
2067: # define EXTERN_C_BEGIN
2068: # define EXTERN_C_END
2069: #endif
2071: /* --------------------------------------------------------------------*/
2073: /*MC
2074: MPI_Comm - the basic object used by MPI to determine which processes are involved in a
2075: communication
2077: Level: beginner
2079: Note: This manual page is a place-holder because MPICH does not have a manual page for MPI_Comm
2081: .seealso: PETSC_COMM_WORLD, PETSC_COMM_SELF
2082: M*/
2084: #if defined(PETSC_HAVE_MPIIO)
2085: PETSC_EXTERN PetscErrorCode MPIU_File_write_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2086: PETSC_EXTERN PetscErrorCode MPIU_File_read_all(MPI_File,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2087: PETSC_EXTERN PetscErrorCode MPIU_File_write_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2088: PETSC_EXTERN PetscErrorCode MPIU_File_read_at(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2089: PETSC_EXTERN PetscErrorCode MPIU_File_write_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2090: PETSC_EXTERN PetscErrorCode MPIU_File_read_at_all(MPI_File,MPI_Offset,void*,PetscMPIInt,MPI_Datatype,MPI_Status*);
2091: #endif
2093: /* the following petsc_static_inline require petscerror.h */
2095: /* Limit MPI to 32-bits */
2096: #define PETSC_MPI_INT_MAX 2147483647
2097: #define PETSC_MPI_INT_MIN -2147483647
2098: /* Limit BLAS to 32-bits */
2099: #define PETSC_BLAS_INT_MAX 2147483647
2100: #define PETSC_BLAS_INT_MIN -2147483647
2102: /*@C
2103: PetscIntCast - casts a PetscInt64 (which is 64 bits in size) to a PetscInt (which may be 32 bits in size), generates an
2104: error if the PetscInt is not large enough to hold the number.
2106: Not Collective
2108: Input Parameter:
2109: . a - the PetscInt64 value
2111: Output Parameter:
2112: . b - the resulting PetscInt value
2114: Level: advanced
2116: Notes: If integers needed for the applications are too large to fit in 32 bit ints you can ./configure using --with-64-bit-indices to make PetscInt use 64 bit ints
2118: Not available from Fortran
2120: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscBLASIntCast()
2121: @*/
2122: PETSC_STATIC_INLINE PetscErrorCode PetscIntCast(PetscInt64 a,PetscInt *b)
2123: {
2125: #if !defined(PETSC_USE_64BIT_INDICES)
2126: if ((a) > PETSC_MAX_INT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Integer too long for PetscInt, you may need to ./configure using --with-64-bit-indices");
2127: #endif
2128: *b = (PetscInt)(a);
2129: return(0);
2130: }
2132: /*@C
2133: PetscBLASIntCast - casts a PetscInt (which may be 64 bits in size) to a PetscBLASInt (which may be 32 bits in size), generates an
2134: error if the PetscBLASInt is not large enough to hold the number.
2136: Not Collective
2138: Input Parameter:
2139: . a - the PetscInt value
2141: Output Parameter:
2142: . b - the resulting PetscBLASInt value
2144: Level: advanced
2146: Notes:
2147: Not available from Fortran
2148: Errors if the integer is negative since PETSc calls to BLAS/LAPACK never need to cast negative integer inputs
2150: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscMPIIntCast(), PetscIntCast()
2151: @*/
2152: PETSC_STATIC_INLINE PetscErrorCode PetscBLASIntCast(PetscInt a,PetscBLASInt *b)
2153: {
2155: #if defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
2156: *b = 0;
2157: if ((a) > PETSC_BLAS_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for BLAS/LAPACK");
2158: #endif
2159: if (a < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Passing negative integer to BLAS/LAPACK routine");
2160: *b = (PetscBLASInt)(a);
2161: return(0);
2162: }
2164: /*@C
2165: PetscMPIIntCast - casts a PetscInt (which may be 64 bits in size) to a PetscMPIInt (which may be 32 bits in size), generates an
2166: error if the PetscMPIInt is not large enough to hold the number.
2168: Not Collective
2170: Input Parameter:
2171: . a - the PetscInt value
2173: Output Parameter:
2174: . b - the resulting PetscMPIInt value
2176: Level: advanced
2178: Not available from Fortran
2180: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntCast()
2181: @*/
2182: PETSC_STATIC_INLINE PetscErrorCode PetscMPIIntCast(PetscInt a,PetscMPIInt *b)
2183: {
2185: #if defined(PETSC_USE_64BIT_INDICES)
2186: *b = 0;
2187: if ((a) > PETSC_MPI_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for MPI");
2188: #endif
2189: *b = (PetscMPIInt)(a);
2190: return(0);
2191: }
2193: #define PetscInt64Mult(a,b) ((PetscInt64)(a))*((PetscInt64)(b))
2195: /*@C
2197: PetscRealIntMultTruncate - Computes the product of a positive PetscReal and a positive PetscInt and truncates the value to slightly less than the maximal possible value
2199: Not Collective
2201: Input Parameter:
2202: + a - the PetscReal value
2203: - b - the second value
2205: Returns:
2206: the result as a PetscInt value
2208: Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2209: Use PetscIntMultTruncate() to compute the product of two positive PetscInt and truncate to fit a PetscInt
2210: 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
2212: Developers Note:
2213: We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks.
2215: This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2217: Not available from Fortran
2219: Level: advanced
2221: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2222: @*/
2223: PETSC_STATIC_INLINE PetscInt PetscRealIntMultTruncate(PetscReal a,PetscInt b)
2224: {
2225: PetscInt64 r;
2227: r = (PetscInt64) (a*(PetscReal)b);
2228: if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2229: return (PetscInt) r;
2230: }
2232: /*@C
2234: PetscIntMultTruncate - Computes the product of two positive PetscInt and truncates the value to slightly less than the maximal possible value
2236: Not Collective
2238: Input Parameter:
2239: + a - the PetscInt value
2240: - b - the second value
2242: Returns:
2243: the result as a PetscInt value
2245: Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2246: Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2247: 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
2249: Not available from Fortran
2251: Developers Note: We currently assume that PetscInt addition can never overflow, this is obviously wrong but requires many more checks.
2253: This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2255: Level: advanced
2257: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2258: @*/
2259: PETSC_STATIC_INLINE PetscInt PetscIntMultTruncate(PetscInt a,PetscInt b)
2260: {
2261: PetscInt64 r;
2263: r = PetscInt64Mult(a,b);
2264: if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2265: return (PetscInt) r;
2266: }
2268: /*@C
2270: PetscIntSumTruncate - Computes the sum of two positive PetscInt and truncates the value to slightly less than the maximal possible value
2272: Not Collective
2274: Input Parameter:
2275: + a - the PetscInt value
2276: - b - the second value
2278: Returns:
2279: the result as a PetscInt value
2281: Use PetscInt64Mult() to compute the product of two PetscInt as a PetscInt64
2282: Use PetscRealIntMultTruncate() to compute the product of a PetscReal and a PetscInt and truncate to fit a PetscInt
2283: 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
2285: This is used where we compute approximate sizes for workspace and need to insure the workspace is index-able.
2287: Not available from Fortran
2289: Level: advanced
2291: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2292: @*/
2293: PETSC_STATIC_INLINE PetscInt PetscIntSumTruncate(PetscInt a,PetscInt b)
2294: {
2295: PetscInt64 r;
2297: r = ((PetscInt64)a) + ((PetscInt64)b);
2298: if (r > PETSC_MAX_INT - 100) r = PETSC_MAX_INT - 100;
2299: return (PetscInt) r;
2300: }
2302: /*@C
2304: PetscIntMultError - Computes the product of two positive PetscInt and generates an error with overflow.
2306: Not Collective
2308: Input Parameter:
2309: + a - the PetscInt value
2310: - b - the second value
2312: Output Parameter:ma
2313: . result - the result as a PetscInt value, or NULL if you do not want the result, you just want to check if it overflows
2315: Use PetscInt64Mult() to compute the product of two 32 bit PetscInt and store in a PetscInt64
2316: Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt
2318: Not available from Fortran
2320: Developers Note: We currently assume that PetscInt addition does not overflow, this is obviously wrong but requires many more checks.
2322: Level: advanced
2324: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscIntMult64(), PetscIntSumError()
2325: @*/
2326: PETSC_STATIC_INLINE PetscErrorCode PetscIntMultError(PetscInt a,PetscInt b,PetscInt *result)
2327: {
2328: PetscInt64 r;
2331: r = PetscInt64Mult(a,b);
2332: #if !defined(PETSC_USE_64BIT_INDICES)
2333: 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);
2334: #endif
2335: if (result) *result = (PetscInt) r;
2336: return(0);
2337: }
2339: /*@C
2341: PetscIntSumError - Computes the sum of two positive PetscInt and generates an error with overflow.
2343: Not Collective
2345: Input Parameter:
2346: + a - the PetscInt value
2347: - b - the second value
2349: Output Parameter:ma
2350: . c - the result as a PetscInt value, or NULL if you do not want the result, you just want to check if it overflows
2352: Use PetscInt64Mult() to compute the product of two 32 bit PetscInt and store in a PetscInt64
2353: Use PetscIntMultTruncate() to compute the product of two PetscInt and truncate it to fit in a PetscInt
2355: Not available from Fortran
2357: Level: advanced
2359: .seealso: PetscBLASInt, PetscMPIInt, PetscInt, PetscBLASIntCast(), PetscInt64Mult()
2360: @*/
2361: PETSC_STATIC_INLINE PetscErrorCode PetscIntSumError(PetscInt a,PetscInt b,PetscInt *result)
2362: {
2363: PetscInt64 r;
2366: r = ((PetscInt64)a) + ((PetscInt64)b);
2367: #if !defined(PETSC_USE_64BIT_INDICES)
2368: 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);
2369: #endif
2370: if (result) *result = (PetscInt) r;
2371: return(0);
2372: }
2374: /*
2375: The IBM include files define hz, here we hide it so that it may be used as a regular user variable.
2376: */
2377: #if defined(hz)
2378: # undef hz
2379: #endif
2381: #include <limits.h>
2383: /* The number of bits in a byte */
2385: #define PETSC_BITS_PER_BYTE CHAR_BIT
2387: #if defined(PETSC_HAVE_SYS_TYPES_H)
2388: # include <sys/types.h>
2389: #endif
2391: /*MC
2393: PETSC_VERSION - This manual page provides information about how PETSc documents and uses its version information. This information is available to both C/C++
2394: and Fortran compilers when petscsys.h is included.
2397: The current PETSc version and the API for accessing it are defined in petscversion.h
2399: The complete version number is given as the triple PETSC_VERSION_MAJOR.PETSC_VERSION_MINOR.PETSC_VERSION_SUBMINOR (in short hand x.y.z)
2401: A change in the minor version number (y) indicates possible/likely changes in the PETSc API. Note this is different than with the semantic versioning convention
2402: where only a change in the major version number (x) indicates a change in the API.
2404: A subminor greater than zero indicates a patch release. Version x.y.z maintains source and binary compatibility with version x.y.w for all z and w
2406: Use the macros PETSC_VERSION_EQ(x,y,z), PETSC_VERSION_LT(x,y,z), PETSC_VERSION_LE(x,y,z), PETSC_VERSION_GT(x,y,z),
2407: PETSC_VERSION_GE(x,y,z) to determine if the current version is equal to, less than, less than or equal to, greater than or greater than or equal to a given
2408: version number (x.y.z).
2410: PETSC_RELEASE_DATE is the date the x.y version was released (i.e. the version before any patch releases)
2412: PETSC_VERSION_DATE is the date the x.y.z version was released
2414: PETSC_VERSION_GIT is the last git commit to the repository given in the form vx.y.z-wwwww
2416: PETSC_VERSION_DATE_GIT is the date of the last git commit to the repository
2418: Level: intermediate
2420: PETSC_VERSION_() and PETSC_VERSION_PATCH are deprecated and will eventually be removed. For several releases PETSC_VERSION_PATCH is always 0
2422: M*/
2424: /*MC
2426: UsingFortran - To use PETSc with Fortran you must use both PETSc include files and modules. At the beginning
2427: of every function and module definition you need something like
2429: $
2430: $#include "petsc/finclude/petscXXX.h"
2431: $ use petscXXX
2433: You can declare PETSc variables using either of the following.
2435: $ XXX variablename
2436: $ type(tXXX) variablename
2438: For example,
2440: $#include "petsc/finclude/petscvec.h"
2441: $ use petscvec
2442: $
2443: $ Vec b
2444: $ type(tVec) x
2446: Level: beginner
2448: M*/
2450: PETSC_EXTERN PetscErrorCode PetscGetArchType(char[],size_t);
2451: PETSC_EXTERN PetscErrorCode PetscGetHostName(char[],size_t);
2452: PETSC_EXTERN PetscErrorCode PetscGetUserName(char[],size_t);
2453: PETSC_EXTERN PetscErrorCode PetscGetProgramName(char[],size_t);
2454: PETSC_EXTERN PetscErrorCode PetscSetProgramName(const char[]);
2455: PETSC_EXTERN PetscErrorCode PetscGetDate(char[],size_t);
2456: PETSC_EXTERN PetscErrorCode PetscGetVersion(char[], size_t);
2457: PETSC_EXTERN PetscErrorCode PetscGetVersionNumber(PetscInt*,PetscInt*,PetscInt*,PetscInt*);
2459: PETSC_EXTERN PetscErrorCode PetscSortedInt(PetscInt,const PetscInt[],PetscBool*);
2460: PETSC_EXTERN PetscErrorCode PetscSortedMPIInt(PetscInt,const PetscMPIInt[],PetscBool*);
2461: PETSC_EXTERN PetscErrorCode PetscSortedReal(PetscInt,const PetscReal[],PetscBool*);
2462: PETSC_EXTERN PetscErrorCode PetscSortInt(PetscInt,PetscInt[]);
2463: PETSC_EXTERN PetscErrorCode PetscSortReverseInt(PetscInt,PetscInt[]);
2464: PETSC_EXTERN PetscErrorCode PetscSortedRemoveDupsInt(PetscInt*,PetscInt[]);
2465: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsInt(PetscInt*,PetscInt[]);
2467: PETSC_EXTERN PetscErrorCode PetscFindInt(PetscInt, PetscInt, const PetscInt[], PetscInt*);
2468: PETSC_EXTERN PetscErrorCode PetscFindMPIInt(PetscMPIInt, PetscInt, const PetscMPIInt[], PetscInt*);
2469: PETSC_EXTERN PetscErrorCode PetscSortIntWithPermutation(PetscInt,const PetscInt[],PetscInt[]);
2470: PETSC_EXTERN PetscErrorCode PetscSortStrWithPermutation(PetscInt,const char*[],PetscInt[]);
2471: PETSC_EXTERN PetscErrorCode PetscSortIntWithArray(PetscInt,PetscInt[],PetscInt[]);
2472: PETSC_EXTERN PetscErrorCode PetscSortIntWithArrayPair(PetscInt,PetscInt[],PetscInt[],PetscInt[]);
2473: PETSC_EXTERN PetscErrorCode PetscSortMPIInt(PetscInt,PetscMPIInt[]);
2474: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsMPIInt(PetscInt*,PetscMPIInt[]);
2475: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithArray(PetscMPIInt,PetscMPIInt[],PetscMPIInt[]);
2476: PETSC_EXTERN PetscErrorCode PetscSortMPIIntWithIntArray(PetscMPIInt,PetscMPIInt[],PetscInt[]);
2477: PETSC_EXTERN PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt[],PetscScalar[]);
2478: PETSC_EXTERN PetscErrorCode PetscSortIntWithDataArray(PetscInt,PetscInt[],void*,size_t,void*);
2479: PETSC_EXTERN PetscErrorCode PetscSortReal(PetscInt,PetscReal[]);
2480: PETSC_EXTERN PetscErrorCode PetscSortRealWithArrayInt(PetscInt,PetscReal[],PetscInt[]);
2481: PETSC_EXTERN PetscErrorCode PetscSortRealWithPermutation(PetscInt,const PetscReal[],PetscInt[]);
2482: PETSC_EXTERN PetscErrorCode PetscSortRemoveDupsReal(PetscInt*,PetscReal[]);
2483: PETSC_EXTERN PetscErrorCode PetscFindReal(PetscReal,PetscInt,const PetscReal[],PetscReal,PetscInt*);
2484: PETSC_EXTERN PetscErrorCode PetscSortSplit(PetscInt,PetscInt,PetscScalar[],PetscInt[]);
2485: PETSC_EXTERN PetscErrorCode PetscSortSplitReal(PetscInt,PetscInt,PetscReal[],PetscInt[]);
2486: PETSC_EXTERN PetscErrorCode PetscProcessTree(PetscInt,const PetscBool [],const PetscInt[],PetscInt*,PetscInt**,PetscInt**,PetscInt**,PetscInt**);
2487: PETSC_EXTERN PetscErrorCode PetscMergeIntArrayPair(PetscInt,const PetscInt[],const PetscInt[],PetscInt,const PetscInt[],const PetscInt[],PetscInt*,PetscInt**,PetscInt**);
2488: PETSC_EXTERN PetscErrorCode PetscMergeIntArray(PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt*,PetscInt**);
2489: PETSC_EXTERN PetscErrorCode PetscMergeMPIIntArray(PetscInt,const PetscMPIInt[],PetscInt,const PetscMPIInt[],PetscInt*,PetscMPIInt**);
2490: PETSC_EXTERN PetscErrorCode PetscParallelSortedInt(MPI_Comm, PetscInt, const PetscInt[], PetscBool *);
2492: PETSC_EXTERN PetscErrorCode PetscTimSort(PetscInt,void*,size_t,int(*)(const void*,const void*,void*),void*);
2493: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrdered(PetscInt,PetscInt[]);
2494: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrdered(PetscInt,PetscMPIInt[]);
2495: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrdered(PetscInt,PetscReal[]);
2496: PETSC_EXTERN PetscErrorCode PetscTimSortWithArray(PetscInt,void*,size_t,void*,size_t,int(*)(const void*,const void*,void*),void*);
2497: PETSC_EXTERN PetscErrorCode PetscIntSortSemiOrderedWithArray(PetscInt,PetscInt[],PetscInt[]);
2498: PETSC_EXTERN PetscErrorCode PetscMPIIntSortSemiOrderedWithArray(PetscInt,PetscMPIInt[],PetscMPIInt[]);
2499: PETSC_EXTERN PetscErrorCode PetscRealSortSemiOrderedWithArrayInt(PetscInt,PetscReal[],PetscInt[]);
2501: PETSC_EXTERN PetscErrorCode PetscSetDisplay(void);
2502: PETSC_EXTERN PetscErrorCode PetscGetDisplay(char[],size_t);
2504: /*J
2505: PetscRandomType - String with the name of a PETSc randomizer
2507: Level: beginner
2509: Notes:
2510: To use SPRNG or RANDOM123 you must have ./configure PETSc
2511: with the option --download-sprng or --download-random123
2513: .seealso: PetscRandomSetType(), PetscRandom, PetscRandomCreate()
2514: J*/
2515: typedef const char* PetscRandomType;
2516: #define PETSCRAND "rand"
2517: #define PETSCRAND48 "rand48"
2518: #define PETSCSPRNG "sprng"
2519: #define PETSCRANDER48 "rander48"
2520: #define PETSCRANDOM123 "random123"
2522: /* Logging support */
2523: PETSC_EXTERN PetscClassId PETSC_RANDOM_CLASSID;
2525: PETSC_EXTERN PetscErrorCode PetscRandomInitializePackage(void);
2527: /* Dynamic creation and loading functions */
2528: PETSC_EXTERN PetscFunctionList PetscRandomList;
2530: PETSC_EXTERN PetscErrorCode PetscRandomRegister(const char[],PetscErrorCode (*)(PetscRandom));
2531: PETSC_EXTERN PetscErrorCode PetscRandomSetType(PetscRandom, PetscRandomType);
2532: PETSC_EXTERN PetscErrorCode PetscRandomSetFromOptions(PetscRandom);
2533: PETSC_EXTERN PetscErrorCode PetscRandomGetType(PetscRandom, PetscRandomType*);
2534: PETSC_EXTERN PetscErrorCode PetscRandomViewFromOptions(PetscRandom,PetscObject,const char[]);
2535: PETSC_EXTERN PetscErrorCode PetscRandomView(PetscRandom,PetscViewer);
2537: PETSC_EXTERN PetscErrorCode PetscRandomCreate(MPI_Comm,PetscRandom*);
2538: PETSC_EXTERN PetscErrorCode PetscRandomGetValue(PetscRandom,PetscScalar*);
2539: PETSC_EXTERN PetscErrorCode PetscRandomGetValueReal(PetscRandom,PetscReal*);
2540: PETSC_EXTERN PetscErrorCode PetscRandomGetInterval(PetscRandom,PetscScalar*,PetscScalar*);
2541: PETSC_EXTERN PetscErrorCode PetscRandomSetInterval(PetscRandom,PetscScalar,PetscScalar);
2542: PETSC_EXTERN PetscErrorCode PetscRandomSetSeed(PetscRandom,unsigned long);
2543: PETSC_EXTERN PetscErrorCode PetscRandomGetSeed(PetscRandom,unsigned long *);
2544: PETSC_EXTERN PetscErrorCode PetscRandomSeed(PetscRandom);
2545: PETSC_EXTERN PetscErrorCode PetscRandomDestroy(PetscRandom*);
2547: PETSC_EXTERN PetscErrorCode PetscGetFullPath(const char[],char[],size_t);
2548: PETSC_EXTERN PetscErrorCode PetscGetRelativePath(const char[],char[],size_t);
2549: PETSC_EXTERN PetscErrorCode PetscGetWorkingDirectory(char[],size_t);
2550: PETSC_EXTERN PetscErrorCode PetscGetRealPath(const char[],char[]);
2551: PETSC_EXTERN PetscErrorCode PetscGetHomeDirectory(char[],size_t);
2552: PETSC_EXTERN PetscErrorCode PetscTestFile(const char[],char,PetscBool *);
2553: PETSC_EXTERN PetscErrorCode PetscTestDirectory(const char[],char,PetscBool *);
2554: PETSC_EXTERN PetscErrorCode PetscMkdir(const char[]);
2555: PETSC_EXTERN PetscErrorCode PetscMkdtemp(char[]);
2556: PETSC_EXTERN PetscErrorCode PetscRMTree(const char[]);
2558: PETSC_STATIC_INLINE PetscBool PetscBinaryBigEndian(void) {long _petsc_v = 1; return ((char*)&_petsc_v)[0] ? PETSC_FALSE : PETSC_TRUE;}
2560: PETSC_EXTERN PetscErrorCode PetscBinaryRead(int,void*,PetscInt,PetscInt*,PetscDataType);
2561: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedRead(MPI_Comm,int,void*,PetscInt,PetscInt*,PetscDataType);
2562: PETSC_EXTERN PetscErrorCode PetscBinaryWrite(int,const void*,PetscInt,PetscDataType);
2563: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedWrite(MPI_Comm,int,const void*,PetscInt,PetscDataType);
2564: PETSC_EXTERN PetscErrorCode PetscBinaryOpen(const char[],PetscFileMode,int *);
2565: PETSC_EXTERN PetscErrorCode PetscBinaryClose(int);
2566: PETSC_EXTERN PetscErrorCode PetscSharedTmp(MPI_Comm,PetscBool *);
2567: PETSC_EXTERN PetscErrorCode PetscSharedWorkingDirectory(MPI_Comm,PetscBool *);
2568: PETSC_EXTERN PetscErrorCode PetscGetTmp(MPI_Comm,char[],size_t);
2569: PETSC_EXTERN PetscErrorCode PetscFileRetrieve(MPI_Comm,const char[],char[],size_t,PetscBool *);
2570: PETSC_EXTERN PetscErrorCode PetscLs(MPI_Comm,const char[],char[],size_t,PetscBool *);
2571: #if defined(PETSC_USE_SOCKET_VIEWER)
2572: PETSC_EXTERN PetscErrorCode PetscOpenSocket(const char[],int,int*);
2573: #endif
2575: PETSC_EXTERN PetscErrorCode PetscBinarySeek(int,off_t,PetscBinarySeekType,off_t*);
2576: PETSC_EXTERN PetscErrorCode PetscBinarySynchronizedSeek(MPI_Comm,int,off_t,PetscBinarySeekType,off_t*);
2577: PETSC_EXTERN PetscErrorCode PetscByteSwap(void *,PetscDataType,PetscInt);
2579: PETSC_EXTERN PetscErrorCode PetscSetDebugTerminal(const char[]);
2580: PETSC_EXTERN PetscErrorCode PetscSetDebugger(const char[],PetscBool);
2581: PETSC_EXTERN PetscErrorCode PetscSetDefaultDebugger(void);
2582: PETSC_EXTERN PetscErrorCode PetscSetDebuggerFromString(const char*);
2583: PETSC_EXTERN PetscErrorCode PetscAttachDebugger(void);
2584: PETSC_EXTERN PetscErrorCode PetscStopForDebugger(void);
2585: PETSC_EXTERN PetscErrorCode PetscWaitOnError(void);
2587: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt*);
2588: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**);
2589: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths2(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscMPIInt**,PetscMPIInt**,PetscMPIInt**);
2590: PETSC_EXTERN PetscErrorCode PetscPostIrecvInt(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscInt***,MPI_Request**);
2591: PETSC_EXTERN PetscErrorCode PetscPostIrecvScalar(MPI_Comm,PetscMPIInt,PetscMPIInt,const PetscMPIInt[],const PetscMPIInt[],PetscScalar***,MPI_Request**);
2592: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSided(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt*,const void*,PetscMPIInt*,PetscMPIInt**,void*)
2593: PetscAttrMPIPointerWithType(6,3);
2594: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2595: PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2596: PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2597: PetscAttrMPIPointerWithType(6,3);
2598: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,
2599: MPI_Request**,MPI_Request**,
2600: PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
2601: PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
2602: PetscAttrMPIPointerWithType(6,3);
2604: PETSC_EXTERN const char *const PetscBuildTwoSidedTypes[];
2605: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm,PetscBuildTwoSidedType);
2606: PETSC_EXTERN PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm,PetscBuildTwoSidedType*);
2608: PETSC_EXTERN PetscErrorCode PetscSSEIsEnabled(MPI_Comm,PetscBool*,PetscBool*);
2610: PETSC_EXTERN MPI_Comm PetscObjectComm(PetscObject);
2612: PETSC_EXTERN const char *const PetscSubcommTypes[];
2614: struct _n_PetscSubcomm {
2615: MPI_Comm parent; /* parent communicator */
2616: MPI_Comm dupparent; /* duplicate parent communicator, under which the processors of this subcomm have contiguous rank */
2617: MPI_Comm child; /* the sub-communicator */
2618: PetscMPIInt n; /* num of subcommunicators under the parent communicator */
2619: PetscMPIInt color; /* color of processors belong to this communicator */
2620: PetscMPIInt *subsize; /* size of subcommunicator[color] */
2621: PetscSubcommType type;
2622: char *subcommprefix;
2623: };
2625: PETSC_STATIC_INLINE MPI_Comm PetscSubcommParent(PetscSubcomm scomm) {return scomm->parent;}
2626: PETSC_STATIC_INLINE MPI_Comm PetscSubcommChild(PetscSubcomm scomm) {return scomm->child;}
2627: PETSC_STATIC_INLINE MPI_Comm PetscSubcommContiguousParent(PetscSubcomm scomm) {return scomm->dupparent;}
2628: PETSC_EXTERN PetscErrorCode PetscSubcommCreate(MPI_Comm,PetscSubcomm*);
2629: PETSC_EXTERN PetscErrorCode PetscSubcommDestroy(PetscSubcomm*);
2630: PETSC_EXTERN PetscErrorCode PetscSubcommSetNumber(PetscSubcomm,PetscInt);
2631: PETSC_EXTERN PetscErrorCode PetscSubcommSetType(PetscSubcomm,PetscSubcommType);
2632: PETSC_EXTERN PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm,PetscMPIInt,PetscMPIInt);
2633: PETSC_EXTERN PetscErrorCode PetscSubcommView(PetscSubcomm,PetscViewer);
2634: PETSC_EXTERN PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm);
2635: PETSC_EXTERN PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm,const char[]);
2636: PETSC_EXTERN PetscErrorCode PetscSubcommGetParent(PetscSubcomm,MPI_Comm*);
2637: PETSC_EXTERN PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm,MPI_Comm*);
2638: PETSC_EXTERN PetscErrorCode PetscSubcommGetChild(PetscSubcomm,MPI_Comm*);
2640: PETSC_EXTERN PetscErrorCode PetscHeapCreate(PetscInt,PetscHeap*);
2641: PETSC_EXTERN PetscErrorCode PetscHeapAdd(PetscHeap,PetscInt,PetscInt);
2642: PETSC_EXTERN PetscErrorCode PetscHeapPop(PetscHeap,PetscInt*,PetscInt*);
2643: PETSC_EXTERN PetscErrorCode PetscHeapPeek(PetscHeap,PetscInt*,PetscInt*);
2644: PETSC_EXTERN PetscErrorCode PetscHeapStash(PetscHeap,PetscInt,PetscInt);
2645: PETSC_EXTERN PetscErrorCode PetscHeapUnstash(PetscHeap);
2646: PETSC_EXTERN PetscErrorCode PetscHeapDestroy(PetscHeap*);
2647: PETSC_EXTERN PetscErrorCode PetscHeapView(PetscHeap,PetscViewer);
2649: PETSC_EXTERN PetscErrorCode PetscProcessPlacementView(PetscViewer);
2650: PETSC_EXTERN PetscErrorCode PetscShmCommGet(MPI_Comm,PetscShmComm*);
2651: PETSC_EXTERN PetscErrorCode PetscShmCommGlobalToLocal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2652: PETSC_EXTERN PetscErrorCode PetscShmCommLocalToGlobal(PetscShmComm,PetscMPIInt,PetscMPIInt*);
2653: PETSC_EXTERN PetscErrorCode PetscShmCommGetMpiShmComm(PetscShmComm,MPI_Comm*);
2655: /* routines to better support OpenMP multithreading needs of some PETSc third party libraries */
2656: PETSC_EXTERN PetscErrorCode PetscOmpCtrlCreate(MPI_Comm,PetscInt,PetscOmpCtrl*);
2657: PETSC_EXTERN PetscErrorCode PetscOmpCtrlGetOmpComms(PetscOmpCtrl,MPI_Comm*,MPI_Comm*,PetscBool*);
2658: PETSC_EXTERN PetscErrorCode PetscOmpCtrlDestroy(PetscOmpCtrl*);
2659: PETSC_EXTERN PetscErrorCode PetscOmpCtrlBarrier(PetscOmpCtrl);
2660: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterBegin(PetscOmpCtrl);
2661: PETSC_EXTERN PetscErrorCode PetscOmpCtrlOmpRegionOnMasterEnd(PetscOmpCtrl);
2663: PETSC_EXTERN PetscErrorCode PetscSegBufferCreate(size_t,size_t,PetscSegBuffer*);
2664: PETSC_EXTERN PetscErrorCode PetscSegBufferDestroy(PetscSegBuffer*);
2665: PETSC_EXTERN PetscErrorCode PetscSegBufferGet(PetscSegBuffer,size_t,void*);
2666: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractAlloc(PetscSegBuffer,void*);
2667: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractTo(PetscSegBuffer,void*);
2668: PETSC_EXTERN PetscErrorCode PetscSegBufferExtractInPlace(PetscSegBuffer,void*);
2669: PETSC_EXTERN PetscErrorCode PetscSegBufferGetSize(PetscSegBuffer,size_t*);
2670: PETSC_EXTERN PetscErrorCode PetscSegBufferUnuse(PetscSegBuffer,size_t);
2674: * prevents the compiler from completely erasing the stub. This is called in inner loops so it has to be as fast as
2675: * possible. */
2676: PETSC_STATIC_INLINE PetscErrorCode PetscSegBufferGetInts(PetscSegBuffer seg,size_t count,PetscInt *PETSC_RESTRICT *slot) {return PetscSegBufferGet(seg,count,(void**)slot);}
2678: extern PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton;
2679: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted*);
2680: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted*);
2681: PETSC_EXTERN PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted,const char*,const char*,PetscBool*);
2683: #include <stdarg.h>
2684: PETSC_EXTERN PetscErrorCode PetscVSNPrintf(char*,size_t,const char[],size_t*,va_list);
2685: PETSC_EXTERN PetscErrorCode (*PetscVFPrintf)(FILE*,const char[],va_list);
2687: PETSC_EXTERN PetscSegBuffer PetscCitationsList;
2689: /*@C
2690: PetscCitationsRegister - Register a bibtex item to obtain credit for an implemented algorithm used in the code.
2692: Not Collective - only what is registered on rank 0 of PETSC_COMM_WORLD will be printed
2694: Input Parameters:
2695: + cite - the bibtex item, formated to displayed on multiple lines nicely
2696: - set - a boolean variable initially set to PETSC_FALSE; this is used to insure only a single registration of the citation
2698: Level: intermediate
2700: Not available from Fortran
2702: Options Database:
2703: . -citations [filename] - print out the bibtex entries for the given computation
2704: @*/
2705: PETSC_STATIC_INLINE PetscErrorCode PetscCitationsRegister(const char cit[],PetscBool *set)
2706: {
2707: size_t len;
2708: char *vstring;
2712: if (set && *set) return(0);
2713: PetscStrlen(cit,&len);
2714: PetscSegBufferGet(PetscCitationsList,len,&vstring);
2715: PetscArraycpy(vstring,cit,len);
2716: if (set) *set = PETSC_TRUE;
2717: return(0);
2718: }
2720: PETSC_EXTERN PetscErrorCode PetscURLShorten(const char[],char[],size_t);
2721: PETSC_EXTERN PetscErrorCode PetscGoogleDriveAuthorize(MPI_Comm,char[],char[],size_t);
2722: PETSC_EXTERN PetscErrorCode PetscGoogleDriveRefresh(MPI_Comm,const char[],char[],size_t);
2723: PETSC_EXTERN PetscErrorCode PetscGoogleDriveUpload(MPI_Comm,const char[],const char []);
2725: PETSC_EXTERN PetscErrorCode PetscBoxAuthorize(MPI_Comm,char[],char[],size_t);
2726: PETSC_EXTERN PetscErrorCode PetscBoxRefresh(MPI_Comm,const char[],char[],char[],size_t);
2728: PETSC_EXTERN PetscErrorCode PetscGlobusGetTransfers(MPI_Comm,const char[],char[],size_t);
2730: PETSC_EXTERN PetscErrorCode PetscTextBelt(MPI_Comm,const char[],const char[],PetscBool*);
2731: PETSC_EXTERN PetscErrorCode PetscTellMyCell(MPI_Comm,const char[],const char[],PetscBool*);
2733: PETSC_EXTERN PetscErrorCode PetscPullJSONValue(const char[],const char[],char[],size_t,PetscBool*);
2734: PETSC_EXTERN PetscErrorCode PetscPushJSONValue(char[],const char[],const char[],size_t);
2737: #if defined(PETSC_USE_DEBUG)
2738: /*
2739: Verify that all processes in the communicator have called this from the same line of code
2740: */
2741: PETSC_EXTERN PetscErrorCode PetscAllreduceBarrierCheck(MPI_Comm,PetscMPIInt,int,const char*,const char *);
2743: /*MC
2744: MPIU_Allreduce - a PETSc replacement for MPI_Allreduce() that tries to determine if the call from all the MPI processes occur from the
2745: same place in the PETSc code. This helps to detect bugs where different MPI processes follow different code paths
2746: resulting in inconsistent and incorrect calls to MPI_Allreduce().
2748: Collective
2750: Synopsis:
2751: #include <petscsys.h>
2752: PetscErrorCode MPIU_Allreduce(void *indata,void *outdata,PetscMPIInt count,MPI_Datatype datatype, MPI_Op op, MPI_Comm comm);
2754: Input Parameters:
2755: + indata - pointer to the input data to be reduced
2756: . count - the number of MPI data items in indata and outdata
2757: . datatype - the MPI datatype, for example MPI_INT
2758: . op - the MPI operation, for example MPI_SUM
2759: - comm - the MPI communicator on which the operation occurs
2761: Output Parameter:
2762: . outdata - the reduced values
2764: Notes:
2765: In optimized mode this directly calls MPI_Allreduce()
2767: Level: developer
2769: .seealso: MPI_Allreduce()
2770: M*/
2771: #define MPIU_Allreduce(a,b,c,d,e,fcomm) (PetscAllreduceBarrierCheck(fcomm,c,__LINE__,PETSC_FUNCTION_NAME,__FILE__) || MPI_Allreduce(a,b,c,d,e,fcomm))
2772: #else
2773: #define MPIU_Allreduce(a,b,c,d,e,fcomm) MPI_Allreduce(a,b,c,d,e,fcomm)
2774: #endif
2776: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
2777: PETSC_EXTERN PetscErrorCode MPIU_Win_allocate_shared(MPI_Aint,PetscMPIInt,MPI_Info,MPI_Comm,void*,MPI_Win*);
2778: PETSC_EXTERN PetscErrorCode MPIU_Win_shared_query(MPI_Win,PetscMPIInt,MPI_Aint*,PetscMPIInt*,void*);
2779: #endif
2781: /*
2782: List of external packages and queries on it
2783: */
2784: PETSC_EXTERN PetscErrorCode PetscHasExternalPackage(const char[],PetscBool*);
2786: #endif