Actual source code: petscmat.h
1: /*
2: Include file for the matrix component of PETSc
3: */
4: #ifndef PETSCMAT_H
5: #define PETSCMAT_H
7: #include <petscvec.h>
9: /* SUBMANSEC = Mat */
11: /*S
12: Mat - Abstract PETSc matrix object used to manage all linear operators in PETSc, even those without
13: an explicit sparse representation (such as matrix-free operators)
15: Level: beginner
17: .seealso: `MatCreate()`, `MatType`, `MatSetType()`, `MatDestroy()`
18: S*/
19: typedef struct _p_Mat *Mat;
21: /*J
22: MatType - String with the name of a PETSc matrix type
24: Level: beginner
26: .seealso: `MatSetType()`, `Mat`, `MatSolverType`, `MatRegister()`
27: J*/
28: typedef const char *MatType;
29: #define MATSAME "same"
30: #define MATMAIJ "maij"
31: #define MATSEQMAIJ "seqmaij"
32: #define MATMPIMAIJ "mpimaij"
33: #define MATKAIJ "kaij"
34: #define MATSEQKAIJ "seqkaij"
35: #define MATMPIKAIJ "mpikaij"
36: #define MATIS "is"
37: #define MATAIJ "aij"
38: #define MATSEQAIJ "seqaij"
39: #define MATMPIAIJ "mpiaij"
40: #define MATAIJCRL "aijcrl"
41: #define MATSEQAIJCRL "seqaijcrl"
42: #define MATMPIAIJCRL "mpiaijcrl"
43: #define MATAIJCUSPARSE "aijcusparse"
44: #define MATSEQAIJCUSPARSE "seqaijcusparse"
45: #define MATMPIAIJCUSPARSE "mpiaijcusparse"
46: #define MATAIJKOKKOS "aijkokkos"
47: #define MATSEQAIJKOKKOS "seqaijkokkos"
48: #define MATMPIAIJKOKKOS "mpiaijkokkos"
49: #define MATAIJVIENNACL "aijviennacl"
50: #define MATSEQAIJVIENNACL "seqaijviennacl"
51: #define MATMPIAIJVIENNACL "mpiaijviennacl"
52: #define MATAIJPERM "aijperm"
53: #define MATSEQAIJPERM "seqaijperm"
54: #define MATMPIAIJPERM "mpiaijperm"
55: #define MATAIJSELL "aijsell"
56: #define MATSEQAIJSELL "seqaijsell"
57: #define MATMPIAIJSELL "mpiaijsell"
58: #define MATAIJMKL "aijmkl"
59: #define MATSEQAIJMKL "seqaijmkl"
60: #define MATMPIAIJMKL "mpiaijmkl"
61: #define MATBAIJMKL "baijmkl"
62: #define MATSEQBAIJMKL "seqbaijmkl"
63: #define MATMPIBAIJMKL "mpibaijmkl"
64: #define MATSHELL "shell"
65: #define MATCENTERING "centering"
66: #define MATDENSE "dense"
67: #define MATDENSECUDA "densecuda"
68: #define MATSEQDENSE "seqdense"
69: #define MATSEQDENSECUDA "seqdensecuda"
70: #define MATMPIDENSE "mpidense"
71: #define MATMPIDENSECUDA "mpidensecuda"
72: #define MATELEMENTAL "elemental"
73: #define MATSCALAPACK "scalapack"
74: #define MATBAIJ "baij"
75: #define MATSEQBAIJ "seqbaij"
76: #define MATMPIBAIJ "mpibaij"
77: #define MATMPIADJ "mpiadj"
78: #define MATSBAIJ "sbaij"
79: #define MATSEQSBAIJ "seqsbaij"
80: #define MATMPISBAIJ "mpisbaij"
81: #define MATMFFD "mffd"
82: #define MATNORMAL "normal"
83: #define MATNORMALHERMITIAN "normalh"
84: #define MATLRC "lrc"
85: #define MATSCATTER "scatter"
86: #define MATBLOCKMAT "blockmat"
87: #define MATCOMPOSITE "composite"
88: #define MATFFT "fft"
89: #define MATFFTW "fftw"
90: #define MATSEQCUFFT "seqcufft"
91: #define MATTRANSPOSEMAT PETSC_DEPRECATED_MACRO("GCC warning \"MATTRANSPOSEMAT macro is deprecated use MATTRANSPOSEVIRTUAL (since version 3.18)\"") "transpose"
92: #define MATTRANSPOSEVIRTUAL "transpose"
93: #define MATHERMITIANTRANSPOSEVIRTUAL "hermitiantranspose"
94: #define MATSCHURCOMPLEMENT "schurcomplement"
95: #define MATPYTHON "python"
96: #define MATHYPRE "hypre"
97: #define MATHYPRESTRUCT "hyprestruct"
98: #define MATHYPRESSTRUCT "hypresstruct"
99: #define MATSUBMATRIX "submatrix"
100: #define MATLOCALREF "localref"
101: #define MATNEST "nest"
102: #define MATPREALLOCATOR "preallocator"
103: #define MATSELL "sell"
104: #define MATSEQSELL "seqsell"
105: #define MATMPISELL "mpisell"
106: #define MATDUMMY "dummy"
107: #define MATLMVM "lmvm"
108: #define MATLMVMDFP "lmvmdfp"
109: #define MATLMVMBFGS "lmvmbfgs"
110: #define MATLMVMSR1 "lmvmsr1"
111: #define MATLMVMBROYDEN "lmvmbroyden"
112: #define MATLMVMBADBROYDEN "lmvmbadbroyden"
113: #define MATLMVMSYMBROYDEN "lmvmsymbroyden"
114: #define MATLMVMSYMBADBROYDEN "lmvmsymbadbroyden"
115: #define MATLMVMDIAGBROYDEN "lmvmdiagbroyden"
116: #define MATCONSTANTDIAGONAL "constantdiagonal"
117: #define MATHTOOL "htool"
118: #define MATH2OPUS "h2opus"
120: /*J
121: MatSolverType - String with the name of a PETSc matrix solver type.
123: For example: "petsc" indicates what PETSc provides, "superlu_dist" the parallel SuperLU_DIST package etc
125: Level: beginner
127: Notes: MATSOLVERUMFPACK, MATSOLVERCHOLMOD, MATSOLVERKLU, MATSOLVERSPQR form the SuiteSparse package for which you can use --download-suitesparse
129: .seealso: `MatGetFactor()`, `PCFactorSetMatSolverType()`, `PCFactorGetMatSolverType()`
130: J*/
131: typedef const char *MatSolverType;
132: #define MATSOLVERSUPERLU "superlu"
133: #define MATSOLVERSUPERLU_DIST "superlu_dist"
134: #define MATSOLVERSTRUMPACK "strumpack"
135: #define MATSOLVERUMFPACK "umfpack"
136: #define MATSOLVERCHOLMOD "cholmod"
137: #define MATSOLVERKLU "klu"
138: #define MATSOLVERSPARSEELEMENTAL "sparseelemental"
139: #define MATSOLVERELEMENTAL "elemental"
140: #define MATSOLVERSCALAPACK "scalapack"
141: #define MATSOLVERESSL "essl"
142: #define MATSOLVERLUSOL "lusol"
143: #define MATSOLVERMUMPS "mumps"
144: #define MATSOLVERMKL_PARDISO "mkl_pardiso"
145: #define MATSOLVERMKL_CPARDISO "mkl_cpardiso"
146: #define MATSOLVERPASTIX "pastix"
147: #define MATSOLVERMATLAB "matlab"
148: #define MATSOLVERPETSC "petsc"
149: #define MATSOLVERBAS "bas"
150: #define MATSOLVERCUSPARSE "cusparse"
151: #define MATSOLVERCUSPARSEBAND "cusparseband"
152: #define MATSOLVERCUDA "cuda"
153: #define MATSOLVERKOKKOS "kokkos"
154: #define MATSOLVERKOKKOSDEVICE "kokkosdevice"
155: #define MATSOLVERSPQR "spqr"
157: /*E
158: MatFactorType - indicates what type of factorization is requested
160: Level: beginner
162: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
164: .seealso: `MatSolverType`, `MatGetFactor()`, `MatGetFactorAvailable()`, `MatSolverTypeRegister()`
165: E*/
166: typedef enum {
167: MAT_FACTOR_NONE,
168: MAT_FACTOR_LU,
169: MAT_FACTOR_CHOLESKY,
170: MAT_FACTOR_ILU,
171: MAT_FACTOR_ICC,
172: MAT_FACTOR_ILUDT,
173: MAT_FACTOR_QR,
174: MAT_FACTOR_NUM_TYPES
175: } MatFactorType;
176: PETSC_EXTERN const char *const MatFactorTypes[];
178: PETSC_EXTERN PetscErrorCode MatGetFactor(Mat, MatSolverType, MatFactorType, Mat *);
179: PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat, MatSolverType, MatFactorType, PetscBool *);
180: PETSC_EXTERN PetscErrorCode MatFactorGetCanUseOrdering(Mat, PetscBool *);
181: PETSC_DEPRECATED_FUNCTION("Use MatFactorGetCanUseOrdering() (since version 3.15)") static inline PetscErrorCode MatFactorGetUseOrdering(Mat A, PetscBool *b)
182: {
183: return MatFactorGetCanUseOrdering(A, b);
184: }
185: PETSC_EXTERN PetscErrorCode MatFactorGetSolverType(Mat, MatSolverType *);
186: PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat, MatFactorType *);
187: PETSC_EXTERN PetscErrorCode MatSetFactorType(Mat, MatFactorType);
188: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*MatSolverFunction)(Mat, MatFactorType, Mat *);
189: PETSC_EXTERN PetscErrorCode MatSolverTypeRegister(MatSolverType, MatType, MatFactorType, MatSolverFunction);
190: PETSC_EXTERN PetscErrorCode MatSolverTypeGet(MatSolverType, MatType, MatFactorType, PetscBool *, PetscBool *, MatSolverFunction *);
191: typedef MatSolverType MatSolverPackage PETSC_DEPRECATED_TYPEDEF("Use MatSolverType (since version 3.9)");
192: PETSC_DEPRECATED_FUNCTION("Use MatSolverTypeRegister() (since version 3.9)") static inline PetscErrorCode MatSolverPackageRegister(MatSolverType stype, MatType mtype, MatFactorType ftype, MatSolverFunction f)
193: {
194: return MatSolverTypeRegister(stype, mtype, ftype, f);
195: }
196: PETSC_DEPRECATED_FUNCTION("Use MatSolverTypeGet() (since version 3.9)") static inline PetscErrorCode MatSolverPackageGet(MatSolverType stype, MatType mtype, MatFactorType ftype, PetscBool *foundmtype, PetscBool *foundstype, MatSolverFunction *f)
197: {
198: return MatSolverTypeGet(stype, mtype, ftype, foundmtype, foundstype, f);
199: }
201: /*E
202: MatProductType - indicates what type of matrix product is requested
204: Level: beginner
206: .seealso: `MatProductSetType()`
207: E*/
208: typedef enum {
209: MATPRODUCT_UNSPECIFIED = 0,
210: MATPRODUCT_AB,
211: MATPRODUCT_AtB,
212: MATPRODUCT_ABt,
213: MATPRODUCT_PtAP,
214: MATPRODUCT_RARt,
215: MATPRODUCT_ABC
216: } MatProductType;
217: PETSC_EXTERN const char *const MatProductTypes[];
219: /*J
220: MatProductAlgorithm - String with the name of an algorithm for a PETSc matrix product implementation
222: Level: beginner
224: .seealso: `MatSetType()`, `Mat`, `MatProductSetAlgorithm()`, `MatProductType`
225: J*/
226: typedef const char *MatProductAlgorithm;
227: #define MATPRODUCTALGORITHMDEFAULT "default"
228: #define MATPRODUCTALGORITHMSORTED "sorted"
229: #define MATPRODUCTALGORITHMSCALABLE "scalable"
230: #define MATPRODUCTALGORITHMSCALABLEFAST "scalable_fast"
231: #define MATPRODUCTALGORITHMHEAP "heap"
232: #define MATPRODUCTALGORITHMBHEAP "btheap"
233: #define MATPRODUCTALGORITHMLLCONDENSED "llcondensed"
234: #define MATPRODUCTALGORITHMROWMERGE "rowmerge"
235: #define MATPRODUCTALGORITHMOUTERPRODUCT "outerproduct"
236: #define MATPRODUCTALGORITHMATB "at*b"
237: #define MATPRODUCTALGORITHMRAP "rap"
238: #define MATPRODUCTALGORITHMNONSCALABLE "nonscalable"
239: #define MATPRODUCTALGORITHMSEQMPI "seqmpi"
240: #define MATPRODUCTALGORITHMBACKEND "backend"
241: #define MATPRODUCTALGORITHMOVERLAPPING "overlapping"
242: #define MATPRODUCTALGORITHMMERGED "merged"
243: #define MATPRODUCTALGORITHMALLATONCE "allatonce"
244: #define MATPRODUCTALGORITHMALLATONCEMERGED "allatonce_merged"
245: #define MATPRODUCTALGORITHMALLGATHERV "allgatherv"
246: #define MATPRODUCTALGORITHMCYCLIC "cyclic"
247: #if defined(PETSC_HAVE_HYPRE)
248: #define MATPRODUCTALGORITHMHYPRE "hypre"
249: #endif
251: PETSC_EXTERN PetscErrorCode MatProductCreate(Mat, Mat, Mat, Mat *);
252: PETSC_EXTERN PetscErrorCode MatProductCreateWithMat(Mat, Mat, Mat, Mat);
253: PETSC_EXTERN PetscErrorCode MatProductSetType(Mat, MatProductType);
254: PETSC_EXTERN PetscErrorCode MatProductSetAlgorithm(Mat, MatProductAlgorithm);
255: PETSC_EXTERN PetscErrorCode MatProductSetFill(Mat, PetscReal);
256: PETSC_EXTERN PetscErrorCode MatProductSetFromOptions(Mat);
257: PETSC_EXTERN PetscErrorCode MatProductSymbolic(Mat);
258: PETSC_EXTERN PetscErrorCode MatProductNumeric(Mat);
259: PETSC_EXTERN PetscErrorCode MatProductReplaceMats(Mat, Mat, Mat, Mat);
260: PETSC_EXTERN PetscErrorCode MatProductClear(Mat);
261: PETSC_EXTERN PetscErrorCode MatProductView(Mat, PetscViewer);
262: PETSC_EXTERN PetscErrorCode MatProductGetType(Mat, MatProductType *);
263: PETSC_EXTERN PetscErrorCode MatProductGetMats(Mat, Mat *, Mat *, Mat *);
265: /* Logging support */
266: #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */
267: PETSC_EXTERN PetscClassId MAT_CLASSID;
268: PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID;
269: PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID;
270: PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
271: PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID;
272: PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID;
273: PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID;
274: PETSC_EXTERN PetscClassId MATMFFD_CLASSID;
276: /*E
277: MatReuse - Indicates if matrices obtained from a previous call to MatCreateSubMatrices(), MatCreateSubMatrix(), MatConvert() or several other functions
278: are to be reused to store the new matrix values.
280: $ MAT_INITIAL_MATRIX - create a new matrix
281: $ MAT_REUSE_MATRIX - reuse the matrix created with a previous call that used MAT_INITIAL_MATRIX
282: $ MAT_INPLACE_MATRIX - replace the first input matrix with the new matrix (not applicable to all functions)
283: $ MAT_IGNORE_MATRIX - do not create a new matrix or reuse a give matrix, just ignore that matrix argument (not applicable to all functions)
285: Level: beginner
287: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
289: .seealso: `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatDestroyMatrices()`, `MatConvert()`
290: E*/
291: typedef enum {
292: MAT_INITIAL_MATRIX,
293: MAT_REUSE_MATRIX,
294: MAT_IGNORE_MATRIX,
295: MAT_INPLACE_MATRIX
296: } MatReuse;
298: /*E
299: MatCreateSubMatrixOption - Indicates if matrices obtained from a call to MatCreateSubMatrices()
300: include the matrix values. Currently it is only used by MatGetSeqNonzeroStructure().
302: Level: beginner
304: .seealso: `MatGetSeqNonzeroStructure()`
305: E*/
306: typedef enum {
307: MAT_DO_NOT_GET_VALUES,
308: MAT_GET_VALUES
309: } MatCreateSubMatrixOption;
311: PETSC_EXTERN PetscErrorCode MatInitializePackage(void);
313: PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm, Mat *);
314: PETSC_EXTERN PetscErrorCode MatSetSizes(Mat, PetscInt, PetscInt, PetscInt, PetscInt);
315: PETSC_EXTERN PetscErrorCode MatSetType(Mat, MatType);
316: PETSC_EXTERN PetscErrorCode MatGetVecType(Mat, VecType *);
317: PETSC_EXTERN PetscErrorCode MatSetVecType(Mat, VecType);
318: PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat);
319: PETSC_EXTERN PetscErrorCode MatViewFromOptions(Mat, PetscObject, const char[]);
320: PETSC_EXTERN PetscErrorCode MatRegister(const char[], PetscErrorCode (*)(Mat));
321: PETSC_EXTERN PetscErrorCode MatRegisterRootName(const char[], const char[], const char[]);
322: PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat, const char[]);
323: PETSC_EXTERN PetscErrorCode MatSetOptionsPrefixFactor(Mat, const char[]);
324: PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefixFactor(Mat, const char[]);
325: PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat, const char[]);
326: PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat, const char *[]);
327: PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat, PetscBool);
329: PETSC_EXTERN PetscFunctionList MatList;
330: PETSC_EXTERN PetscFunctionList MatColoringList;
331: PETSC_EXTERN PetscFunctionList MatPartitioningList;
333: /*E
334: MatStructure - Indicates if two matrices have the same nonzero structure
336: Level: beginner
338: $ SAME_NONZERO_PATTERN - the two matrices have identical nonzero patterns
339: $ DIFFERENT_NONZERO_PATTERN - the two matrices may have different nonzero patterns
340: $ SUBSET_NONZERO_PATTERN - the nonzero pattern of the second matrix is a subset of the nonzero pattern of the first matrix
341: $ UNKNOWN_NONZERO_PATTERN - there is no known relationship between the nonzero patterns. In this case the implementations may try to detect a relationship to optimize the operation
343: Developer Notes:
344: Any additions/changes here MUST also be made in src/mat/f90-mod/petscmat.h
346: .seealso: `MatCopy()`, `MatAXPY()`, `MatAYPX()`
347: E*/
348: typedef enum {
349: DIFFERENT_NONZERO_PATTERN,
350: SUBSET_NONZERO_PATTERN,
351: SAME_NONZERO_PATTERN,
352: UNKNOWN_NONZERO_PATTERN
353: } MatStructure;
354: PETSC_EXTERN const char *const MatStructures[];
356: #if defined PETSC_HAVE_MKL_SPARSE
357: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJMKL(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *);
358: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJMKL(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *);
359: PETSC_EXTERN PetscErrorCode MatCreateBAIJMKL(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *);
360: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJMKL(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *);
361: #endif
363: PETSC_EXTERN PetscErrorCode MatCreateSeqSELL(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *);
364: PETSC_EXTERN PetscErrorCode MatCreateSELL(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *);
365: PETSC_EXTERN PetscErrorCode MatSeqSELLSetPreallocation(Mat, PetscInt, const PetscInt[]);
366: PETSC_EXTERN PetscErrorCode MatMPISELLSetPreallocation(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[]);
368: PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm, PetscInt, PetscInt, PetscScalar[], Mat *);
369: PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar[], Mat *);
370: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *);
371: PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *);
372: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], Mat *);
373: PETSC_EXTERN PetscErrorCode MatUpdateMPIAIJWithArrays(Mat, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[]);
374: PETSC_EXTERN PetscErrorCode MatUpdateMPIAIJWithArray(Mat, const PetscScalar[]);
375: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], PetscInt[], PetscInt[], PetscScalar[], Mat *);
376: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSeqAIJ(MPI_Comm, Mat, Mat, const PetscInt[], Mat *);
378: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *);
379: PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *);
380: PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], Mat *);
382: PETSC_EXTERN PetscErrorCode MatSetPreallocationCOO(Mat, PetscCount, PetscInt[], PetscInt[]);
383: PETSC_EXTERN PetscErrorCode MatSetPreallocationCOOLocal(Mat, PetscCount, PetscInt[], PetscInt[]);
384: PETSC_EXTERN PetscErrorCode MatSetValuesCOO(Mat, const PetscScalar[], InsertMode);
386: PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscInt[], Mat *);
387: PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *);
389: PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *);
390: PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], Mat *);
391: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[]);
392: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[]);
393: PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat, PetscInt, const PetscInt[], const PetscInt[], const PetscInt[], const PetscInt[]);
395: PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, void *, Mat *);
396: PETSC_EXTERN PetscErrorCode MatCreateCentering(MPI_Comm, PetscInt, PetscInt, Mat *);
397: PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat, Mat *);
398: PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat, Mat *);
399: PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat, Mat, Vec, Mat, Mat *);
400: PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat, Mat *, Mat *, Vec *, Mat *);
401: PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, ISLocalToGlobalMapping, ISLocalToGlobalMapping, Mat *);
402: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *);
403: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *);
405: PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm, VecScatter, Mat *);
406: PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat, VecScatter);
407: PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat, VecScatter *);
408: PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, Mat *);
409: PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat, Mat);
410: PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
411: typedef enum {
412: MAT_COMPOSITE_MERGE_RIGHT,
413: MAT_COMPOSITE_MERGE_LEFT
414: } MatCompositeMergeType;
415: PETSC_EXTERN PetscErrorCode MatCompositeSetMergeType(Mat, MatCompositeMergeType);
416: PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm, PetscInt, const Mat *, Mat *);
417: typedef enum {
418: MAT_COMPOSITE_ADDITIVE,
419: MAT_COMPOSITE_MULTIPLICATIVE
420: } MatCompositeType;
421: PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat, MatCompositeType);
422: PETSC_EXTERN PetscErrorCode MatCompositeGetType(Mat, MatCompositeType *);
423: PETSC_EXTERN PetscErrorCode MatCompositeSetMatStructure(Mat, MatStructure);
424: PETSC_EXTERN PetscErrorCode MatCompositeGetMatStructure(Mat, MatStructure *);
425: PETSC_EXTERN PetscErrorCode MatCompositeGetNumberMat(Mat, PetscInt *);
426: PETSC_EXTERN PetscErrorCode MatCompositeGetMat(Mat, PetscInt, Mat *);
427: PETSC_EXTERN PetscErrorCode MatCompositeSetScalings(Mat, const PetscScalar *);
429: PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm, PetscInt, const PetscInt[], MatType, Mat *);
430: PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm, PetscInt, const PetscInt[], Mat *);
432: PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat, Mat *);
433: PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat, Mat *);
434: PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat, Mat *);
435: PETSC_EXTERN PetscErrorCode MatHermitianTransposeGetMat(Mat, Mat *);
436: PETSC_EXTERN PetscErrorCode MatNormalGetMat(Mat, Mat *);
437: PETSC_EXTERN PetscErrorCode MatNormalHermitianGetMat(Mat, Mat *);
438: PETSC_EXTERN PetscErrorCode MatCreateSubMatrixVirtual(Mat, IS, IS, Mat *);
439: PETSC_EXTERN PetscErrorCode MatSubMatrixVirtualUpdate(Mat, Mat, IS, IS);
440: PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat, IS, IS, Mat *);
441: PETSC_EXTERN PetscErrorCode MatCreateConstantDiagonal(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar, Mat *);
443: #if defined(PETSC_HAVE_HYPRE)
444: PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[]);
445: #endif
447: PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat, const char[]);
448: PETSC_EXTERN PetscErrorCode MatPythonGetType(Mat, const char *[]);
450: PETSC_EXTERN PetscErrorCode MatResetPreallocation(Mat);
451: PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
452: PETSC_EXTERN PetscErrorCode MatDestroy(Mat *);
453: PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat, PetscObjectState *);
455: PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
456: PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
457: PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
458: PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat, Mat *);
459: PETSC_EXTERN PetscErrorCode MatGetTrace(Mat, PetscScalar *);
460: PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat, const PetscScalar **);
461: PETSC_EXTERN PetscErrorCode MatInvertVariableBlockDiagonal(Mat, PetscInt, const PetscInt *, PetscScalar *);
462: PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonalMat(Mat, Mat);
463: PETSC_EXTERN PetscErrorCode MatInvertVariableBlockEnvelope(Mat, MatReuse, Mat *);
465: /* ------------------------------------------------------------*/
466: PETSC_EXTERN PetscErrorCode MatSetValues(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
467: PETSC_EXTERN PetscErrorCode MatSetValuesIS(Mat, IS, IS, const PetscScalar[], InsertMode);
468: PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
469: PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat, PetscInt, const PetscScalar[]);
470: PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat, PetscInt, const PetscScalar[]);
471: PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat, PetscInt, PetscInt, PetscInt[], const PetscScalar[]);
472: PETSC_EXTERN PetscErrorCode MatSetRandom(Mat, PetscRandom);
474: /*S
475: MatStencil - Data structure (C struct) for storing information about a single row or
476: column of a matrix as indexed on an associated grid. These are arguments to `MatSetStencil()` and `MatSetBlockStencil()`
478: The i,j, and k represent the logical coordinates over the entire grid (for 2 and 1 dimensional problems the k and j entries are ignored).
479: The c represents the the degrees of freedom at each grid point (the dof argument to `DMDASetDOF()`). If dof is 1 then this entry is ignored.
481: For stencil access to vectors see `DMDAVecGetArray()`, `DMDAVecGetArrayF90()`.
483: Fortran usage is different, see `MatSetValuesStencil()` for details.
485: Level: beginner
487: .seealso: `MatSetValuesStencil()`, `MatSetStencil()`, `MatSetValuesBlockedStencil()`, `DMDAVecGetArray()`, `DMDAVecGetArrayF90()`
488: S*/
489: typedef struct {
490: PetscInt k, j, i, c;
491: } MatStencil;
493: PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat, PetscInt, const MatStencil[], PetscInt, const MatStencil[], const PetscScalar[], InsertMode);
494: PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat, PetscInt, const MatStencil[], PetscInt, const MatStencil[], const PetscScalar[], InsertMode);
495: PETSC_EXTERN PetscErrorCode MatSetStencil(Mat, PetscInt, const PetscInt[], const PetscInt[], PetscInt);
497: /*E
498: MatAssemblyType - Indicates if the matrix is now to be used, for example in a solver, or if you plan
499: to continue to add or insert values to it
501: Level: beginner
503: .seealso: `MatAssemblyBegin()`, `MatAssemblyEnd()`
504: E*/
505: typedef enum {
506: MAT_FLUSH_ASSEMBLY = 1,
507: MAT_FINAL_ASSEMBLY = 0
508: } MatAssemblyType;
509: PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat, MatAssemblyType);
510: PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat, MatAssemblyType);
511: PETSC_EXTERN PetscErrorCode MatAssembled(Mat, PetscBool *);
513: /*E
514: MatOption - Options that may be set for a matrix and its behavior or storage
516: Level: beginner
518: Developer Notes:
519: Entries that are negative need not be called collectively by all processes.
521: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
523: Any additions/changes here must also be made in src/mat/interface/dlregismat.c in MatOptions[]
525: .seealso: `MatSetOption()`
526: E*/
527: typedef enum {
528: MAT_OPTION_MIN = -3,
529: MAT_UNUSED_NONZERO_LOCATION_ERR = -2,
530: MAT_ROW_ORIENTED = -1,
531: MAT_SYMMETRIC = 1,
532: MAT_STRUCTURALLY_SYMMETRIC = 2,
533: MAT_FORCE_DIAGONAL_ENTRIES = 3,
534: MAT_IGNORE_OFF_PROC_ENTRIES = 4,
535: MAT_USE_HASH_TABLE = 5,
536: MAT_KEEP_NONZERO_PATTERN = 6,
537: MAT_IGNORE_ZERO_ENTRIES = 7,
538: MAT_USE_INODES = 8,
539: MAT_HERMITIAN = 9,
540: MAT_SYMMETRY_ETERNAL = 10,
541: MAT_NEW_NONZERO_LOCATION_ERR = 11,
542: MAT_IGNORE_LOWER_TRIANGULAR = 12,
543: MAT_ERROR_LOWER_TRIANGULAR = 13,
544: MAT_GETROW_UPPERTRIANGULAR = 14,
545: MAT_SPD = 15,
546: MAT_NO_OFF_PROC_ZERO_ROWS = 16,
547: MAT_NO_OFF_PROC_ENTRIES = 17,
548: MAT_NEW_NONZERO_LOCATIONS = 18,
549: MAT_NEW_NONZERO_ALLOCATION_ERR = 19,
550: MAT_SUBSET_OFF_PROC_ENTRIES = 20,
551: MAT_SUBMAT_SINGLEIS = 21,
552: MAT_STRUCTURE_ONLY = 22,
553: MAT_SORTED_FULL = 23,
554: MAT_FORM_EXPLICIT_TRANSPOSE = 24,
555: MAT_STRUCTURAL_SYMMETRY_ETERNAL = 25,
556: MAT_SPD_ETERNAL = 26,
557: MAT_OPTION_MAX = 27
558: } MatOption;
560: PETSC_EXTERN const char *const *MatOptions;
561: PETSC_EXTERN PetscErrorCode MatSetOption(Mat, MatOption, PetscBool);
562: PETSC_EXTERN PetscErrorCode MatGetOption(Mat, MatOption, PetscBool *);
563: PETSC_EXTERN PetscErrorCode MatPropagateSymmetryOptions(Mat, Mat);
564: PETSC_EXTERN PetscErrorCode MatGetType(Mat, MatType *);
566: PETSC_EXTERN PetscErrorCode MatGetValues(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscScalar[]);
567: PETSC_EXTERN PetscErrorCode MatGetRow(Mat, PetscInt, PetscInt *, const PetscInt *[], const PetscScalar *[]);
568: PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat, PetscInt, PetscInt *, const PetscInt *[], const PetscScalar *[]);
569: PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
570: PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
571: PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat, Vec, PetscInt);
572: PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat, PetscScalar *[]);
573: PETSC_EXTERN PetscErrorCode MatSeqAIJGetArrayRead(Mat, const PetscScalar *[]);
574: PETSC_EXTERN PetscErrorCode MatSeqAIJGetArrayWrite(Mat, PetscScalar *[]);
575: PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat, PetscScalar *[]);
576: PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArrayRead(Mat, const PetscScalar *[]);
577: PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArrayWrite(Mat, PetscScalar *[]);
578: PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat, PetscInt *);
579: PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
580: PETSC_EXTERN PetscErrorCode MatSeqAIJSetType(Mat, MatType);
581: PETSC_EXTERN PetscErrorCode MatSeqAIJKron(Mat, Mat, MatReuse, Mat *);
582: PETSC_EXTERN PetscErrorCode MatSeqAIJRegister(const char[], PetscErrorCode (*)(Mat, MatType, MatReuse, Mat *));
583: PETSC_EXTERN PetscFunctionList MatSeqAIJList;
584: PETSC_EXTERN PetscErrorCode MatSeqBAIJGetArray(Mat, PetscScalar *[]);
585: PETSC_EXTERN PetscErrorCode MatSeqBAIJRestoreArray(Mat, PetscScalar *[]);
586: PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat, PetscScalar *[]);
587: PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat, PetscScalar *[]);
588: PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat, PetscScalar *[]);
589: PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat, PetscScalar *[]);
590: PETSC_EXTERN PetscErrorCode MatDensePlaceArray(Mat, const PetscScalar[]);
591: PETSC_EXTERN PetscErrorCode MatDenseReplaceArray(Mat, const PetscScalar[]);
592: PETSC_EXTERN PetscErrorCode MatDenseResetArray(Mat);
593: PETSC_EXTERN PetscErrorCode MatDenseGetArrayRead(Mat, const PetscScalar *[]);
594: PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayRead(Mat, const PetscScalar *[]);
595: PETSC_EXTERN PetscErrorCode MatDenseGetArrayWrite(Mat, PetscScalar *[]);
596: PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayWrite(Mat, PetscScalar *[]);
597: PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat, PetscInt *);
598: PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat, PetscInt);
599: PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat, PetscInt *, PetscInt *);
600: PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat, PetscInt, PetscInt);
601: PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat, Mat, Mat);
602: PETSC_EXTERN PetscErrorCode MatSetVariableBlockSizes(Mat, PetscInt, PetscInt *);
603: PETSC_EXTERN PetscErrorCode MatGetVariableBlockSizes(Mat, PetscInt *, const PetscInt **);
605: PETSC_EXTERN PetscErrorCode MatDenseGetColumn(Mat, PetscInt, PetscScalar *[]);
606: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumn(Mat, PetscScalar *[]);
607: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVec(Mat, PetscInt, Vec *);
608: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVec(Mat, PetscInt, Vec *);
609: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecRead(Mat, PetscInt, Vec *);
610: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecRead(Mat, PetscInt, Vec *);
611: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecWrite(Mat, PetscInt, Vec *);
612: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecWrite(Mat, PetscInt, Vec *);
613: PETSC_EXTERN PetscErrorCode MatDenseGetSubMatrix(Mat, PetscInt, PetscInt, PetscInt, PetscInt, Mat *);
614: PETSC_EXTERN PetscErrorCode MatDenseRestoreSubMatrix(Mat, Mat *);
616: PETSC_EXTERN PetscErrorCode MatMult(Mat, Vec, Vec);
617: PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat, Vec, Vec);
618: PETSC_EXTERN PetscErrorCode MatMultAdd(Mat, Vec, Vec, Vec);
619: PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat, Vec, Vec);
620: PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat, Vec, Vec);
621: PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat, Mat, PetscReal, PetscBool *);
622: PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat, Mat, PetscReal, PetscBool *);
623: PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat, Vec, Vec, Vec);
624: PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat, Vec, Vec, Vec);
625: PETSC_EXTERN PetscErrorCode MatMatSolve(Mat, Mat, Mat);
626: PETSC_EXTERN PetscErrorCode MatMatSolveTranspose(Mat, Mat, Mat);
627: PETSC_EXTERN PetscErrorCode MatMatTransposeSolve(Mat, Mat, Mat);
628: PETSC_EXTERN PetscErrorCode MatResidual(Mat, Vec, Vec, Vec);
630: /*E
631: MatDuplicateOption - Indicates if a duplicated sparse matrix should have
632: its numerical values copied over or just its nonzero structure.
634: Level: beginner
636: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
638: $ `MAT_DO_NOT_COPY_VALUES` - Create a matrix using the same nonzero pattern as the original matrix,
639: $ with zeros for the numerical values.
640: $ `MAT_COPY_VALUES` - Create a matrix with the same nonzero pattern as the original matrix
641: $ and with the same numerical values.
642: $ `MAT_SHARE_NONZERO_PATTERN` - Create a matrix that shares the nonzero structure with the previous matrix
643: $ and does not copy it, using zeros for the numerical values. The parent and
644: $ child matrices will share their index (i and j) arrays, and you cannot
645: $ insert new nonzero entries into either matrix.
647: Note:
648: Many matrix types (including `MATSEQAIJ`) do not support the `MAT_SHARE_NONZERO_PATTERN` optimization; in
649: this case the behavior is as if `MAT_DO_NOT_COPY_VALUES` has been specified.
651: .seealso: `MatDuplicate()`
652: E*/
653: typedef enum {
654: MAT_DO_NOT_COPY_VALUES,
655: MAT_COPY_VALUES,
656: MAT_SHARE_NONZERO_PATTERN
657: } MatDuplicateOption;
659: PETSC_EXTERN PetscErrorCode MatConvert(Mat, MatType, MatReuse, Mat *);
660: PETSC_EXTERN PetscErrorCode MatDuplicate(Mat, MatDuplicateOption, Mat *);
662: PETSC_EXTERN PetscErrorCode MatCopy(Mat, Mat, MatStructure);
663: PETSC_EXTERN PetscErrorCode MatView(Mat, PetscViewer);
664: PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat, PetscReal, PetscBool *);
665: PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat, PetscBool *);
666: PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat, PetscReal, PetscBool *);
667: PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat, PetscBool *, PetscBool *);
668: PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat, PetscBool *, PetscBool *);
669: PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetricKnown(Mat, PetscBool *, PetscBool *);
670: PETSC_EXTERN PetscErrorCode MatIsSPDKnown(Mat, PetscBool *, PetscBool *);
671: PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat, PetscBool *, PetscInt *);
672: PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
674: PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
675: PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
676: PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
677: PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
679: /*S
680: MatInfo - Context of matrix information, used with `MatGetInfo()`
682: In Fortran this is simply a double precision array of dimension `MAT_INFO_SIZE`
684: Level: intermediate
686: .seealso: `MatGetInfo()`, `MatInfoType`
687: S*/
688: typedef struct {
689: PetscLogDouble block_size; /* block size */
690: PetscLogDouble nz_allocated, nz_used, nz_unneeded; /* number of nonzeros */
691: PetscLogDouble memory; /* memory allocated */
692: PetscLogDouble assemblies; /* number of matrix assemblies called */
693: PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */
694: PetscLogDouble fill_ratio_given, fill_ratio_needed; /* fill ratio for LU/ILU */
695: PetscLogDouble factor_mallocs; /* number of mallocs during factorization */
696: } MatInfo;
698: /*E
699: MatInfoType - Indicates if you want information about the local part of the matrix,
700: the entire parallel matrix or the maximum over all the local parts.
702: Level: beginner
704: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
706: .seealso: `MatGetInfo()`, `MatInfo`
707: E*/
708: typedef enum {
709: MAT_LOCAL = 1,
710: MAT_GLOBAL_MAX = 2,
711: MAT_GLOBAL_SUM = 3
712: } MatInfoType;
713: PETSC_EXTERN PetscErrorCode MatGetInfo(Mat, MatInfoType, MatInfo *);
714: PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat, Vec);
715: PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat, Vec, PetscInt[]);
716: PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat, Vec, PetscInt[]);
717: PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat, Vec, PetscInt[]);
718: PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat, Vec, PetscInt[]);
719: PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat, Vec);
720: PETSC_EXTERN PetscErrorCode MatTranspose(Mat, MatReuse, Mat *);
721: PETSC_EXTERN PetscErrorCode MatTransposeSymbolic(Mat, Mat *);
722: PETSC_EXTERN PetscErrorCode MatTransposeSetPrecursor(Mat, Mat);
723: PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat, MatReuse, Mat *);
724: PETSC_EXTERN PetscErrorCode MatPermute(Mat, IS, IS, Mat *);
725: PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat, Vec, Vec);
726: PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat, Vec, InsertMode);
728: PETSC_EXTERN PetscErrorCode MatEqual(Mat, Mat, PetscBool *);
729: PETSC_EXTERN PetscErrorCode MatMultEqual(Mat, Mat, PetscInt, PetscBool *);
730: PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat, Mat, PetscInt, PetscBool *);
731: PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat, Mat, PetscInt, PetscBool *);
732: PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat, Mat, PetscInt, PetscBool *);
733: PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeEqual(Mat, Mat, PetscInt, PetscBool *);
734: PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAddEqual(Mat, Mat, PetscInt, PetscBool *);
735: PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat, Mat, Mat, PetscInt, PetscBool *);
736: PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat, Mat, Mat, PetscInt, PetscBool *);
737: PETSC_EXTERN PetscErrorCode MatMatTransposeMultEqual(Mat, Mat, Mat, PetscInt, PetscBool *);
738: PETSC_EXTERN PetscErrorCode MatPtAPMultEqual(Mat, Mat, Mat, PetscInt, PetscBool *);
739: PETSC_EXTERN PetscErrorCode MatRARtMultEqual(Mat, Mat, Mat, PetscInt, PetscBool *);
740: PETSC_EXTERN PetscErrorCode MatIsLinear(Mat, PetscInt, PetscBool *);
742: PETSC_EXTERN PetscErrorCode MatNorm(Mat, NormType, PetscReal *);
743: PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat, NormType, PetscReal *);
744: PETSC_EXTERN PetscErrorCode MatGetColumnSums(Mat, PetscScalar *);
745: PETSC_EXTERN PetscErrorCode MatGetColumnSumsRealPart(Mat, PetscReal *);
746: PETSC_EXTERN PetscErrorCode MatGetColumnSumsImaginaryPart(Mat, PetscReal *);
747: PETSC_EXTERN PetscErrorCode MatGetColumnMeans(Mat, PetscScalar *);
748: PETSC_EXTERN PetscErrorCode MatGetColumnMeansRealPart(Mat, PetscReal *);
749: PETSC_EXTERN PetscErrorCode MatGetColumnMeansImaginaryPart(Mat, PetscReal *);
750: PETSC_EXTERN PetscErrorCode MatGetColumnReductions(Mat, PetscInt, PetscReal *);
751: PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
752: PETSC_EXTERN PetscErrorCode MatSetInf(Mat);
753: PETSC_EXTERN PetscErrorCode MatZeroRows(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec);
754: PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat, IS, PetscScalar, Vec, Vec);
755: PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat, PetscInt, const MatStencil[], PetscScalar, Vec, Vec);
756: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat, PetscInt, const MatStencil[], PetscScalar, Vec, Vec);
757: PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec);
758: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat, IS, PetscScalar, Vec, Vec);
760: PETSC_EXTERN PetscErrorCode MatGetSize(Mat, PetscInt *, PetscInt *);
761: PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat, PetscInt *, PetscInt *);
762: PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat, PetscInt *, PetscInt *);
763: PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat, const PetscInt **);
764: PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat, PetscInt *, PetscInt *);
765: PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat, const PetscInt **);
766: PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat, IS *, IS *);
768: PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *[]);
769: PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatrices() (since version 3.8)") static inline PetscErrorCode MatGetSubMatrices(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
770: {
771: return MatCreateSubMatrices(mat, n, irow, icol, scall, submat);
772: }
773: PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *[]);
774: PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatricesMPI() (since version 3.8)") static inline PetscErrorCode MatGetSubMatricesMPI(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
775: {
776: return MatCreateSubMatricesMPI(mat, n, irow, icol, scall, submat);
777: }
778: PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt, Mat *[]);
779: PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt, Mat *[]);
780: PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat, IS, IS, MatReuse, Mat *);
781: PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatrix() (since version 3.8)") static inline PetscErrorCode MatGetSubMatrix(Mat mat, IS isrow, IS iscol, MatReuse cll, Mat *newmat)
782: {
783: return MatCreateSubMatrix(mat, isrow, iscol, cll, newmat);
784: }
785: PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat, IS, IS, Mat *);
786: PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat, IS, IS, Mat *);
787: PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat, Mat *);
788: PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat *);
790: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm, Mat, PetscInt, PetscInt, MatReuse, Mat *);
791: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm, Mat, PetscInt, PetscInt, Mat *);
792: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat, Mat);
793: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat, MatReuse, Mat *);
794: PETSC_EXTERN PetscErrorCode MatAIJGetLocalMat(Mat, Mat *);
795: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat, MatReuse, IS *, IS *, Mat *);
796: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatMerge(Mat, MatReuse, IS *, Mat *);
797: PETSC_EXTERN PetscErrorCode MatMPIAIJGetNumberNonzeros(Mat, PetscCount *);
798: PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat, Mat, MatReuse, IS *, IS *, Mat *);
799: PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *, const PetscInt *[]);
801: PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat, PetscInt, IS[], PetscInt);
802: PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat, PetscInt n, IS is[], PetscInt ov);
803: PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat, PetscBool);
805: PETSC_EXTERN PetscErrorCode MatMatMult(Mat, Mat, MatReuse, PetscReal, Mat *);
807: PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat, Mat, Mat, MatReuse, PetscReal, Mat *);
808: PETSC_EXTERN PetscErrorCode MatGalerkin(Mat, Mat, Mat, MatReuse, PetscReal, Mat *);
810: PETSC_EXTERN PetscErrorCode MatPtAP(Mat, Mat, MatReuse, PetscReal, Mat *);
811: PETSC_EXTERN PetscErrorCode MatRARt(Mat, Mat, MatReuse, PetscReal, Mat *);
813: PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat, Mat, MatReuse, PetscReal, Mat *);
814: PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat, Mat, MatReuse, PetscReal, Mat *);
816: PETSC_EXTERN PetscErrorCode MatAXPY(Mat, PetscScalar, Mat, MatStructure);
817: PETSC_EXTERN PetscErrorCode MatAYPX(Mat, PetscScalar, Mat, MatStructure);
819: PETSC_EXTERN PetscErrorCode MatScale(Mat, PetscScalar);
820: PETSC_EXTERN PetscErrorCode MatShift(Mat, PetscScalar);
822: PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat, ISLocalToGlobalMapping, ISLocalToGlobalMapping);
823: PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat, ISLocalToGlobalMapping *, ISLocalToGlobalMapping *);
824: PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat, PetscLayout *, PetscLayout *);
825: PETSC_EXTERN PetscErrorCode MatSetLayouts(Mat, PetscLayout, PetscLayout);
826: PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec);
827: PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat, IS, PetscScalar, Vec, Vec);
828: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec);
829: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat, IS, PetscScalar, Vec, Vec);
830: PETSC_EXTERN PetscErrorCode MatGetValuesLocal(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscScalar[]);
831: PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
832: PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
834: PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat, PetscInt, PetscInt);
835: PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
837: PETSC_EXTERN PetscErrorCode MatInterpolate(Mat, Vec, Vec);
838: PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat, Vec, Vec, Vec);
839: PETSC_EXTERN PetscErrorCode MatRestrict(Mat, Vec, Vec);
840: PETSC_EXTERN PetscErrorCode MatMatInterpolate(Mat, Mat, Mat *);
841: PETSC_EXTERN PetscErrorCode MatMatInterpolateAdd(Mat, Mat, Mat, Mat *);
842: PETSC_EXTERN PetscErrorCode MatMatRestrict(Mat, Mat, Mat *);
843: PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat, Vec *, Vec *);
844: PETSC_DEPRECATED_FUNCTION("Use MatCreateVecs() (since version 3.6)") static inline PetscErrorCode MatGetVecs(Mat mat, Vec *x, Vec *y)
845: {
846: return MatCreateVecs(mat, x, y);
847: }
848: PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat, PetscInt, MPI_Comm, MatReuse, Mat *);
849: PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat, MPI_Comm, MatReuse, Mat *);
850: PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat, IS *);
851: PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat, IS *);
852: PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);
854: /*MC
855: MatSetValue - Set a single entry into a matrix.
857: Not collective
859: Synopsis:
860: #include <petscmat.h>
861: PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode)
863: Input Parameters:
864: + m - the matrix
865: . row - the row location of the entry
866: . col - the column location of the entry
867: . value - the value to insert
868: - mode - either `INSERT_VALUES` or `ADD_VALUES`
870: Notes:
871: For efficiency one should use `MatSetValues()` and set several values simultaneously.
873: Level: beginner
875: .seealso: `MatGetValue()`, `MatSetValues()`, `MatSetValueLocal()`, `MatSetValuesLocal()`
876: M*/
877: static inline PetscErrorCode MatSetValue(Mat v, PetscInt i, PetscInt j, PetscScalar va, InsertMode mode)
878: {
879: return MatSetValues(v, 1, &i, 1, &j, &va, mode);
880: }
882: /*@C
883: MatGetValue - Gets a single value from a matrix
885: Not Collective; can only return a value owned by the given process
887: Input Parameters:
888: + mat - the matrix
889: . row - the row location of the entry
890: - col - the column location of the entry
892: Output Parameter:
893: . va - the value
895: Notes:
896: For efficiency one should use `MatGetValues()` and get several values simultaneously.
898: See notes for `MatGetValues()`.
900: Level: advanced
902: .seealso: `MatSetValue()`, `MatGetValueLocal()`, `MatGetValues()`
903: @*/
904: static inline PetscErrorCode MatGetValue(Mat mat, PetscInt row, PetscInt col, PetscScalar *va)
905: {
906: return MatGetValues(mat, 1, &row, 1, &col, va);
907: }
909: /*MC
910: MatSetValueLocal - Inserts or adds a single value into a matrix,
911: using a local numbering of the nodes.
913: Not Collective
915: Input Parameters:
916: + m - the matrix
917: . row - the row location of the entry
918: . col - the column location of the entry
919: . value - the value to insert
920: - mode - either `INSERT_VALUES` or `ADD_VALUES`
922: Notes:
923: For efficiency one should use `MatSetValuesLocal()` and set several values simultaneously.
925: See notes for `MatSetValuesLocal()` for additional information on when and how this function can be used.
927: Level: intermediate
929: .seealso: `MatSetValue()`, `MatSetValuesLocal()`
930: M*/
931: static inline PetscErrorCode MatSetValueLocal(Mat v, PetscInt i, PetscInt j, PetscScalar va, InsertMode mode)
932: {
933: return MatSetValuesLocal(v, 1, &i, 1, &j, &va, mode);
934: }
936: /*MC
937: MatPreallocateBegin - Begins the block of code that will count the number of nonzeros per
938: row in a matrix providing the data that one can use to correctly preallocate the matrix.
940: Synopsis:
941: #include <petscmat.h>
942: PetscErrorCode MatPreallocateBegin(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
944: Collective
946: Input Parameters:
947: + comm - the communicator that will share the eventually allocated matrix
948: . nrows - the number of LOCAL rows in the matrix
949: - ncols - the number of LOCAL columns in the matrix
951: Output Parameters:
952: + dnz - the array that will be passed to the matrix preallocation routines
953: - onz - the other array passed to the matrix preallocation routines
955: Level: intermediate
957: Notes:
958: This is a macro that handles its own error checking, it does not return an error code.
960: See Users-Manual: ch_performance for more details.
962: Do not malloc or free dnz and onz, that is handled internally by these routines
964: Developer Notes:
965: This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
967: .seealso: `MatPreallocateEnd()`, `MatPreallocateSet()`, `MatPreallocateSymmetricSetBlock()`, `MatPreallocateSetLocal()`,
968: `MatPreallocateSymmetricSetLocalBlock()`
969: M*/
970: #define MatPreallocateBegin(comm, nrows, ncols, dnz, onz) \
971: do { \
972: PetscInt __nrows = (nrows), __ncols = (ncols), __rstart, __start, __end = 0; \
973: PetscCalloc2(__nrows, &(dnz), __nrows, &(onz)); \
974: MPI_Scan(&__ncols, &__end, 1, MPIU_INT, MPI_SUM, comm); \
975: __start = __end - __ncols; \
976: (void)__start; \
977: MPI_Scan(&__nrows, &__rstart, 1, MPIU_INT, MPI_SUM, comm); \
978: __rstart -= __nrows
980: #define MatPreallocateInitialize(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use MatPreallocateBegin() (since version 3.18)\"") MatPreallocateBegin(__VA_ARGS__)
982: /*MC
983: MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
984: inserted using a local number of the rows and columns
986: Synopsis:
987: #include <petscmat.h>
988: PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
990: Not Collective
992: Input Parameters:
993: + map - the row mapping from local numbering to global numbering
994: . nrows - the number of rows indicated
995: . rows - the indices of the rows
996: . cmap - the column mapping from local to global numbering
997: . ncols - the number of columns in the matrix
998: . cols - the columns indicated
999: . dnz - the array that will be passed to the matrix preallocation routines
1000: - onz - the other array passed to the matrix preallocation routines
1002: Level: intermediate
1004: Notes:
1005: See Users-Manual: ch_performance for more details.
1007: Do not malloc or free dnz and onz, that is handled internally by these routines
1009: .seealso: `MatPreallocateEnd()`, `MatPreallocateSet()`, `MatPreallocateSymmetricSetBlock()`
1010: `MatPreallocateBegin()`, `MatPreallocateSymmetricSetLocalBlock()`, `MatPreallocateSetLocalRemoveDups()`
1011: M*/
1012: #define MatPreallocateSetLocal(rmap, nrows, rows, cmap, ncols, cols, dnz, onz) \
1013: PetscMacroReturnStandard(ISLocalToGlobalMappingApply(rmap, nrows, rows, rows)); PetscCall(ISLocalToGlobalMappingApply(cmap, ncols, cols, cols)); for (PetscInt __l = 0; __l < nrows; __l++) PetscCall(MatPreallocateSet((rows)[__l], ncols, cols, dnz, onz);)
1015: /*MC
1016: MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be
1017: inserted using a local number of the rows and columns. This version removes any duplicate columns in cols
1019: Synopsis:
1020: #include <petscmat.h>
1021: PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
1023: Not Collective
1025: Input Parameters:
1026: + map - the row mapping from local numbering to global numbering
1027: . nrows - the number of rows indicated
1028: . rows - the indices of the rows (these values are mapped to the global values)
1029: . cmap - the column mapping from local to global numbering
1030: . ncols - the number of columns in the matrix (this value will be changed if duplicate columns are found)
1031: . cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed)
1032: . dnz - the array that will be passed to the matrix preallocation routines
1033: - onz - the other array passed to the matrix preallocation routines
1035: Level: intermediate
1037: Notes:
1038: See Users-Manual: ch_performance for more details.
1040: Do not malloc or free dnz and onz, that is handled internally by these routines
1042: .seealso: `MatPreallocateEnd()`, `MatPreallocateSet()`, `MatPreallocateSymmetricSetBlock()`
1043: `MatPreallocateBegin()`, `MatPreallocateSymmetricSetLocalBlock()`, `MatPreallocateSetLocal()`
1044: M*/
1045: #define MatPreallocateSetLocalRemoveDups(rmap, nrows, rows, cmap, ncols, cols, dnz, onz) \
1046: PetscMacroReturnStandard(ISLocalToGlobalMappingApply(rmap, nrows, rows, rows)); PetscCall(ISLocalToGlobalMappingApply(cmap, ncols, cols, cols)); PetscCall(PetscSortRemoveDupsInt(&ncols, cols)); for (PetscInt __l = 0; __l < nrows; __l++) PetscCall(MatPreallocateSet((rows)[__l], ncols, cols, dnz, onz);)
1048: /*MC
1049: MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
1050: inserted using a local number of the rows and columns
1052: Synopsis:
1053: #include <petscmat.h>
1054: PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
1056: Not Collective
1058: Input Parameters:
1059: + map - the row mapping from local numbering to global numbering
1060: . nrows - the number of rows indicated
1061: . rows - the indices of the rows
1062: . cmap - the column mapping from local to global numbering
1063: . ncols - the number of columns in the matrix
1064: . cols - the columns indicated
1065: . dnz - the array that will be passed to the matrix preallocation routines
1066: - onz - the other array passed to the matrix preallocation routines
1068: Level: intermediate
1070: Notes:
1071: See Users-Manual: ch_performance for more details.
1073: Do not malloc or free dnz and onz, that is handled internally by these routines
1075: .seealso: `MatPreallocateEnd()`, `MatPreallocateSet()`, `MatPreallocateSymmetricSetBlock()`
1076: `MatPreallocateBegin()`, `MatPreallocateSymmetricSetLocalBlock()`
1077: M*/
1078: #define MatPreallocateSetLocalBlock(rmap, nrows, rows, cmap, ncols, cols, dnz, onz) \
1079: PetscMacroReturnStandard(ISLocalToGlobalMappingApplyBlock(rmap, nrows, rows, rows)); PetscCall(ISLocalToGlobalMappingApplyBlock(cmap, ncols, cols, cols)); for (PetscInt __l = 0; __l < nrows; __l++) PetscCall(MatPreallocateSet((rows)[__l], ncols, cols, dnz, onz);)
1081: /*MC
1082: MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
1083: inserted using a local number of the rows and columns
1085: Synopsis:
1086: #include <petscmat.h>
1087: PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
1089: Not Collective
1091: Input Parameters:
1092: + map - the mapping between local numbering and global numbering
1093: . nrows - the number of rows indicated
1094: . rows - the indices of the rows
1095: . ncols - the number of columns in the matrix
1096: . cols - the columns indicated
1097: . dnz - the array that will be passed to the matrix preallocation routines
1098: - onz - the other array passed to the matrix preallocation routines
1100: Level: intermediate
1102: Notes:
1103: See Users-Manual: ch_performance for more details.
1105: Do not malloc or free dnz and onz that is handled internally by these routines
1107: .seealso: `MatPreallocateEnd()`, `MatPreallocateSet()`
1108: `MatPreallocateBegin()`, `MatPreallocateSetLocal()`
1109: M*/
1110: #define MatPreallocateSymmetricSetLocalBlock(map, nrows, rows, ncols, cols, dnz, onz) \
1111: PetscMacroReturnStandard(ISLocalToGlobalMappingApplyBlock(map, nrows, rows, rows)); PetscCall(ISLocalToGlobalMappingApplyBlock(map, ncols, cols, cols)); for (PetscInt __l = 0; __l < nrows; __l++) PetscCall(MatPreallocateSymmetricSetBlock((rows)[__l], ncols, cols, dnz, onz);)
1113: /*MC
1114: MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
1115: inserted using a local number of the rows and columns
1117: Synopsis:
1118: #include <petscmat.h>
1119: PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
1121: Not Collective
1123: Input Parameters:
1124: + row - the row
1125: . ncols - the number of columns in the matrix
1126: - cols - the columns indicated
1128: Output Parameters:
1129: + dnz - the array that will be passed to the matrix preallocation routines
1130: - onz - the other array passed to the matrix preallocation routines
1132: Level: intermediate
1134: Notes:
1135: See Users-Manual: ch_performance for more details.
1137: Do not malloc or free dnz and onz that is handled internally by these routines
1139: This is a MACRO not a function because it uses variables declared in MatPreallocateBegin().
1141: .seealso: `MatPreallocateEnd()`, `MatPreallocateSet()`, `MatPreallocateSymmetricSetBlock()`
1142: `MatPreallocateBegin()`, `MatPreallocateSetLocal()`
1143: M*/
1144: #define MatPreallocateSet(row, nc, cols, dnz, onz) \
1146: if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
1147: else if (dnz[row - __rstart] < __ncols) dnz[row - __rstart]++; \
1148: })
1150: /*MC
1151: MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
1152: inserted using a local number of the rows and columns
1154: Synopsis:
1155: #include <petscmat.h>
1156: PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
1158: Not Collective
1160: Input Parameters:
1161: + nrows - the number of rows indicated
1162: . rows - the indices of the rows
1163: . ncols - the number of columns in the matrix
1164: . cols - the columns indicated
1165: . dnz - the array that will be passed to the matrix preallocation routines
1166: - onz - the other array passed to the matrix preallocation routines
1168: Level: intermediate
1170: Notes:
1171: See Users-Manual: ch_performance for more details.
1173: Do not malloc or free dnz and onz that is handled internally by these routines
1175: This is a MACRO not a function because it uses variables declared in MatPreallocateBegin().
1177: .seealso: `MatPreallocateEnd()`, `MatPreallocateSet()`, `MatPreallocateBegin()`,
1178: `MatPreallocateSymmetricSetLocalBlock()`, `MatPreallocateSetLocal()`
1179: M*/
1180: #define MatPreallocateSymmetricSetBlock(row, nc, cols, dnz, onz) \
1181: PetscMacroReturnStandard(for (PetscInt __i = 0; __i < nc; __i++) { \
1182: if (cols[__i] >= __end) onz[row - __rstart]++; \
1183: else if (cols[__i] >= row && dnz[row - __rstart] < __ncols) dnz[row - __rstart]++; \
1184: })
1186: /*MC
1187: MatPreallocateLocation - An alternative to MatPreallocateSet() that puts the nonzero locations into the matrix if it exists
1189: Synopsis:
1190: #include <petscmat.h>
1191: PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
1193: Not Collective
1195: Input Parameters:
1196: + A - matrix
1197: . row - row where values exist (must be local to this process)
1198: . ncols - number of columns
1199: . cols - columns with nonzeros
1200: . dnz - the array that will be passed to the matrix preallocation routines
1201: - onz - the other array passed to the matrix preallocation routines
1203: Level: intermediate
1205: Notes:
1206: See Users-Manual: ch_performance for more details.
1208: Do not malloc or free dnz and onz that is handled internally by these routines
1210: This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
1212: .seealso: `MatPreallocateBegin()`, `MatPreallocateSet()`, `MatPreallocateSymmetricSetBlock()`, `MatPreallocateSetLocal()`,
1213: `MatPreallocateSymmetricSetLocalBlock()`
1214: M*/
1215: #define MatPreallocateLocation(A, row, ncols, cols, dnz, onz) (A ? MatSetValues(A, 1, &row, ncols, cols, NULL, INSERT_VALUES) : MatPreallocateSet(row, ncols, cols, dnz, onz))
1217: /*MC
1218: MatPreallocateEnd - Ends the block of code that will count the number of nonzeros per
1219: row in a matrix providing the data that one can use to correctly preallocate the matrix.
1221: Synopsis:
1222: #include <petscmat.h>
1223: PetscErrorCode MatPreallocateEnd(PetscInt *dnz, PetscInt *onz)
1225: Collective
1227: Input Parameters:
1228: + dnz - the array that was be passed to the matrix preallocation routines
1229: - onz - the other array passed to the matrix preallocation routines
1231: Level: intermediate
1233: Notes:
1234: This is a macro that handles its own error checking, it does not return an error code.
1236: See Users-Manual: ch_performance for more details.
1238: Do not malloc or free dnz and onz, that is handled internally by these routines
1240: Developer Notes:
1241: This is a MACRO not a function because it closes the { started in MatPreallocateBegin().
1243: .seealso: `MatPreallocateBegin()`, `MatPreallocateSet()`, `MatPreallocateSymmetricSetBlock()`, `MatPreallocateSetLocal()`,
1244: `MatPreallocateSymmetricSetLocalBlock()`
1245: M*/
1246: #define MatPreallocateEnd(dnz, onz) \
1247: PetscFree2(dnz, onz); \
1248: } \
1249: while (0)
1251: #define MatPreallocateFinalize(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use MatPreallocateEnd() (since version 3.18)\"") MatPreallocateEnd(__VA_ARGS__)
1253: /* Routines unique to particular data structures */
1254: PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat, void *);
1256: PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat, IS *, IS *);
1257: PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat, PetscInt *, PetscInt *[], PetscInt *);
1259: PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat, PetscInt[]);
1260: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat, PetscInt[]);
1261: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], Mat *);
1262: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], Mat *);
1263: PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], Mat *);
1264: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], Mat *, PetscInt, PetscBool);
1266: #define MAT_SKIP_ALLOCATION -4
1268: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat, PetscInt, PetscInt, const PetscInt[]);
1269: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat, PetscInt, PetscInt, const PetscInt[]);
1270: PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat, PetscInt, const PetscInt[]);
1271: PETSC_EXTERN PetscErrorCode MatSeqAIJSetTotalPreallocation(Mat, PetscInt);
1273: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[]);
1274: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[]);
1275: PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[]);
1276: PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat, const PetscInt[], const PetscInt[], const PetscScalar[]);
1277: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[]);
1278: PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat, const PetscInt[], const PetscInt[], const PetscScalar[]);
1279: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[]);
1280: PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat, PetscInt[], PetscInt[], PetscInt[]);
1281: PETSC_EXTERN PetscErrorCode MatMPIAdjToSeq(Mat, Mat *);
1282: PETSC_EXTERN PetscErrorCode MatMPIAdjToSeqRankZero(Mat, Mat *);
1283: PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat, PetscScalar[]);
1284: PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat, PetscScalar[]);
1285: PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat, Mat *, Mat *, const PetscInt *[]);
1286: PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat, Mat *, Mat *, const PetscInt *[]);
1287: PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat, Mat *);
1289: PETSC_EXTERN PetscErrorCode MatDenseGetLDA(Mat, PetscInt *);
1290: PETSC_EXTERN PetscErrorCode MatDenseSetLDA(Mat, PetscInt);
1291: PETSC_DEPRECATED_FUNCTION("Use MatDenseSetLDA() (since version 3.14)") static inline PetscErrorCode MatSeqDenseSetLDA(Mat A, PetscInt lda)
1292: {
1293: return MatDenseSetLDA(A, lda);
1294: }
1295: PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat, Mat *);
1297: PETSC_EXTERN PetscErrorCode MatBlockMatSetPreallocation(Mat, PetscInt, PetscInt, const PetscInt[]);
1299: PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
1300: PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
1302: PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat, IS *);
1303: PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat, IS *);
1304: /*
1305: These routines are not usually accessed directly, rather solving is
1306: done through the KSP and PC interfaces.
1307: */
1309: /*J
1310: MatOrderingType - String with the name of a PETSc matrix ordering
1312: Level: beginner
1314: Notes:
1315: If `MATORDERINGEXTERNAL` is used then PETSc does not compute an ordering and utilizes one built into the factorization package
1317: .seealso: `MatGetOrdering()`
1318: J*/
1319: typedef const char *MatOrderingType;
1320: #define MATORDERINGNATURAL "natural"
1321: #define MATORDERINGND "nd"
1322: #define MATORDERING1WD "1wd"
1323: #define MATORDERINGRCM "rcm"
1324: #define MATORDERINGQMD "qmd"
1325: #define MATORDERINGROWLENGTH "rowlength"
1326: #define MATORDERINGWBM "wbm"
1327: #define MATORDERINGSPECTRAL "spectral"
1328: #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */
1329: #define MATORDERINGMETISND "metisnd" /* only works if METIS is installed with PETSc */
1330: #define MATORDERINGNATURAL_OR_ND "natural_or_nd" /* special coase used for Cholesky and ICC, allows ND when AIJ matrix is used but Natural when SBAIJ is used */
1331: #define MATORDERINGEXTERNAL "external" /* uses an ordering type internal to the factorization package */
1333: PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat, MatOrderingType, IS *, IS *);
1334: PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList *);
1335: PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[], PetscErrorCode (*)(Mat, MatOrderingType, IS *, IS *));
1336: PETSC_EXTERN PetscFunctionList MatOrderingList;
1338: #include "petscmatcoarsen.h"
1340: PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat, PetscReal, IS, IS);
1341: PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat, PetscReal, PetscBool, Mat *);
1343: PETSC_EXTERN PetscErrorCode MatFactorGetPreferredOrdering(Mat, MatFactorType, MatOrderingType *);
1345: /*S
1346: MatFactorShiftType - Numeric Shift for factorizations
1348: Level: beginner
1350: .seealso: `MatGetFactor()`
1351: S*/
1352: typedef enum {
1353: MAT_SHIFT_NONE,
1354: MAT_SHIFT_NONZERO,
1355: MAT_SHIFT_POSITIVE_DEFINITE,
1356: MAT_SHIFT_INBLOCKS
1357: } MatFactorShiftType;
1358: PETSC_EXTERN const char *const MatFactorShiftTypes[];
1359: PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];
1361: /*S
1362: MatFactorError - indicates what type of error was generated in a matrix factorization
1364: Level: beginner
1366: Developer Note:
1367: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1369: .seealso: `MatGetFactor()`
1370: S*/
1371: typedef enum {
1372: MAT_FACTOR_NOERROR,
1373: MAT_FACTOR_STRUCT_ZEROPIVOT,
1374: MAT_FACTOR_NUMERIC_ZEROPIVOT,
1375: MAT_FACTOR_OUTMEMORY,
1376: MAT_FACTOR_OTHER
1377: } MatFactorError;
1379: PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat, MatFactorError *);
1380: PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat);
1381: PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat, PetscReal *, PetscInt *);
1383: /*S
1384: MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization
1386: In Fortran these are simply double precision arrays of size `MAT_FACTORINFO_SIZE`, that is use
1387: $ MatFactorInfo info(MAT_FACTORINFO_SIZE)
1389: Notes:
1390: These are not usually directly used by users, instead use `PC` type of `PCLU`, `PCILU`, `PCCHOLESKY` or `PCICC`.
1392: You can use `MatFactorInfoInitialize()` to set default values.
1394: Level: developer
1396: .seealso: `MatLUFactorSymbolic()`, `MatILUFactorSymbolic()`, `MatCholeskyFactorSymbolic()`, `MatICCFactorSymbolic()`, `MatICCFactor()`,
1397: `MatFactorInfoInitialize()`
1399: S*/
1400: typedef struct {
1401: PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */
1402: PetscReal usedt;
1403: PetscReal dt; /* drop tolerance */
1404: PetscReal dtcol; /* tolerance for pivoting */
1405: PetscReal dtcount; /* maximum nonzeros to be allowed per row */
1406: PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1407: PetscReal levels; /* ICC/ILU(levels) */
1408: PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1409: factorization may be faster if do not pivot */
1410: PetscReal zeropivot; /* pivot is called zero if less than this */
1411: PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */
1412: PetscReal shiftamount; /* how large the shift is */
1413: } MatFactorInfo;
1415: PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo *);
1416: PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat, IS, const MatFactorInfo *);
1417: PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat, Mat, IS, const MatFactorInfo *);
1418: PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat, Mat, const MatFactorInfo *);
1419: PETSC_EXTERN PetscErrorCode MatLUFactor(Mat, IS, IS, const MatFactorInfo *);
1420: PETSC_EXTERN PetscErrorCode MatILUFactor(Mat, IS, IS, const MatFactorInfo *);
1421: PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat, Mat, IS, IS, const MatFactorInfo *);
1422: PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat, Mat, IS, IS, const MatFactorInfo *);
1423: PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat, Mat, IS, const MatFactorInfo *);
1424: PETSC_EXTERN PetscErrorCode MatICCFactor(Mat, IS, const MatFactorInfo *);
1425: PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat, Mat, const MatFactorInfo *);
1426: PETSC_EXTERN PetscErrorCode MatQRFactor(Mat, IS, const MatFactorInfo *);
1427: PETSC_EXTERN PetscErrorCode MatQRFactorSymbolic(Mat, Mat, IS, const MatFactorInfo *);
1428: PETSC_EXTERN PetscErrorCode MatQRFactorNumeric(Mat, Mat, const MatFactorInfo *);
1429: PETSC_EXTERN PetscErrorCode MatGetInertia(Mat, PetscInt *, PetscInt *, PetscInt *);
1430: PETSC_EXTERN PetscErrorCode MatSolve(Mat, Vec, Vec);
1431: PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat, Vec, Vec);
1432: PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat, Vec, Vec);
1433: PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat, Vec, Vec, Vec);
1434: PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat, Vec, Vec);
1435: PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat, Vec, Vec, Vec);
1436: PETSC_EXTERN PetscErrorCode MatSolves(Mat, Vecs, Vecs);
1437: PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
1439: typedef enum {
1440: MAT_FACTOR_SCHUR_UNFACTORED,
1441: MAT_FACTOR_SCHUR_FACTORED,
1442: MAT_FACTOR_SCHUR_INVERTED
1443: } MatFactorSchurStatus;
1444: PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat, IS);
1445: PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat, Mat *, MatFactorSchurStatus *);
1446: PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat, Mat *, MatFactorSchurStatus);
1447: PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat);
1448: PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat, Mat *, MatFactorSchurStatus *);
1449: PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat, Vec, Vec);
1450: PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat, Vec, Vec);
1451: PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat);
1453: PETSC_EXTERN PetscErrorCode MatSeqDenseInvert(Mat);
1454: /*E
1455: MatSORType - What type of (S)SOR to perform
1457: Level: beginner
1459: May be bitwise ORd together
1461: Developer Notes:
1462: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1464: `MatSORType` may be bitwise ORd together, so do not change the numerical values
1466: .seealso: `MatSOR()`
1467: E*/
1468: typedef enum {
1469: SOR_FORWARD_SWEEP = 1,
1470: SOR_BACKWARD_SWEEP = 2,
1471: SOR_SYMMETRIC_SWEEP = 3,
1472: SOR_LOCAL_FORWARD_SWEEP = 4,
1473: SOR_LOCAL_BACKWARD_SWEEP = 8,
1474: SOR_LOCAL_SYMMETRIC_SWEEP = 12,
1475: SOR_ZERO_INITIAL_GUESS = 16,
1476: SOR_EISENSTAT = 32,
1477: SOR_APPLY_UPPER = 64,
1478: SOR_APPLY_LOWER = 128
1479: } MatSORType;
1480: PETSC_EXTERN PetscErrorCode MatSOR(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);
1482: /*S
1483: MatColoring - Object for managing the coloring of matrices.
1485: Level: beginner
1487: Notes:
1488: Coloring of matrices can be computed directly from the sparse matrix nonzero structure via the `MatColoring` object or from the mesh from which the
1489: matrix comes from via `DMCreateColoring()`. In general using the mesh produces a more optimal coloring (fewer colors).
1491: Once a coloring is available `MatFDColoringCreate()` creates an object that can be used to efficiently compute Jacobians using that coloring. This
1492: same object can also be used to efficiently convert data created by Automatic Differentation tools to PETSc sparse matrices.
1494: .seealso: `MatFDColoringCreate()`, `MatColoringWeightType`, `ISColoring`, `MatFDColoring`, `DMCreateColoring()`, `MatColoringCreate()`, `MatOrdering`, `MatPartitioning`, `MatColoringType`
1495: S*/
1496: typedef struct _p_MatColoring *MatColoring;
1498: /*J
1499: MatColoringType - String with the name of a PETSc matrix coloring
1501: Level: beginner
1503: .seealso: `MatColoringSetType()`, `MatColoring`
1504: J*/
1505: typedef const char *MatColoringType;
1506: #define MATCOLORINGJP "jp"
1507: #define MATCOLORINGPOWER "power"
1508: #define MATCOLORINGNATURAL "natural"
1509: #define MATCOLORINGSL "sl"
1510: #define MATCOLORINGLF "lf"
1511: #define MATCOLORINGID "id"
1512: #define MATCOLORINGGREEDY "greedy"
1514: /*E
1515: MatColoringWeightType - Type of weight scheme
1517: Not Collective
1519: + `MAT_COLORING_RANDOM` - Random weights
1520: . `MAT_COLORING_LEXICAL` - Lexical weighting based upon global numbering.
1521: - `MAT_COLORING_LF` - Last-first weighting.
1523: Level: intermediate
1525: Developer Note:
1526: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1528: .seealso: `MatColoring`, `MatColoringCreate()`
1529: E*/
1530: typedef enum {
1531: MAT_COLORING_WEIGHT_RANDOM,
1532: MAT_COLORING_WEIGHT_LEXICAL,
1533: MAT_COLORING_WEIGHT_LF,
1534: MAT_COLORING_WEIGHT_SL
1535: } MatColoringWeightType;
1537: PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat, MatColoring *);
1538: PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat, PetscInt, PetscInt *);
1539: PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring *);
1540: PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring, PetscViewer);
1541: PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring, MatColoringType);
1542: PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1543: PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring, PetscInt);
1544: PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring, PetscInt *);
1545: PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring, PetscInt);
1546: PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring, PetscInt *);
1547: PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring, ISColoring *);
1548: PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[], PetscErrorCode (*)(MatColoring));
1549: PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat, PetscInt, PetscInt, ISColoringValue[], ISColoring *);
1550: PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring, MatColoringWeightType);
1551: PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring, PetscReal *, PetscInt *);
1552: PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring, PetscReal **, PetscInt **lperm);
1553: PETSC_EXTERN PetscErrorCode MatColoringTest(MatColoring, ISColoring);
1554: PETSC_DEPRECATED_FUNCTION("Use MatColoringTest() (since version 3.10)") static inline PetscErrorCode MatColoringTestValid(MatColoring matcoloring, ISColoring iscoloring)
1555: {
1556: return MatColoringTest(matcoloring, iscoloring);
1557: }
1558: PETSC_EXTERN PetscErrorCode MatISColoringTest(Mat, ISColoring);
1560: /*S
1561: MatFDColoring - Object for computing a sparse Jacobian via finite differences with coloring
1563: Level: beginner
1565: Notes:
1566: This object is creating utilizing a coloring provided by the `MatColoring` object or `DMCreateColoring()`
1568: .seealso: `MatFDColoringCreate()`, `MatFDColoringSetFunction()`, `MatColoring`, `DMCreateColoring()`
1569: S*/
1570: typedef struct _p_MatFDColoring *MatFDColoring;
1572: PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat, ISColoring, MatFDColoring *);
1573: PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring *);
1574: PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring, PetscViewer);
1575: PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring, PetscErrorCode (*)(void), void *);
1576: PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring, PetscErrorCode (**)(void), void **);
1577: PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring, PetscReal, PetscReal);
1578: PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1579: PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat, MatFDColoring, Vec, void *);
1580: PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring, Vec);
1581: PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring, PetscInt *, const PetscInt *[]);
1582: PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat, ISColoring, MatFDColoring);
1583: PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring, PetscInt, PetscInt);
1584: PETSC_EXTERN PetscErrorCode MatFDColoringSetValues(Mat, MatFDColoring, const PetscScalar *);
1586: /*S
1587: MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1589: Level: beginner
1591: .seealso: `MatTransposeColoringCreate()`
1592: S*/
1593: typedef struct _p_MatTransposeColoring *MatTransposeColoring;
1595: PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat, ISColoring, MatTransposeColoring *);
1596: PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring, Mat, Mat);
1597: PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring, Mat, Mat);
1598: PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring *);
1600: /*S
1601: MatPartitioning - Object for managing the partitioning of a matrix or graph
1603: Level: beginner
1605: Note:
1606: There is also a `PetscPartitioner` object that provides the same functionality. It can utilize the `MatPartitioning` operations
1607: via `PetscPartitionerSetType`(p,`PETSCPARTITIONERMATPARTITIONING`)
1609: Developers Note:
1610: It is an extra maintainance and documentation cost to have two objects with the same functionality.
1612: .seealso: `MatPartitioningCreate()`, `MatPartitioningType`, `MatColoring`, `MatGetOrdering()`
1613: S*/
1614: typedef struct _p_MatPartitioning *MatPartitioning;
1616: /*J
1617: MatPartitioningType - String with the name of a PETSc matrix partitioning
1619: Level: beginner
1620: dm
1621: .seealso: `MatPartitioningCreate()`, `MatPartitioning`, `MatPartitioningSetType()`
1622: J*/
1623: typedef const char *MatPartitioningType;
1624: #define MATPARTITIONINGCURRENT "current"
1625: #define MATPARTITIONINGAVERAGE "average"
1626: #define MATPARTITIONINGSQUARE "square"
1627: #define MATPARTITIONINGPARMETIS "parmetis"
1628: #define MATPARTITIONINGCHACO "chaco"
1629: #define MATPARTITIONINGPARTY "party"
1630: #define MATPARTITIONINGPTSCOTCH "ptscotch"
1631: #define MATPARTITIONINGHIERARCH "hierarch"
1633: PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm, MatPartitioning *);
1634: PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning, MatPartitioningType);
1635: PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning, PetscInt);
1636: PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning, Mat);
1637: PETSC_EXTERN PetscErrorCode MatPartitioningSetNumberVertexWeights(MatPartitioning, PetscInt);
1638: PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning, const PetscInt[]);
1639: PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning, const PetscReal[]);
1640: PETSC_EXTERN PetscErrorCode MatPartitioningSetUseEdgeWeights(MatPartitioning, PetscBool);
1641: PETSC_EXTERN PetscErrorCode MatPartitioningGetUseEdgeWeights(MatPartitioning, PetscBool *);
1642: PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning, IS *);
1643: PETSC_EXTERN PetscErrorCode MatPartitioningImprove(MatPartitioning, IS *);
1644: PETSC_EXTERN PetscErrorCode MatPartitioningViewImbalance(MatPartitioning, IS);
1645: PETSC_EXTERN PetscErrorCode MatPartitioningApplyND(MatPartitioning, IS *);
1646: PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning *);
1647: PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[], PetscErrorCode (*)(MatPartitioning));
1648: PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning, PetscViewer);
1649: PETSC_EXTERN PetscErrorCode MatPartitioningViewFromOptions(MatPartitioning, PetscObject, const char[]);
1650: PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
1651: PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning, MatPartitioningType *);
1653: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1654: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1655: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
1657: typedef enum {
1658: MP_CHACO_MULTILEVEL = 1,
1659: MP_CHACO_SPECTRAL = 2,
1660: MP_CHACO_LINEAR = 4,
1661: MP_CHACO_RANDOM = 5,
1662: MP_CHACO_SCATTERED = 6
1663: } MPChacoGlobalType;
1664: PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1665: typedef enum {
1666: MP_CHACO_KERNIGHAN = 1,
1667: MP_CHACO_NONE = 2
1668: } MPChacoLocalType;
1669: PETSC_EXTERN const char *const MPChacoLocalTypes[];
1670: typedef enum {
1671: MP_CHACO_LANCZOS = 0,
1672: MP_CHACO_RQI = 1
1673: } MPChacoEigenType;
1674: PETSC_EXTERN const char *const MPChacoEigenTypes[];
1676: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
1677: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning, MPChacoGlobalType *);
1678: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1679: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning, MPChacoLocalType *);
1680: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning, PetscReal);
1681: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning, MPChacoEigenType);
1682: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning, MPChacoEigenType *);
1683: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1684: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning, PetscReal *);
1685: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
1686: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning, PetscInt *);
1688: #define MP_PARTY_OPT "opt"
1689: #define MP_PARTY_LIN "lin"
1690: #define MP_PARTY_SCA "sca"
1691: #define MP_PARTY_RAN "ran"
1692: #define MP_PARTY_GBF "gbf"
1693: #define MP_PARTY_GCF "gcf"
1694: #define MP_PARTY_BUB "bub"
1695: #define MP_PARTY_DEF "def"
1696: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning, const char *);
1697: #define MP_PARTY_HELPFUL_SETS "hs"
1698: #define MP_PARTY_KERNIGHAN_LIN "kl"
1699: #define MP_PARTY_NONE "no"
1700: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning, const char *);
1701: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning, PetscReal);
1702: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning, PetscBool);
1703: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning, PetscBool);
1705: typedef enum {
1706: MP_PTSCOTCH_DEFAULT,
1707: MP_PTSCOTCH_QUALITY,
1708: MP_PTSCOTCH_SPEED,
1709: MP_PTSCOTCH_BALANCE,
1710: MP_PTSCOTCH_SAFETY,
1711: MP_PTSCOTCH_SCALABILITY
1712: } MPPTScotchStrategyType;
1713: PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1715: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning, PetscReal);
1716: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning, PetscReal *);
1717: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning, MPPTScotchStrategyType);
1718: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning, MPPTScotchStrategyType *);
1720: /*
1721: * hierarchical partitioning
1722: */
1723: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning, IS *);
1724: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning, IS *);
1725: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning, PetscInt);
1726: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);
1728: PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat, PetscInt, Mat *);
1730: /*
1731: If you add entries here you must also add them to src/mat/f90-mod/petscmat.h
1732: If any of the enum values are changed, also update dMatOps dict at src/binding/petsc4py/src/libpetsc4py/libpetsc4py.pyx
1733: */
1734: typedef enum {
1735: MATOP_SET_VALUES = 0,
1736: MATOP_GET_ROW = 1,
1737: MATOP_RESTORE_ROW = 2,
1738: MATOP_MULT = 3,
1739: MATOP_MULT_ADD = 4,
1740: MATOP_MULT_TRANSPOSE = 5,
1741: MATOP_MULT_TRANSPOSE_ADD = 6,
1742: MATOP_SOLVE = 7,
1743: MATOP_SOLVE_ADD = 8,
1744: MATOP_SOLVE_TRANSPOSE = 9,
1745: MATOP_SOLVE_TRANSPOSE_ADD = 10,
1746: MATOP_LUFACTOR = 11,
1747: MATOP_CHOLESKYFACTOR = 12,
1748: MATOP_SOR = 13,
1749: MATOP_TRANSPOSE = 14,
1750: MATOP_GETINFO = 15,
1751: MATOP_EQUAL = 16,
1752: MATOP_GET_DIAGONAL = 17,
1753: MATOP_DIAGONAL_SCALE = 18,
1754: MATOP_NORM = 19,
1755: MATOP_ASSEMBLY_BEGIN = 20,
1756: MATOP_ASSEMBLY_END = 21,
1757: MATOP_SET_OPTION = 22,
1758: MATOP_ZERO_ENTRIES = 23,
1759: MATOP_ZERO_ROWS = 24,
1760: MATOP_LUFACTOR_SYMBOLIC = 25,
1761: MATOP_LUFACTOR_NUMERIC = 26,
1762: MATOP_CHOLESKY_FACTOR_SYMBOLIC = 27,
1763: MATOP_CHOLESKY_FACTOR_NUMERIC = 28,
1764: MATOP_SETUP = 29,
1765: MATOP_ILUFACTOR_SYMBOLIC = 30,
1766: MATOP_ICCFACTOR_SYMBOLIC = 31,
1767: MATOP_GET_DIAGONAL_BLOCK = 32,
1768: MATOP_SET_INF = 33,
1769: MATOP_DUPLICATE = 34,
1770: MATOP_FORWARD_SOLVE = 35,
1771: MATOP_BACKWARD_SOLVE = 36,
1772: MATOP_ILUFACTOR = 37,
1773: MATOP_ICCFACTOR = 38,
1774: MATOP_AXPY = 39,
1775: MATOP_CREATE_SUBMATRICES = 40,
1776: MATOP_INCREASE_OVERLAP = 41,
1777: MATOP_GET_VALUES = 42,
1778: MATOP_COPY = 43,
1779: MATOP_GET_ROW_MAX = 44,
1780: MATOP_SCALE = 45,
1781: MATOP_SHIFT = 46,
1782: MATOP_DIAGONAL_SET = 47,
1783: MATOP_ZERO_ROWS_COLUMNS = 48,
1784: MATOP_SET_RANDOM = 49,
1785: MATOP_GET_ROW_IJ = 50,
1786: MATOP_RESTORE_ROW_IJ = 51,
1787: MATOP_GET_COLUMN_IJ = 52,
1788: MATOP_RESTORE_COLUMN_IJ = 53,
1789: MATOP_FDCOLORING_CREATE = 54,
1790: MATOP_COLORING_PATCH = 55,
1791: MATOP_SET_UNFACTORED = 56,
1792: MATOP_PERMUTE = 57,
1793: MATOP_SET_VALUES_BLOCKED = 58,
1794: MATOP_CREATE_SUBMATRIX = 59,
1795: MATOP_DESTROY = 60,
1796: MATOP_VIEW = 61,
1797: MATOP_CONVERT_FROM = 62,
1798: /* MATOP_PLACEHOLDER_63=63 */
1799: MATOP_MATMAT_MULT_SYMBOLIC = 64,
1800: MATOP_MATMAT_MULT_NUMERIC = 65,
1801: MATOP_SET_LOCAL_TO_GLOBAL_MAP = 66,
1802: MATOP_SET_VALUES_LOCAL = 67,
1803: MATOP_ZERO_ROWS_LOCAL = 68,
1804: MATOP_GET_ROW_MAX_ABS = 69,
1805: MATOP_GET_ROW_MIN_ABS = 70,
1806: MATOP_CONVERT = 71,
1807: MATOP_HAS_OPERATION = 72,
1808: /* MATOP_PLACEHOLDER_73=73, */
1809: MATOP_SET_VALUES_ADIFOR = 74,
1810: MATOP_FD_COLORING_APPLY = 75,
1811: MATOP_SET_FROM_OPTIONS = 76,
1812: /* MATOP_PLACEHOLDER_77=77, */
1813: /* MATOP_PLACEHOLDER_78=78, */
1814: MATOP_FIND_ZERO_DIAGONALS = 79,
1815: MATOP_MULT_MULTIPLE = 80,
1816: MATOP_SOLVE_MULTIPLE = 81,
1817: MATOP_GET_INERTIA = 82,
1818: MATOP_LOAD = 83,
1819: MATOP_IS_SYMMETRIC = 84,
1820: MATOP_IS_HERMITIAN = 85,
1821: MATOP_IS_STRUCTURALLY_SYMMETRIC = 86,
1822: MATOP_SET_VALUES_BLOCKEDLOCAL = 87,
1823: MATOP_CREATE_VECS = 88,
1824: /* MATOP_PLACEHOLDER_89=89, */
1825: MATOP_MAT_MULT_SYMBOLIC = 90,
1826: MATOP_MAT_MULT_NUMERIC = 91,
1827: /* MATOP_PLACEHOLDER_92=92, */
1828: MATOP_PTAP_SYMBOLIC = 93,
1829: MATOP_PTAP_NUMERIC = 94,
1830: /* MATOP_PLACEHOLDER_95=95, */
1831: MATOP_MAT_TRANSPOSE_MULT_SYMBO = 96,
1832: MATOP_MAT_TRANSPOSE_MULT_NUMER = 97,
1833: MATOP_BIND_TO_CPU = 98,
1834: MATOP_PRODUCTSETFROMOPTIONS = 99,
1835: MATOP_PRODUCTSYMBOLIC = 100,
1836: MATOP_PRODUCTNUMERIC = 101,
1837: MATOP_CONJUGATE = 102,
1838: MATOP_VIEW_NATIVE = 103,
1839: MATOP_SET_VALUES_ROW = 104,
1840: MATOP_REAL_PART = 105,
1841: MATOP_IMAGINARY_PART = 106,
1842: MATOP_GET_ROW_UPPER_TRIANGULAR = 107,
1843: MATOP_RESTORE_ROW_UPPER_TRIANG = 108,
1844: MATOP_MAT_SOLVE = 109,
1845: MATOP_MAT_SOLVE_TRANSPOSE = 110,
1846: MATOP_GET_ROW_MIN = 111,
1847: MATOP_GET_COLUMN_VECTOR = 112,
1848: MATOP_MISSING_DIAGONAL = 113,
1849: MATOP_GET_SEQ_NONZERO_STRUCTUR = 114,
1850: MATOP_CREATE = 115,
1851: MATOP_GET_GHOSTS = 116,
1852: MATOP_GET_LOCAL_SUB_MATRIX = 117,
1853: MATOP_RESTORE_LOCALSUB_MATRIX = 118,
1854: MATOP_MULT_DIAGONAL_BLOCK = 119,
1855: MATOP_HERMITIAN_TRANSPOSE = 120,
1856: MATOP_MULT_HERMITIAN_TRANSPOSE = 121,
1857: MATOP_MULT_HERMITIAN_TRANS_ADD = 122,
1858: MATOP_GET_MULTI_PROC_BLOCK = 123,
1859: MATOP_FIND_NONZERO_ROWS = 124,
1860: MATOP_GET_COLUMN_NORMS = 125,
1861: MATOP_INVERT_BLOCK_DIAGONAL = 126,
1862: MATOP_INVERT_VBLOCK_DIAGONAL = 127,
1863: MATOP_CREATE_SUB_MATRICES_MPI = 128,
1864: MATOP_SET_VALUES_BATCH = 129,
1865: /* MATOP_PLACEHOLDER_130=130, */
1866: MATOP_TRANSPOSE_MAT_MULT_SYMBO = 131,
1867: MATOP_TRANSPOSE_MAT_MULT_NUMER = 132,
1868: MATOP_TRANSPOSE_COLORING_CREAT = 133,
1869: MATOP_TRANS_COLORING_APPLY_SPT = 134,
1870: MATOP_TRANS_COLORING_APPLY_DEN = 135,
1871: /* MATOP_PLACEHOLDER_136=136, */
1872: MATOP_RART_SYMBOLIC = 137,
1873: MATOP_RART_NUMERIC = 138,
1874: MATOP_SET_BLOCK_SIZES = 139,
1875: MATOP_AYPX = 140,
1876: MATOP_RESIDUAL = 141,
1877: MATOP_FDCOLORING_SETUP = 142,
1878: MATOP_FIND_OFFBLOCK_ENTRIES = 143,
1879: MATOP_MPICONCATENATESEQ = 144,
1880: MATOP_DESTROYSUBMATRICES = 145,
1881: MATOP_TRANSPOSE_SOLVE = 146,
1882: MATOP_GET_VALUES_LOCAL = 147
1883: } MatOperation;
1884: PETSC_EXTERN PetscErrorCode MatSetOperation(Mat, MatOperation, void (*)(void));
1885: PETSC_EXTERN PetscErrorCode MatGetOperation(Mat, MatOperation, void (**)(void));
1886: PETSC_EXTERN PetscErrorCode MatHasOperation(Mat, MatOperation, PetscBool *);
1887: PETSC_EXTERN PetscErrorCode MatHasCongruentLayouts(Mat, PetscBool *);
1888: PETSC_DEPRECATED_FUNCTION("Use MatProductClear() (since version 3.14)") static inline PetscErrorCode MatFreeIntermediateDataStructures(Mat A)
1889: {
1890: return MatProductClear(A);
1891: }
1892: PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat, MatOperation, void (*)(void));
1893: PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat, MatOperation, void (**)(void));
1894: PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat, void *);
1895: PETSC_EXTERN PetscErrorCode MatShellSetContextDestroy(Mat, PetscErrorCode (*)(void *));
1896: PETSC_EXTERN PetscErrorCode MatShellSetVecType(Mat, VecType);
1897: PETSC_EXTERN PetscErrorCode MatShellTestMult(Mat, PetscErrorCode (*)(void *, Vec, Vec), Vec, void *, PetscBool *);
1898: PETSC_EXTERN PetscErrorCode MatShellTestMultTranspose(Mat, PetscErrorCode (*)(void *, Vec, Vec), Vec, void *, PetscBool *);
1899: PETSC_EXTERN PetscErrorCode MatShellSetManageScalingShifts(Mat);
1900: PETSC_EXTERN PetscErrorCode MatShellSetMatProductOperation(Mat, MatProductType, PetscErrorCode (*)(Mat, Mat, Mat, void **), PetscErrorCode (*)(Mat, Mat, Mat, void *), PetscErrorCode (*)(void *), MatType, MatType);
1901: PETSC_EXTERN PetscErrorCode MatIsShell(Mat, PetscBool *);
1903: /*
1904: Codes for matrices stored on disk. By default they are
1905: stored in a universal format. By changing the format with
1906: PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
1907: be stored in a way natural for the matrix, for example dense matrices
1908: would be stored as dense. Matrices stored this way may only be
1909: read into matrices of the same type.
1910: */
1911: #define MATRIX_BINARY_FORMAT_DENSE -1
1913: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat, PetscReal);
1915: PETSC_EXTERN PetscErrorCode MatISSetLocalMatType(Mat, MatType);
1916: PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[]);
1917: PETSC_EXTERN PetscErrorCode MatISStoreL2L(Mat, PetscBool);
1918: PETSC_EXTERN PetscErrorCode MatISFixLocalEmpty(Mat, PetscBool);
1919: PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat, Mat *);
1920: PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat, Mat *);
1921: PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat, Mat);
1922: PETSC_EXTERN PetscErrorCode MatISGetLocalToGlobalMapping(Mat, ISLocalToGlobalMapping *, ISLocalToGlobalMapping *);
1923: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use the MatConvert() interface (since version 3.10)") PetscErrorCode MatISGetMPIXAIJ(Mat, MatReuse, Mat *);
1925: /*S
1926: MatNullSpace - Object that removes a null space from a vector, i.e.
1927: orthogonalizes the vector to a subspace
1929: Level: advanced
1931: .seealso: `MatNullSpaceCreate()`, `MatNullSpaceSetFunction()`, `MatGetNullSpace()`, `MatSetNullSpace()`
1932: S*/
1933: typedef struct _p_MatNullSpace *MatNullSpace;
1935: PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm, PetscBool, PetscInt, const Vec[], MatNullSpace *);
1936: PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace, PetscErrorCode (*)(MatNullSpace, Vec, void *), void *);
1937: PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace *);
1938: PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace, Vec);
1939: PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1940: PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *);
1941: PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat, MatNullSpace);
1942: PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat, MatNullSpace);
1943: PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat, MatNullSpace);
1944: PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat, MatNullSpace *);
1945: PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace, Mat, PetscBool *);
1946: PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace, PetscViewer);
1947: PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace, PetscBool *, PetscInt *, const Vec **);
1948: PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec, MatNullSpace *);
1950: PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat, IS);
1951: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat, PetscReal);
1952: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat, PetscInt *);
1954: PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat, PetscInt, Mat *);
1955: PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat, PetscInt, Mat *);
1956: PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat, Mat *);
1958: PETSC_EXTERN PetscErrorCode MatComputeOperator(Mat, MatType, Mat *);
1959: PETSC_EXTERN PetscErrorCode MatComputeOperatorTranspose(Mat, MatType, Mat *);
1961: PETSC_DEPRECATED_FUNCTION("Use MatComputeOperator() (since version 3.12)") static inline PetscErrorCode MatComputeExplicitOperator(Mat A, Mat *B)
1962: {
1963: return MatComputeOperator(A, NULL, B);
1964: }
1965: PETSC_DEPRECATED_FUNCTION("Use MatComputeOperatorTranspose() (since version 3.12)") static inline PetscErrorCode MatComputeExplicitOperatorTranspose(Mat A, Mat *B)
1966: {
1967: return MatComputeOperatorTranspose(A, NULL, B);
1968: }
1970: PETSC_EXTERN PetscErrorCode MatCreateKAIJ(Mat, PetscInt, PetscInt, const PetscScalar[], const PetscScalar[], Mat *);
1971: PETSC_EXTERN PetscErrorCode MatKAIJGetAIJ(Mat, Mat *);
1972: PETSC_EXTERN PetscErrorCode MatKAIJGetS(Mat, PetscInt *, PetscInt *, PetscScalar **);
1973: PETSC_EXTERN PetscErrorCode MatKAIJGetSRead(Mat, PetscInt *, PetscInt *, const PetscScalar **);
1974: PETSC_EXTERN PetscErrorCode MatKAIJRestoreS(Mat, PetscScalar **);
1975: PETSC_EXTERN PetscErrorCode MatKAIJRestoreSRead(Mat, const PetscScalar **);
1976: PETSC_EXTERN PetscErrorCode MatKAIJGetT(Mat, PetscInt *, PetscInt *, PetscScalar **);
1977: PETSC_EXTERN PetscErrorCode MatKAIJGetTRead(Mat, PetscInt *, PetscInt *, const PetscScalar **);
1978: PETSC_EXTERN PetscErrorCode MatKAIJRestoreT(Mat, PetscScalar **);
1979: PETSC_EXTERN PetscErrorCode MatKAIJRestoreTRead(Mat, const PetscScalar **);
1980: PETSC_EXTERN PetscErrorCode MatKAIJSetAIJ(Mat, Mat);
1981: PETSC_EXTERN PetscErrorCode MatKAIJSetS(Mat, PetscInt, PetscInt, const PetscScalar[]);
1982: PETSC_EXTERN PetscErrorCode MatKAIJSetT(Mat, PetscInt, PetscInt, const PetscScalar[]);
1983: PETSC_EXTERN PetscErrorCode MatKAIJGetScaledIdentity(Mat, PetscBool *);
1985: PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat, Vec);
1987: PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, Mat *);
1988: PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat, Vec, Vec);
1989: PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat, PetscErrorCode (*)(void *, Vec, Vec), void *);
1990: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat, PetscErrorCode (*)(void *, PetscInt, Vec, PetscScalar *));
1991: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat, PetscErrorCode (*)(void *, Vec));
1992: PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat, PetscScalar[], PetscInt);
1993: PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1994: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat, PetscReal);
1995: PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat, PetscInt);
1996: PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat, PetscScalar *);
1997: PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat, const char[]);
1998: PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void *, Vec, Vec, PetscScalar *);
1999: PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat, PetscErrorCode (*)(void *, Vec, Vec, PetscScalar *), void *);
2001: /*S
2002: MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
2003: Jacobian vector products
2005: Notes:
2006: `MATMFFD` is a specific `MatType` which uses the `MatMFFD` data structure
2008: MatMFFD*() methods actually take the Mat as their first argument. Not a `MatMFFD` data structure
2010: Level: developer
2012: .seealso: `MATMFFD`, `MatCreateMFFD()`, `MatMFFDSetFuction()`, `MatMFFDSetType()`, `MatMFFDRegister()`
2013: S*/
2014: typedef struct _p_MatMFFD *MatMFFD;
2016: /*J
2017: MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
2019: Level: beginner
2021: .seealso: `MatMFFDSetType()`, `MatMFFDRegister()`
2022: J*/
2023: typedef const char *MatMFFDType;
2024: #define MATMFFD_DS "ds"
2025: #define MATMFFD_WP "wp"
2027: PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat, MatMFFDType);
2028: PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[], PetscErrorCode (*)(MatMFFD));
2030: PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat, PetscReal);
2031: PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat, PetscBool);
2033: PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring, MatMFFDType);
2035: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
2036: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
2038: #ifdef PETSC_HAVE_H2OPUS
2039: PETSC_EXTERN_TYPEDEF typedef PetscScalar (*MatH2OpusKernel)(PetscInt, PetscReal[], PetscReal[], void *);
2040: PETSC_EXTERN PetscErrorCode MatCreateH2OpusFromKernel(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscReal[], PetscBool, MatH2OpusKernel, void *, PetscReal, PetscInt, PetscInt, Mat *);
2041: PETSC_EXTERN PetscErrorCode MatCreateH2OpusFromMat(Mat, PetscInt, const PetscReal[], PetscBool, PetscReal, PetscInt, PetscInt, PetscInt, PetscReal, Mat *);
2042: PETSC_EXTERN PetscErrorCode MatH2OpusSetSamplingMat(Mat, Mat, PetscInt, PetscReal);
2043: PETSC_EXTERN PetscErrorCode MatH2OpusOrthogonalize(Mat);
2044: PETSC_EXTERN PetscErrorCode MatH2OpusCompress(Mat, PetscReal);
2045: PETSC_EXTERN PetscErrorCode MatH2OpusSetNativeMult(Mat, PetscBool);
2046: PETSC_EXTERN PetscErrorCode MatH2OpusGetNativeMult(Mat, PetscBool *);
2047: PETSC_EXTERN PetscErrorCode MatH2OpusGetIndexMap(Mat, IS *);
2048: PETSC_EXTERN PetscErrorCode MatH2OpusMapVec(Mat, PetscBool, Vec, Vec *);
2049: PETSC_EXTERN PetscErrorCode MatH2OpusLowRankUpdate(Mat, Mat, Mat, PetscScalar);
2050: #endif
2052: #ifdef PETSC_HAVE_HTOOL
2053: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*MatHtoolKernel)(PetscInt, PetscInt, PetscInt, const PetscInt *, const PetscInt *, PetscScalar *, void *);
2054: PETSC_EXTERN PetscErrorCode MatCreateHtoolFromKernel(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscReal[], const PetscReal[], MatHtoolKernel, void *, Mat *);
2055: PETSC_EXTERN PetscErrorCode MatHtoolSetKernel(Mat, MatHtoolKernel, void *);
2056: PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationSource(Mat, IS *);
2057: PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationTarget(Mat, IS *);
2058: PETSC_EXTERN PetscErrorCode MatHtoolUsePermutation(Mat, PetscBool);
2059: /*E
2060: MatHtoolCompressorType - Indicates the type of compressor used by a `MATHTOOL`
2062: Level: beginner
2064: Values:
2065: + `MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA` (default) - symmetric partial adaptive cross approximation
2066: . `MAT_HTOOL_COMPRESSOR_FULL_ACA` - full adaptive cross approximation
2067: - `MAT_HTOOL_COMPRESSOR_SVD` - singular value decomposition
2069: .seealso: `MatCreateHtoolFromKernel()`, `MATHTOOL`, `MatHtoolClusteringType`
2070: E*/
2071: typedef enum {
2072: MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA,
2073: MAT_HTOOL_COMPRESSOR_FULL_ACA,
2074: MAT_HTOOL_COMPRESSOR_SVD
2075: } MatHtoolCompressorType;
2076: /*E
2077: MatHtoolClusteringType - Indicates the type of clustering used by a `MATHTOOL`
2079: Level: beginner
2081: Values:
2082: + `MAT_HTOOL_CLUSTERING_PCA_REGULAR` (default) - axis computed via principle component analysis, split uniformly
2083: . `MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC` - axis computed via principle component analysis, split barycentrically
2084: . `MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR` - axis along the largest extent of the bounding box, split uniformly
2085: - `MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC` - axis along the largest extent of the bounding box, split barycentrically
2087: Notes: higher-dimensional clustering is not yet supported in Htool, but once it is, one should add BOUNDING_BOX_{2,3} types
2089: .seealso: `MatCreateHtoolFromKernel()`, `MATHTOOL`, `MatHtoolCompressorType`
2090: E*/
2091: typedef enum {
2092: MAT_HTOOL_CLUSTERING_PCA_REGULAR,
2093: MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC,
2094: MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR,
2095: MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC
2096: } MatHtoolClusteringType;
2097: #endif
2099: #ifdef PETSC_HAVE_MUMPS
2100: PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat, PetscInt, PetscInt);
2101: PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat, PetscInt, PetscInt *);
2102: PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat, PetscInt, PetscReal);
2103: PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat, PetscInt, PetscReal *);
2105: PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat, PetscInt, PetscInt *);
2106: PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat, PetscInt, PetscInt *);
2107: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat, PetscInt, PetscReal *);
2108: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat, PetscInt, PetscReal *);
2109: PETSC_EXTERN PetscErrorCode MatMumpsGetInverse(Mat, Mat);
2110: PETSC_EXTERN PetscErrorCode MatMumpsGetInverseTranspose(Mat, Mat);
2111: #endif
2113: #ifdef PETSC_HAVE_MKL_PARDISO
2114: PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat, PetscInt, PetscInt);
2115: #endif
2117: #ifdef PETSC_HAVE_MKL_CPARDISO
2118: PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat, PetscInt, PetscInt);
2119: #endif
2121: #ifdef PETSC_HAVE_SUPERLU
2122: PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat, PetscReal);
2123: #endif
2125: #ifdef PETSC_HAVE_SUPERLU_DIST
2126: PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat, PetscScalar *);
2127: #endif
2129: #ifdef PETSC_HAVE_STRUMPACK
2130: /*E
2131: MatSTRUMPACKReordering - sparsity reducing ordering to be used in `STRUMPACK`
2133: Level: intermediate
2134: E*/
2135: typedef enum {
2136: MAT_STRUMPACK_NATURAL = 0,
2137: MAT_STRUMPACK_METIS = 1,
2138: MAT_STRUMPACK_PARMETIS = 2,
2139: MAT_STRUMPACK_SCOTCH = 3,
2140: MAT_STRUMPACK_PTSCOTCH = 4,
2141: MAT_STRUMPACK_RCM = 5
2142: } MatSTRUMPACKReordering;
2143: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetReordering(Mat, MatSTRUMPACKReordering);
2144: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat, PetscBool);
2145: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat, PetscReal);
2146: PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSRelTol() (since version 3.9)") static inline PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat mat, PetscReal rtol)
2147: {
2148: return MatSTRUMPACKSetHSSRelTol(mat, rtol);
2149: }
2150: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat, PetscReal);
2151: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat, PetscInt);
2152: PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSMinSepSize() (since version 3.9)") static inline PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat mat, PetscInt hssminsize)
2153: {
2154: return MatSTRUMPACKSetHSSMinSepSize(mat, hssminsize);
2155: }
2156: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat, PetscInt);
2157: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat, PetscInt);
2158: #endif
2160: PETSC_EXTERN PetscErrorCode MatBindToCPU(Mat, PetscBool);
2161: PETSC_EXTERN PetscErrorCode MatBoundToCPU(Mat, PetscBool *);
2162: PETSC_DEPRECATED_FUNCTION("Use MatBindToCPU (since v3.13)") static inline PetscErrorCode MatPinToCPU(Mat A, PetscBool flg)
2163: {
2164: return MatBindToCPU(A, flg);
2165: }
2166: PETSC_EXTERN PetscErrorCode MatSetBindingPropagates(Mat, PetscBool);
2167: PETSC_EXTERN PetscErrorCode MatGetBindingPropagates(Mat, PetscBool *);
2169: typedef struct _n_SplitCSRMat *PetscSplitCSRDataStructure;
2170: PETSC_EXTERN PetscErrorCode MatCUSPARSEGetDeviceMatWrite(Mat, PetscSplitCSRDataStructure *);
2172: #ifdef PETSC_HAVE_KOKKOS_KERNELS
2173: PETSC_EXTERN PetscErrorCode MatKokkosGetDeviceMatWrite(Mat, PetscSplitCSRDataStructure *);
2174: PETSC_EXTERN PetscErrorCode MatSeqAIJKokkosSetDeviceMat(Mat, PetscSplitCSRDataStructure);
2175: PETSC_EXTERN PetscErrorCode MatSeqAIJKokkosGetDeviceMat(Mat, PetscSplitCSRDataStructure *);
2176: #endif
2178: #ifdef PETSC_HAVE_CUDA
2179: /*E
2180: MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
2181: matrices.
2183: Not Collective
2185: + `MAT_CUSPARSE_CSR` - Compressed Sparse Row
2186: . `MAT_CUSPARSE_ELL` - Ellpack (requires CUDA 4.2 or later).
2187: - `MAT_CUSPARSE_HYB` - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).
2189: Level: intermediate
2191: Developer Note:
2192: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
2194: .seealso: `MatCUSPARSESetFormat()`, `MatCUSPARSEFormatOperation`
2195: E*/
2197: typedef enum {
2198: MAT_CUSPARSE_CSR,
2199: MAT_CUSPARSE_ELL,
2200: MAT_CUSPARSE_HYB
2201: } MatCUSPARSEStorageFormat;
2203: /* these will be strings associated with enumerated type defined above */
2204: PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
2206: /*E
2207: MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
2208: matrices whose operation should use a particular storage format.
2210: Not Collective
2212: + `MAT_CUSPARSE_MULT_DIAG` - sets the storage format for the diagonal matrix in the parallel MatMult
2213: . `MAT_CUSPARSE_MULT_OFFDIAG` - sets the storage format for the offdiagonal matrix in the parallel MatMult
2214: . `MAT_CUSPARSE_MULT` - sets the storage format for the entire matrix in the serial (single GPU) MatMult
2215: - `MAT_CUSPARSE_ALL` - sets the storage format for all CUSPARSE (GPU) matrices
2217: Level: intermediate
2219: .seealso: `MatCUSPARSESetFormat()`, `MatCUSPARSEStorageFormat`
2220: E*/
2221: typedef enum {
2222: MAT_CUSPARSE_MULT_DIAG,
2223: MAT_CUSPARSE_MULT_OFFDIAG,
2224: MAT_CUSPARSE_MULT,
2225: MAT_CUSPARSE_ALL
2226: } MatCUSPARSEFormatOperation;
2228: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *);
2229: PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *);
2230: PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat, MatCUSPARSEFormatOperation, MatCUSPARSEStorageFormat);
2231: PETSC_EXTERN PetscErrorCode MatCUSPARSESetUseCPUSolve(Mat, PetscBool);
2232: typedef struct Mat_SeqAIJCUSPARSETriFactors *Mat_SeqAIJCUSPARSETriFactors_p;
2233: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetIJ(Mat, PetscBool, const int **, const int **);
2234: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreIJ(Mat, PetscBool, const int **, const int **);
2235: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetArrayRead(Mat, const PetscScalar **);
2236: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArrayRead(Mat, const PetscScalar **);
2237: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetArrayWrite(Mat, PetscScalar **);
2238: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArrayWrite(Mat, PetscScalar **);
2239: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetArray(Mat, PetscScalar **);
2240: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArray(Mat, PetscScalar **);
2242: PETSC_EXTERN PetscErrorCode MatCreateDenseCUDA(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar[], Mat *);
2243: PETSC_EXTERN PetscErrorCode MatCreateSeqDenseCUDA(MPI_Comm, PetscInt, PetscInt, PetscScalar[], Mat *);
2244: PETSC_EXTERN PetscErrorCode MatMPIDenseCUDASetPreallocation(Mat, PetscScalar[]);
2245: PETSC_EXTERN PetscErrorCode MatSeqDenseCUDASetPreallocation(Mat, PetscScalar[]);
2246: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayWrite(Mat, PetscScalar **);
2247: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayRead(Mat, const PetscScalar **);
2248: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArray(Mat, PetscScalar **);
2249: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayWrite(Mat, PetscScalar **);
2250: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayRead(Mat, const PetscScalar **);
2251: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArray(Mat, PetscScalar **);
2252: PETSC_EXTERN PetscErrorCode MatDenseCUDAPlaceArray(Mat, const PetscScalar *);
2253: PETSC_EXTERN PetscErrorCode MatDenseCUDAReplaceArray(Mat, const PetscScalar *);
2254: PETSC_EXTERN PetscErrorCode MatDenseCUDAResetArray(Mat);
2256: #endif
2258: #if defined(PETSC_HAVE_VIENNACL)
2259: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *);
2260: PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *);
2261: #endif
2263: #if defined(PETSC_HAVE_FFTW)
2264: PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat, Vec, Vec);
2265: PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat, Vec, Vec);
2266: PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat, Vec *, Vec *, Vec *);
2267: #endif
2269: #if defined(PETSC_HAVE_SCALAPACK)
2270: PETSC_EXTERN PetscErrorCode MatCreateScaLAPACK(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, Mat *);
2271: PETSC_EXTERN PetscErrorCode MatScaLAPACKSetBlockSizes(Mat, PetscInt, PetscInt);
2272: PETSC_EXTERN PetscErrorCode MatScaLAPACKGetBlockSizes(Mat, PetscInt *, PetscInt *);
2273: #endif
2275: PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm, PetscInt, const IS[], PetscInt, const IS[], const Mat[], Mat *);
2276: PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat, PetscInt *, PetscInt *);
2277: PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat, IS[], IS[]);
2278: PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat, IS[], IS[]);
2279: PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat, PetscInt *, PetscInt *, Mat ***);
2280: PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat, PetscInt, PetscInt, Mat *);
2281: PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat, VecType);
2282: PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat, PetscInt, const IS[], PetscInt, const IS[], const Mat[]);
2283: PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat, PetscInt, PetscInt, Mat);
2285: PETSC_EXTERN PetscErrorCode MatChop(Mat, PetscReal);
2286: PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat, PetscReal, PetscInt *);
2288: PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat, PetscInt, PetscInt *, IS **);
2290: PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat, PetscBool, Mat);
2292: PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat, Mat *);
2293: PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat, Mat *);
2295: PETSC_EXTERN PetscErrorCode MatSeqAIJGetCSRAndMemType(Mat, const PetscInt **, const PetscInt **, PetscScalar **, PetscMemType *);
2297: PETSC_EXTERN PetscErrorCode MatCreateGraph(Mat, PetscBool, PetscBool, PetscReal, Mat *);
2298: #endif