Actual source code: petscmat.h
petsc-3.14.1 2020-11-03
1: /*
2: Include file for the matrix component of PETSc
3: */
4: #ifndef PETSCMAT_H
5: #define PETSCMAT_H
6: #include <petscvec.h>
8: /*S
9: Mat - Abstract PETSc matrix object used to manage all linear operators in PETSc, even those without
10: an explicit sparse representation (such as matrix-free operators)
12: Level: beginner
14: .seealso: MatCreate(), MatType, MatSetType(), MatDestroy()
15: S*/
16: typedef struct _p_Mat* Mat;
18: /*J
19: MatType - String with the name of a PETSc matrix type
21: Level: beginner
23: .seealso: MatSetType(), Mat, MatSolverType, MatRegister()
24: J*/
25: typedef const char* MatType;
26: #define MATSAME "same"
27: #define MATMAIJ "maij"
28: #define MATSEQMAIJ "seqmaij"
29: #define MATMPIMAIJ "mpimaij"
30: #define MATKAIJ "kaij"
31: #define MATSEQKAIJ "seqkaij"
32: #define MATMPIKAIJ "mpikaij"
33: #define MATIS "is"
34: #define MATAIJ "aij"
35: #define MATSEQAIJ "seqaij"
36: #define MATMPIAIJ "mpiaij"
37: #define MATAIJCRL "aijcrl"
38: #define MATSEQAIJCRL "seqaijcrl"
39: #define MATMPIAIJCRL "mpiaijcrl"
40: #define MATAIJCUSPARSE "aijcusparse"
41: #define MATSEQAIJCUSPARSE "seqaijcusparse"
42: #define MATMPIAIJCUSPARSE "mpiaijcusparse"
43: #define MATAIJKOKKOS "aijkokkos"
44: #define MATSEQAIJKOKKOS "seqaijkokkos"
45: #define MATMPIAIJKOKKOS "mpiaijkokkos"
46: #define MATAIJVIENNACL "aijviennacl"
47: #define MATSEQAIJVIENNACL "seqaijviennacl"
48: #define MATMPIAIJVIENNACL "mpiaijviennacl"
49: #define MATAIJPERM "aijperm"
50: #define MATSEQAIJPERM "seqaijperm"
51: #define MATMPIAIJPERM "mpiaijperm"
52: #define MATAIJSELL "aijsell"
53: #define MATSEQAIJSELL "seqaijsell"
54: #define MATMPIAIJSELL "mpiaijsell"
55: #define MATAIJMKL "aijmkl"
56: #define MATSEQAIJMKL "seqaijmkl"
57: #define MATMPIAIJMKL "mpiaijmkl"
58: #define MATBAIJMKL "baijmkl"
59: #define MATSEQBAIJMKL "seqbaijmkl"
60: #define MATMPIBAIJMKL "mpibaijmkl"
61: #define MATSHELL "shell"
62: #define MATDENSE "dense"
63: #define MATDENSECUDA "densecuda"
64: #define MATSEQDENSE "seqdense"
65: #define MATSEQDENSECUDA "seqdensecuda"
66: #define MATMPIDENSE "mpidense"
67: #define MATMPIDENSECUDA "mpidensecuda"
68: #define MATELEMENTAL "elemental"
69: #define MATSCALAPACK "scalapack"
70: #define MATBAIJ "baij"
71: #define MATSEQBAIJ "seqbaij"
72: #define MATMPIBAIJ "mpibaij"
73: #define MATMPIADJ "mpiadj"
74: #define MATSBAIJ "sbaij"
75: #define MATSEQSBAIJ "seqsbaij"
76: #define MATMPISBAIJ "mpisbaij"
77: #define MATDAAD "daad"
78: #define MATMFFD "mffd"
79: #define MATNORMAL "normal"
80: #define MATNORMALHERMITIAN "normalh"
81: #define MATLRC "lrc"
82: #define MATSCATTER "scatter"
83: #define MATBLOCKMAT "blockmat"
84: #define MATCOMPOSITE "composite"
85: #define MATFFT "fft"
86: #define MATFFTW "fftw"
87: #define MATSEQCUFFT "seqcufft"
88: #define MATTRANSPOSEMAT "transpose"
89: #define MATSCHURCOMPLEMENT "schurcomplement"
90: #define MATPYTHON "python"
91: #define MATHYPRE "hypre"
92: #define MATHYPRESTRUCT "hyprestruct"
93: #define MATHYPRESSTRUCT "hypresstruct"
94: #define MATSUBMATRIX "submatrix"
95: #define MATLOCALREF "localref"
96: #define MATNEST "nest"
97: #define MATPREALLOCATOR "preallocator"
98: #define MATSELL "sell"
99: #define MATSEQSELL "seqsell"
100: #define MATMPISELL "mpisell"
101: #define MATDUMMY "dummy"
102: #define MATLMVM "lmvm"
103: #define MATLMVMDFP "lmvmdfp"
104: #define MATLMVMBFGS "lmvmbfgs"
105: #define MATLMVMSR1 "lmvmsr1"
106: #define MATLMVMBROYDEN "lmvmbroyden"
107: #define MATLMVMBADBROYDEN "lmvmbadbroyden"
108: #define MATLMVMSYMBROYDEN "lmvmsymbroyden"
109: #define MATLMVMSYMBADBROYDEN "lmvmsymbadbroyden"
110: #define MATLMVMDIAGBROYDEN "lmvmdiagbroyden"
111: #define MATCONSTANTDIAGONAL "constantdiagonal"
112: #define MATHARA "hara"
114: /*J
115: MatSolverType - String with the name of a PETSc matrix solver type.
117: For example: "petsc" indicates what PETSc provides, "superlu_dist" the parallel SuperLU_DIST package etc
119: Level: beginner
121: Notes: MATSOLVERUMFPACK, MATSOLVERCHOLMOD, MATSOLVERKLU form the SuiteSparse package for which you can use --download-suitesparse
123: .seealso: MatGetFactor(), PCFactorSetMatSolverType(), PCFactorGetMatSolverType()
124: J*/
125: typedef const char* MatSolverType;
126: #define MATSOLVERSUPERLU "superlu"
127: #define MATSOLVERSUPERLU_DIST "superlu_dist"
128: #define MATSOLVERSTRUMPACK "strumpack"
129: #define MATSOLVERUMFPACK "umfpack"
130: #define MATSOLVERCHOLMOD "cholmod"
131: #define MATSOLVERKLU "klu"
132: #define MATSOLVERSPARSEELEMENTAL "sparseelemental"
133: #define MATSOLVERELEMENTAL "elemental"
134: #define MATSOLVERSCALAPACK "scalapack"
135: #define MATSOLVERESSL "essl"
136: #define MATSOLVERLUSOL "lusol"
137: #define MATSOLVERMUMPS "mumps"
138: #define MATSOLVERMKL_PARDISO "mkl_pardiso"
139: #define MATSOLVERMKL_CPARDISO "mkl_cpardiso"
140: #define MATSOLVERPASTIX "pastix"
141: #define MATSOLVERMATLAB "matlab"
142: #define MATSOLVERPETSC "petsc"
143: #define MATSOLVERBAS "bas"
144: #define MATSOLVERCUSPARSE "cusparse"
145: #define MATSOLVERCUDA "cuda"
147: /*E
148: MatFactorType - indicates what type of factorization is requested
150: Level: beginner
152: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
154: .seealso: MatSolverType, MatGetFactor(), MatGetFactorAvailable(), MatSolverTypeRegister()
155: E*/
156: typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
157: PETSC_EXTERN const char *const MatFactorTypes[];
159: PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,MatSolverType,MatFactorType,Mat*);
160: PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,MatSolverType,MatFactorType,PetscBool *);
161: PETSC_EXTERN PetscErrorCode MatFactorGetUseOrdering(Mat, PetscBool*);
162: PETSC_EXTERN PetscErrorCode MatFactorGetSolverType(Mat,MatSolverType*);
163: PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*);
164: PETSC_EXTERN PetscErrorCode MatSetFactorType(Mat,MatFactorType);
165: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*MatSolverFunction)(Mat,MatFactorType,Mat*);
166: PETSC_EXTERN PetscErrorCode MatSolverTypeRegister(MatSolverType,MatType,MatFactorType,MatSolverFunction);
167: PETSC_EXTERN PetscErrorCode MatSolverTypeGet(MatSolverType,MatType,MatFactorType,PetscBool*,PetscBool*,MatSolverFunction*);
168: typedef MatSolverType MatSolverPackage PETSC_DEPRECATED_TYPEDEF("Use MatSolverType (since version 3.9)");
169: PETSC_DEPRECATED_FUNCTION("Use MatSolverTypeRegister() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSolverPackageRegister(MatSolverType stype,MatType mtype,MatFactorType ftype,MatSolverFunction f)
170: { return MatSolverTypeRegister(stype,mtype,ftype,f); }
171: PETSC_DEPRECATED_FUNCTION("Use MatSolverTypeGet() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSolverPackageGet(MatSolverType stype,MatType mtype,MatFactorType ftype,PetscBool *foundmtype,PetscBool *foundstype,MatSolverFunction *f)
172: { return MatSolverTypeGet(stype,mtype,ftype,foundmtype,foundstype,f); }
174: /*E
175: MatProductType - indicates what type of matrix product is requested
177: Level: beginner
179: .seealso: MatSolverType, MatProductSetType()
180: E*/
181: typedef enum {MATPRODUCT_UNSPECIFIED=0,MATPRODUCT_AB,MATPRODUCT_AtB,MATPRODUCT_ABt,MATPRODUCT_PtAP,MATPRODUCT_RARt,MATPRODUCT_ABC} MatProductType;
182: PETSC_EXTERN const char *const MatProductTypes[];
184: /*J
185: MatProductAlgorithm - String with the name of an algorithm for a PETSc matrix product implementation
187: Level: beginner
189: .seealso: MatSetType(), Mat, MatSolverType, MatRegister(), MatProductSetAlgorithm(), MatProductType
190: J*/
191: typedef const char* MatProductAlgorithm;
192: #define MATPRODUCTALGORITHM_DEFAULT "default"
194: PETSC_EXTERN PetscErrorCode MatProductCreate(Mat,Mat,Mat,Mat*);
195: PETSC_EXTERN PetscErrorCode MatProductCreateWithMat(Mat,Mat,Mat,Mat);
196: PETSC_EXTERN PetscErrorCode MatProductSetType(Mat,MatProductType);
197: PETSC_EXTERN PetscErrorCode MatProductSetAlgorithm(Mat,MatProductAlgorithm);
198: PETSC_EXTERN PetscErrorCode MatProductSetFill(Mat,PetscReal);
199: PETSC_EXTERN PetscErrorCode MatProductSetFromOptions(Mat);
200: PETSC_EXTERN PetscErrorCode MatProductSymbolic(Mat);
201: PETSC_EXTERN PetscErrorCode MatProductNumeric(Mat);
202: PETSC_EXTERN PetscErrorCode MatProductReplaceMats(Mat,Mat,Mat,Mat);
203: PETSC_EXTERN PetscErrorCode MatProductClear(Mat);
204: PETSC_EXTERN PetscErrorCode MatProductView(Mat,PetscViewer);
206: /* Logging support */
207: #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */
208: PETSC_EXTERN PetscClassId MAT_CLASSID;
209: PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID;
210: PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID;
211: PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
212: PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID;
213: PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID;
214: PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID;
215: PETSC_EXTERN PetscClassId MATMFFD_CLASSID;
217: /*E
218: MatReuse - Indicates if matrices obtained from a previous call to MatCreateSubMatrices(), MatCreateSubMatrix(), MatConvert() or several other functions
219: are to be reused to store the new matrix values.
221: $ MAT_INITIAL_MATRIX - create a new matrix
222: $ MAT_REUSE_MATRIX - reuse the matrix created with a previous call that used MAT_INITIAL_MATRIX
223: $ MAT_INPLACE_MATRIX - replace the first input matrix with the new matrix (not applicable to all functions)
224: $ MAT_IGNORE_MATRIX - do not create a new matrix or reuse a give matrix, just ignore that matrix argument (not applicable to all functions)
226: Level: beginner
228: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
230: .seealso: MatCreateSubMatrices(), MatCreateSubMatrix(), MatDestroyMatrices(), MatConvert()
231: E*/
232: typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse;
234: /*E
235: MatCreateSubMatrixOption - Indicates if matrices obtained from a call to MatCreateSubMatrices()
236: include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
238: Level: beginner
240: .seealso: MatGetSeqNonzerostructure()
241: E*/
242: typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatCreateSubMatrixOption;
244: PETSC_EXTERN PetscErrorCode MatInitializePackage(void);
246: PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*);
247: PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
248: PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType);
249: PETSC_EXTERN PetscErrorCode MatGetVecType(Mat,VecType*);
250: PETSC_EXTERN PetscErrorCode MatSetVecType(Mat,VecType);
251: PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat);
252: PETSC_EXTERN PetscErrorCode MatViewFromOptions(Mat,PetscObject,const char[]);
253: PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat));
254: PETSC_EXTERN PetscErrorCode MatRegisterRootName(const char[],const char[],const char[]);
255: PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]);
256: PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]);
257: PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]);
258: PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat,PetscBool);
260: PETSC_EXTERN PetscFunctionList MatList;
261: PETSC_EXTERN PetscFunctionList MatColoringList;
262: PETSC_EXTERN PetscFunctionList MatPartitioningList;
264: /*E
265: MatStructure - Indicates if two matrices have the same nonzero structure
267: Level: beginner
269: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
271: .seealso: MatCopy(), MatAXPY()
272: E*/
273: typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure;
275: #if defined PETSC_HAVE_MKL_SPARSE
276: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
277: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
278: PETSC_EXTERN PetscErrorCode MatCreateBAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
279: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
280: #endif
282: PETSC_EXTERN PetscErrorCode MatCreateSeqSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
283: PETSC_EXTERN PetscErrorCode MatCreateSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
284: PETSC_EXTERN PetscErrorCode MatSeqSELLSetPreallocation(Mat,PetscInt,const PetscInt[]);
285: PETSC_EXTERN PetscErrorCode MatMPISELLSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
287: PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
288: PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
289: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
290: PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
291: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
292: PETSC_EXTERN PetscErrorCode MatUpdateMPIAIJWithArrays(Mat,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
293: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
294: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSeqAIJ(MPI_Comm,Mat,Mat,const PetscInt[],Mat*);
296: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
297: PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
298: PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
300: PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
301: PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
303: PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
304: PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
305: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
306: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
307: PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]);
309: PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
310: PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
311: PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*);
312: PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Vec,Mat,Mat*);
313: PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat,Mat*,Mat*,Vec*,Mat*);
314: PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*);
315: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
316: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
318: PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
319: PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
320: PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
321: PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
322: PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
323: PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
324: typedef enum {MAT_COMPOSITE_MERGE_RIGHT,MAT_COMPOSITE_MERGE_LEFT} MatCompositeMergeType;
325: PETSC_EXTERN PetscErrorCode MatCompositeSetMergeType(Mat,MatCompositeMergeType);
326: PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
327: typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
328: PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
329: PETSC_EXTERN PetscErrorCode MatCompositeGetType(Mat,MatCompositeType*);
330: PETSC_EXTERN PetscErrorCode MatCompositeSetMatStructure(Mat,MatStructure);
331: PETSC_EXTERN PetscErrorCode MatCompositeGetMatStructure(Mat,MatStructure*);
332: PETSC_EXTERN PetscErrorCode MatCompositeGetNumberMat(Mat,PetscInt*);
333: PETSC_EXTERN PetscErrorCode MatCompositeGetMat(Mat,PetscInt,Mat*);
334: PETSC_EXTERN PetscErrorCode MatCompositeSetScalings(Mat,const PetscScalar*);
336: PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*);
337: PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
339: PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
340: PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat,Mat*);
341: PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*);
342: PETSC_EXTERN PetscErrorCode MatHermitianTransposeGetMat(Mat,Mat*);
343: PETSC_EXTERN PetscErrorCode MatCreateSubMatrixVirtual(Mat,IS,IS,Mat*);
344: PETSC_EXTERN PetscErrorCode MatSubMatrixVirtualUpdate(Mat,Mat,IS,IS);
345: PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
346: PETSC_EXTERN PetscErrorCode MatCreateConstantDiagonal(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar,Mat*);
348: #if defined(PETSC_HAVE_HYPRE)
349: PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
350: #endif
352: PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]);
354: PETSC_EXTERN PetscErrorCode MatResetPreallocation(Mat);
355: PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
356: PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
357: PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*);
359: PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
360: PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
361: PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
362: PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
363: PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
364: PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
365: PETSC_EXTERN PetscErrorCode MatInvertVariableBlockDiagonal(Mat,PetscInt,const PetscInt*,PetscScalar*);
366: PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonalMat(Mat,Mat);
368: /* ------------------------------------------------------------*/
369: PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
370: PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
371: PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
372: PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
373: PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
374: PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);
376: /*S
377: MatStencil - Data structure (C struct) for storing information about a single row or
378: column of a matrix as indexed on an associated grid. These are arguments to MatSetStencil() and MatSetBlockStencil()
380: 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).
381: 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.
383: For stencil access to vectors see DMDAVecGetArray(), DMDAVecGetArrayF90().
385: Fortran usage is different, see MatSetValuesStencil() for details.
387: Level: beginner
389: .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90()
390: S*/
391: typedef struct {
392: PetscInt k,j,i,c;
393: } MatStencil;
395: PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
396: PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
397: PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
399: /*E
400: MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
401: to continue to add or insert values to it
403: Level: beginner
405: .seealso: MatAssemblyBegin(), MatAssemblyEnd()
406: E*/
407: typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
408: PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
409: PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
410: PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);
414: /*E
415: MatOption - Options that may be set for a matrix and its behavior or storage
417: Level: beginner
419: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
420: Any additions/changes here must also be made in src/mat/interface/dlregismat.c in MatOptions[]
422: Developer Notes:
423: Entries that are negative need not be called collectively by all processes.
425: .seealso: MatSetOption()
426: E*/
427: typedef enum {MAT_OPTION_MIN = -3,
428: MAT_UNUSED_NONZERO_LOCATION_ERR = -2,
429: MAT_ROW_ORIENTED = -1,
430: MAT_SYMMETRIC = 1,
431: MAT_STRUCTURALLY_SYMMETRIC = 2,
432: MAT_NEW_DIAGONALS = 3,
433: MAT_IGNORE_OFF_PROC_ENTRIES = 4,
434: MAT_USE_HASH_TABLE = 5,
435: MAT_KEEP_NONZERO_PATTERN = 6,
436: MAT_IGNORE_ZERO_ENTRIES = 7,
437: MAT_USE_INODES = 8,
438: MAT_HERMITIAN = 9,
439: MAT_SYMMETRY_ETERNAL = 10,
440: MAT_NEW_NONZERO_LOCATION_ERR = 11,
441: MAT_IGNORE_LOWER_TRIANGULAR = 12,
442: MAT_ERROR_LOWER_TRIANGULAR = 13,
443: MAT_GETROW_UPPERTRIANGULAR = 14,
444: MAT_SPD = 15,
445: MAT_NO_OFF_PROC_ZERO_ROWS = 16,
446: MAT_NO_OFF_PROC_ENTRIES = 17,
447: MAT_NEW_NONZERO_LOCATIONS = 18,
448: MAT_NEW_NONZERO_ALLOCATION_ERR = 19,
449: MAT_SUBSET_OFF_PROC_ENTRIES = 20,
450: MAT_SUBMAT_SINGLEIS = 21,
451: MAT_STRUCTURE_ONLY = 22,
452: MAT_SORTED_FULL = 23,
453: MAT_OPTION_MAX = 24} MatOption;
455: PETSC_EXTERN const char *const *MatOptions;
456: PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool);
457: PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*);
458: PETSC_EXTERN PetscErrorCode MatPropagateSymmetryOptions(Mat,Mat);
459: PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);
461: PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
462: PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
463: PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
464: PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
465: PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
466: PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
467: PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
468: PETSC_EXTERN PetscErrorCode MatSeqAIJGetArrayRead(Mat,const PetscScalar *[]);
469: PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
470: PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArrayRead(Mat,const PetscScalar *[]);
471: PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*);
472: PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
473: PETSC_EXTERN PetscErrorCode MatSeqAIJSetType(Mat,MatType);
474: PETSC_EXTERN PetscErrorCode MatSeqAIJRegister(const char[],PetscErrorCode (*)(Mat,MatType,MatReuse,Mat *));
475: PETSC_EXTERN PetscFunctionList MatSeqAIJList;
476: PETSC_EXTERN PetscErrorCode MatSeqBAIJGetArray(Mat,PetscScalar *[]);
477: PETSC_EXTERN PetscErrorCode MatSeqBAIJRestoreArray(Mat,PetscScalar *[]);
478: PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat,PetscScalar *[]);
479: PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat,PetscScalar *[]);
480: PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
481: PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
482: PETSC_EXTERN PetscErrorCode MatDensePlaceArray(Mat,const PetscScalar[]);
483: PETSC_EXTERN PetscErrorCode MatDenseReplaceArray(Mat,const PetscScalar[]);
484: PETSC_EXTERN PetscErrorCode MatDenseResetArray(Mat);
485: PETSC_EXTERN PetscErrorCode MatDenseGetArrayRead(Mat,const PetscScalar *[]);
486: PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayRead(Mat,const PetscScalar *[]);
487: PETSC_EXTERN PetscErrorCode MatDenseGetArrayWrite(Mat,PetscScalar *[]);
488: PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayWrite(Mat,PetscScalar *[]);
489: PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
490: PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
491: PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
492: PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
493: PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat);
494: PETSC_EXTERN PetscErrorCode MatSetVariableBlockSizes(Mat,PetscInt,PetscInt*);
495: PETSC_EXTERN PetscErrorCode MatGetVariableBlockSizes(Mat,PetscInt*,const PetscInt**);
497: PETSC_EXTERN PetscErrorCode MatDenseGetColumn(Mat,PetscInt,PetscScalar*[]);
498: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumn(Mat,PetscScalar*[]);
499: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVec(Mat,PetscInt,Vec*);
500: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVec(Mat,PetscInt,Vec*);
501: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecRead(Mat,PetscInt,Vec*);
502: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecRead(Mat,PetscInt,Vec*);
503: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecWrite(Mat,PetscInt,Vec*);
504: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecWrite(Mat,PetscInt,Vec*);
505: PETSC_EXTERN PetscErrorCode MatDenseGetSubMatrix(Mat,PetscInt,PetscInt,Mat*);
506: PETSC_EXTERN PetscErrorCode MatDenseRestoreSubMatrix(Mat,Mat*);
508: PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
509: PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
510: PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
511: PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
512: PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
513: PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
514: PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
515: PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
516: PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
517: PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
518: PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
519: PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
520: PETSC_EXTERN PetscErrorCode MatMatSolveTranspose(Mat,Mat,Mat);
521: PETSC_EXTERN PetscErrorCode MatMatTransposeSolve(Mat,Mat,Mat);
522: PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);
524: /*E
525: MatDuplicateOption - Indicates if a duplicated sparse matrix should have
526: its numerical values copied over or just its nonzero structure.
528: Level: beginner
530: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
532: $ MAT_DO_NOT_COPY_VALUES - Create a matrix using the same nonzero pattern as the original matrix,
533: $ with zeros for the numerical values.
534: $ MAT_COPY_VALUES - Create a matrix with the same nonzero pattern as the original matrix
535: $ and with the same numerical values.
536: $ MAT_SHARE_NONZERO_PATTERN - Create a matrix that shares the nonzero structure with the previous matrix
537: $ and does not copy it, using zeros for the numerical values. The parent and
538: $ child matrices will share their index (i and j) arrays, and you cannot
539: $ insert new nonzero entries into either matrix.
541: Notes:
542: Many matrix types (including SeqAIJ) do not support the MAT_SHARE_NONZERO_PATTERN optimization; in
543: this case the behavior is as if MAT_DO_NOT_COPY_VALUES has been specified.
545: .seealso: MatDuplicate()
546: E*/
547: typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
549: PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
550: PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);
553: PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
554: PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
555: PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
556: PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
557: PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
558: PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
559: PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
560: PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *);
561: PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
563: PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *);
564: PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *);
565: PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *);
566: PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *);
568: /*S
569: MatInfo - Context of matrix information, used with MatGetInfo()
571: In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
573: Level: intermediate
575: .seealso: MatGetInfo(), MatInfoType
576: S*/
577: typedef struct {
578: PetscLogDouble block_size; /* block size */
579: PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */
580: PetscLogDouble memory; /* memory allocated */
581: PetscLogDouble assemblies; /* number of matrix assemblies called */
582: PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */
583: PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
584: PetscLogDouble factor_mallocs; /* number of mallocs during factorization */
585: } MatInfo;
587: /*E
588: MatInfoType - Indicates if you want information about the local part of the matrix,
589: the entire parallel matrix or the maximum over all the local parts.
591: Level: beginner
593: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
595: .seealso: MatGetInfo(), MatInfo
596: E*/
597: typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
598: PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
599: PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
600: PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
601: PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
602: PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
603: PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
604: PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
605: PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
606: PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
607: PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*);
608: PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
609: PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);
611: PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*);
612: PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*);
613: PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*);
614: PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*);
615: PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*);
616: PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
617: PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
618: PETSC_EXTERN PetscErrorCode MatMatTransposeMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
619: PETSC_EXTERN PetscErrorCode MatPtAPMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
620: PETSC_EXTERN PetscErrorCode MatRARtMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
621: PETSC_EXTERN PetscErrorCode MatIsLinear(Mat,PetscInt,PetscBool*);
623: PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*);
624: PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*);
625: PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
626: PETSC_EXTERN PetscErrorCode MatSetInf(Mat);
627: PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
628: PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
629: PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
630: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
631: PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
632: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
634: PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
635: PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
636: PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
637: PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
638: PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
639: PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
640: PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);
642: PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
643: PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatrices() (since version 3.8)") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatrices(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) {return MatCreateSubMatrices(mat,n,irow,icol,scall,submat);}
644: PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
645: PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatricesMPI() (since version 3.8)") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatricesMPI(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) {return MatCreateSubMatricesMPI(mat,n,irow,icol,scall,submat);}
646: PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
647: PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]);
648: PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *);
649: PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatrix() (since version 3.8)") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatrix(Mat mat,IS isrow,IS iscol,MatReuse cll,Mat *newmat) {return MatCreateSubMatrix(mat,isrow,iscol,cll,newmat);}
650: PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
651: PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
652: PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
653: PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);
655: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
656: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
657: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
658: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
659: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
660: PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
661: PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
663: PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
664: PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov);
665: PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool);
667: PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
669: PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
670: PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
672: PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
673: PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
675: PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
676: PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
678: PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
679: PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);
681: PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
682: PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);
684: PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
685: PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
686: PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*);
687: PETSC_EXTERN PetscErrorCode MatSetLayouts(Mat,PetscLayout,PetscLayout);
688: PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
689: PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
690: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
691: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
692: PETSC_EXTERN PetscErrorCode MatGetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
693: PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
694: PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
696: PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
697: PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
699: PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
700: PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
701: PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
702: PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*);
703: PETSC_DEPRECATED_FUNCTION("Use MatCreateVecs() (since version 3.6)") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);}
704: PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
705: PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
706: PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
707: PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
708: PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
710: /*MC
711: MatSetValue - Set a single entry into a matrix.
713: Not collective
715: Synopsis:
716: #include <petscmat.h>
717: PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode)
719: Input Parameters:
720: + m - the matrix
721: . row - the row location of the entry
722: . col - the column location of the entry
723: . value - the value to insert
724: - mode - either INSERT_VALUES or ADD_VALUES
726: Notes:
727: For efficiency one should use MatSetValues() and set several or many
728: values simultaneously if possible.
730: Level: beginner
732: .seealso: MatSetValues(), MatSetValueLocal()
733: M*/
734: PETSC_STATIC_INLINE PetscErrorCode MatSetValue(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValues(v,1,&i,1,&j,&va,mode);}
736: PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
738: PETSC_STATIC_INLINE PetscErrorCode MatSetValueLocal(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValuesLocal(v,1,&i,1,&j,&va,mode);}
740: /*MC
741: MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
742: row in a matrix providing the data that one can use to correctly preallocate the matrix.
744: Synopsis:
745: #include <petscmat.h>
746: PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
748: Collective
750: Input Parameters:
751: + comm - the communicator that will share the eventually allocated matrix
752: . nrows - the number of LOCAL rows in the matrix
753: - ncols - the number of LOCAL columns in the matrix
755: Output Parameters:
756: + dnz - the array that will be passed to the matrix preallocation routines
757: - onz - the other array passed to the matrix preallocation routines
759: Level: intermediate
761: Notes:
762: See Users-Manual: ch_performance for more details.
764: Do not malloc or free dnz and onz, that is handled internally by these routines
766: This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
768: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
769: MatPreallocateSymmetricSetLocalBlock()
770: M*/
771: #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
772: do { \
773: PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ncols = (ncols),__rstart,__start,__end; \
774: _4_PetscCalloc2((size_t)__nrows,&dnz,(size_t)__nrows,&onz);CHKERRQ(_4_ierr); \
775: __start = 0; __end = __start; \
776: _4_MPI_Scan(&__ncols,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ncols;\
777: _4_MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; \
778: do { } while (0)
780: /*MC
781: MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
782: inserted using a local number of the rows and columns
784: Synopsis:
785: #include <petscmat.h>
786: PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
788: Not Collective
790: Input Parameters:
791: + map - the row mapping from local numbering to global numbering
792: . nrows - the number of rows indicated
793: . rows - the indices of the rows
794: . cmap - the column mapping from local to global numbering
795: . ncols - the number of columns in the matrix
796: . cols - the columns indicated
797: . dnz - the array that will be passed to the matrix preallocation routines
798: - onz - the other array passed to the matrix preallocation routines
800: Level: intermediate
802: Notes:
803: See Users-Manual: ch_performance for more details.
805: Do not malloc or free dnz and onz, that is handled internally by these routines
807: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
808: MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups()
809: M*/
810: #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
811: do {\
812: PetscInt __l;\
813: _4_ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
814: _4_ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
815: for (__l=0;__l<nrows;__l++) {\
816: _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
817: }\
818: } while (0)
820: /*MC
821: MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be
822: inserted using a local number of the rows and columns. This version removes any duplicate columns in cols
824: Synopsis:
825: #include <petscmat.h>
826: PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
828: Not Collective
830: Input Parameters:
831: + map - the row mapping from local numbering to global numbering
832: . nrows - the number of rows indicated
833: . rows - the indices of the rows (these values are mapped to the global values)
834: . cmap - the column mapping from local to global numbering
835: . ncols - the number of columns in the matrix (this value will be changed if duplicate columns are found)
836: . cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed)
837: . dnz - the array that will be passed to the matrix preallocation routines
838: - onz - the other array passed to the matrix preallocation routines
840: Level: intermediate
842: Notes:
843: See Users-Manual: ch_performance for more details.
845: Do not malloc or free dnz and onz, that is handled internally by these routines
847: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
848: MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
849: M*/
850: #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
851: do {\
852: PetscInt __l;\
853: _4_ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
854: _4_ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
855: _4_PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\
856: for (__l=0;__l<nrows;__l++) {\
857: _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
858: }\
859: } while (0)
861: /*MC
862: MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
863: inserted using a local number of the rows and columns
865: Synopsis:
866: #include <petscmat.h>
867: PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
869: Not Collective
871: Input Parameters:
872: + map - the row mapping from local numbering to global numbering
873: . nrows - the number of rows indicated
874: . rows - the indices of the rows
875: . cmap - the column mapping from local to global numbering
876: . ncols - the number of columns in the matrix
877: . cols - the columns indicated
878: . dnz - the array that will be passed to the matrix preallocation routines
879: - onz - the other array passed to the matrix preallocation routines
881: Level: intermediate
883: Notes:
884: See Users-Manual: ch_performance for more details.
886: Do not malloc or free dnz and onz, that is handled internally by these routines
888: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
889: MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
890: M*/
891: #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
892: do {\
893: PetscInt __l;\
894: _4_ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
895: _4_ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
896: for (__l=0;__l<nrows;__l++) {\
897: _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
898: }\
899: } while (0)
901: /*MC
902: MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
903: inserted using a local number of the rows and columns
905: Synopsis:
906: #include <petscmat.h>
907: PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
909: Not Collective
911: Input Parameters:
912: + map - the mapping between local numbering and global numbering
913: . nrows - the number of rows indicated
914: . rows - the indices of the rows
915: . ncols - the number of columns in the matrix
916: . cols - the columns indicated
917: . dnz - the array that will be passed to the matrix preallocation routines
918: - onz - the other array passed to the matrix preallocation routines
920: Level: intermediate
922: Notes:
923: See Users-Manual: ch_performance for more details.
925: Do not malloc or free dnz and onz that is handled internally by these routines
927: .seealso: MatPreallocateFinalize(), MatPreallocateSet()
928: MatPreallocateInitialize(), MatPreallocateSetLocal()
929: M*/
930: #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\
931: do {\
932: PetscInt __l;\
933: _4_ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
934: _4_ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
935: for (__l=0;__l<nrows;__l++) {\
936: _4_MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
937: }\
938: } while (0)
940: /*MC
941: MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
942: inserted using a local number of the rows and columns
944: Synopsis:
945: #include <petscmat.h>
946: PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
948: Not Collective
950: Input Parameters:
951: + row - the row
952: . ncols - the number of columns in the matrix
953: - cols - the columns indicated
955: Output Parameters:
956: + dnz - the array that will be passed to the matrix preallocation routines
957: - onz - the other array passed to the matrix preallocation routines
959: Level: intermediate
961: Notes:
962: See Users-Manual: ch_performance for more details.
964: Do not malloc or free dnz and onz that is handled internally by these routines
966: This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
968: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
969: MatPreallocateInitialize(), MatPreallocateSetLocal()
970: M*/
971: #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
972: do { PetscInt __i; \
973: if (row < __rstart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
974: if (row >= __rstart+__nrows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D greater than last local row %D",row,__rstart+__nrows-1);\
975: for (__i=0; __i<nc; __i++) {\
976: if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
977: else if (dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\
978: }\
979: } while (0)
981: /*MC
982: MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
983: inserted using a local number of the rows and columns
985: Synopsis:
986: #include <petscmat.h>
987: PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
989: Not Collective
991: Input Parameters:
992: + nrows - the number of rows indicated
993: . rows - the indices of the rows
994: . ncols - the number of columns in the matrix
995: . cols - the columns indicated
996: . dnz - the array that will be passed to the matrix preallocation routines
997: - onz - the other array passed to the matrix preallocation routines
999: Level: intermediate
1001: Notes:
1002: See Users-Manual: ch_performance for more details.
1004: Do not malloc or free dnz and onz that is handled internally by these routines
1006: This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
1008: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(),
1009: MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
1010: M*/
1011: #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\
1012: do { PetscInt __i; \
1013: for (__i=0; __i<nc; __i++) {\
1014: if (cols[__i] >= __end) onz[row - __rstart]++; \
1015: else if (cols[__i] >= row && dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\
1016: }\
1017: } while (0)
1019: /*MC
1020: MatPreallocateLocation - An alternative to MatPreallocateSet() that puts the nonzero locations into the matrix if it exists
1022: Synopsis:
1023: #include <petscmat.h>
1024: PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
1026: Not Collective
1028: Input Parameters:
1029: + A - matrix
1030: . row - row where values exist (must be local to this process)
1031: . ncols - number of columns
1032: . cols - columns with nonzeros
1033: . dnz - the array that will be passed to the matrix preallocation routines
1034: - onz - the other array passed to the matrix preallocation routines
1036: Level: intermediate
1038: Notes:
1039: See Users-Manual: ch_performance for more details.
1041: Do not malloc or free dnz and onz that is handled internally by these routines
1043: This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
1045: .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
1046: MatPreallocateSymmetricSetLocalBlock()
1047: M*/
1048: #define MatPreallocateLocation(A,row,ncols,cols,dnz,onz) 0; do {if (A) {MatSetValues(A,1,&row,ncols,cols,NULL,INSERT_VALUES);} else { MatPreallocateSet(row,ncols,cols,dnz,onz);}} while (0)
1051: /*MC
1052: MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
1053: row in a matrix providing the data that one can use to correctly preallocate the matrix.
1055: Synopsis:
1056: #include <petscmat.h>
1057: PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
1059: Collective
1061: Input Parameters:
1062: + dnz - the array that was be passed to the matrix preallocation routines
1063: - onz - the other array passed to the matrix preallocation routines
1065: Level: intermediate
1067: Notes:
1068: See Users-Manual: ch_performance for more details.
1070: Do not malloc or free dnz and onz that is handled internally by these routines
1072: This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
1074: .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
1075: MatPreallocateSymmetricSetLocalBlock()
1076: M*/
1077: #define MatPreallocateFinalize(dnz,onz) 0;_4_PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} while (0)
1079: /* Routines unique to particular data structures */
1080: PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
1082: PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
1083: PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
1085: PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
1086: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
1087: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1088: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1089: PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1090: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
1092: #define MAT_SKIP_ALLOCATION -4
1094: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1095: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1096: PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
1097: PETSC_EXTERN PetscErrorCode MatSeqAIJSetTotalPreallocation(Mat,PetscInt);
1099: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1100: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1101: PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1102: PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
1103: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1104: PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
1105: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1106: PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
1107: PETSC_EXTERN PetscErrorCode MatMPIAdjToSeq(Mat,Mat*);
1108: PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
1109: PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
1110: PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1111: PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1112: PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
1114: PETSC_EXTERN PetscErrorCode MatDenseGetLDA(Mat,PetscInt*);
1115: PETSC_EXTERN PetscErrorCode MatDenseSetLDA(Mat,PetscInt);
1116: PETSC_DEPRECATED_FUNCTION("Use MatDenseSetLDA() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode MatSeqDenseSetLDA(Mat A,PetscInt lda) {return MatDenseSetLDA(A,lda);}
1117: PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
1119: PETSC_EXTERN PetscErrorCode MatBlockMatSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1121: PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
1122: PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
1124: PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
1126: PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
1127: PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*);
1128: /*
1129: These routines are not usually accessed directly, rather solving is
1130: done through the KSP and PC interfaces.
1131: */
1133: /*J
1134: MatOrderingType - String with the name of a PETSc matrix ordering
1136: Level: beginner
1138: Notes:
1139: If MATORDERINGEXTERNAL is used then PETSc does not compute an ordering and utilizes one built into the factorization package
1141: .seealso: MatGetOrdering()
1142: J*/
1143: typedef const char* MatOrderingType;
1144: #define MATORDERINGNATURAL "natural"
1145: #define MATORDERINGND "nd"
1146: #define MATORDERING1WD "1wd"
1147: #define MATORDERINGRCM "rcm"
1148: #define MATORDERINGQMD "qmd"
1149: #define MATORDERINGROWLENGTH "rowlength"
1150: #define MATORDERINGWBM "wbm"
1151: #define MATORDERINGSPECTRAL "spectral"
1152: #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */
1153: #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 */
1154: #define MATORDERINGEXTERNAL "external" /* uses an ordering type internal to the factorization package */
1156: PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
1157: PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
1158: PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
1159: PETSC_EXTERN PetscFunctionList MatOrderingList;
1161: PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1162: PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);
1164: /*S
1165: MatFactorShiftType - Numeric Shift for factorizations
1167: Level: beginner
1169: .seealso: MatGetFactor()
1170: S*/
1171: typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1172: PETSC_EXTERN const char *const MatFactorShiftTypes[];
1173: PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];
1175: /*S
1176: MatFactorError - indicates what type of error was generated in a matrix factorization
1178: Level: beginner
1180: Developer Notes:
1181: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1183: .seealso: MatGetFactor()
1184: S*/
1185: typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError;
1187: PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*);
1188: PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat);
1189: PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*);
1191: /*S
1192: MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization
1194: In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
1195: $ MatFactorInfo info(MAT_FACTORINFO_SIZE)
1197: Notes:
1198: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1200: You can use MatFactorInfoInitialize() to set default values.
1202: Level: developer
1204: .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1205: MatFactorInfoInitialize()
1207: S*/
1208: typedef struct {
1209: PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */
1210: PetscReal usedt;
1211: PetscReal dt; /* drop tolerance */
1212: PetscReal dtcol; /* tolerance for pivoting */
1213: PetscReal dtcount; /* maximum nonzeros to be allowed per row */
1214: PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1215: PetscReal levels; /* ICC/ILU(levels) */
1216: PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1217: factorization may be faster if do not pivot */
1218: PetscReal zeropivot; /* pivot is called zero if less than this */
1219: PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */
1220: PetscReal shiftamount; /* how large the shift is */
1221: } MatFactorInfo;
1223: PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
1224: PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
1225: PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1226: PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
1227: PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
1228: PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
1229: PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1230: PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1231: PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1232: PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
1233: PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1234: PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1235: PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1236: PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1237: PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1238: PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1239: PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1240: PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1241: PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
1242: PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
1244: typedef enum {MAT_FACTOR_SCHUR_UNFACTORED, MAT_FACTOR_SCHUR_FACTORED, MAT_FACTOR_SCHUR_INVERTED} MatFactorSchurStatus;
1245: PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS);
1246: PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1247: PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*,MatFactorSchurStatus);
1248: PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat);
1249: PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1250: PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec);
1251: PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec);
1252: PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat);
1254: /*E
1255: MatSORType - What type of (S)SOR to perform
1257: Level: beginner
1259: May be bitwise ORd together
1261: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1263: MatSORType may be bitwise ORd together, so do not change the numbers
1265: .seealso: MatSOR()
1266: E*/
1267: typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1268: SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1269: SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
1270: SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1271: PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1273: /*
1274: These routines are for efficiently computing Jacobians via finite differences.
1275: */
1277: /*S
1278: MatColoring - Object for managing the coloring of matrices.
1280: Level: beginner
1282: Notes:
1283: Coloring of matrices can be computed directly from the sparse matrix nonzero structure via the MatColoring object or from the mesh from which the
1284: matrix comes from via DMCreateColoring(). In general using the mesh produces a more optimal coloring (fewer colors).
1286: Once a coloring is available MatFDColoringCreate() creates an object that can be used to efficiently compute Jacobians using that coloring. This
1287: same object can also be used to efficiently convert data created by Automatic Differentation tools to PETSc sparse matrices.
1289: .seealso: MatFDColoringCreate(), MatColoringWeightType, ISColoring, MatFDColoring, DMCreateColoring(), MatColoringCreate(), MatOrdering, MatPartitioning
1290: S*/
1291: typedef struct _p_MatColoring* MatColoring;
1293: /*J
1294: MatColoringType - String with the name of a PETSc matrix coloring
1296: Level: beginner
1298: .seealso: MatColoringSetType(), MatColoring
1299: J*/
1300: typedef const char* MatColoringType;
1301: #define MATCOLORINGJP "jp"
1302: #define MATCOLORINGPOWER "power"
1303: #define MATCOLORINGNATURAL "natural"
1304: #define MATCOLORINGSL "sl"
1305: #define MATCOLORINGLF "lf"
1306: #define MATCOLORINGID "id"
1307: #define MATCOLORINGGREEDY "greedy"
1309: /*E
1310: MatColoringWeightType - Type of weight scheme
1312: Not Collective
1314: + MAT_COLORING_RANDOM - Random weights
1315: . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
1316: - MAT_COLORING_LF - Last-first weighting.
1318: Level: intermediate
1320: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1322: .seealso: MatColoring, MatColoringCreate()
1323: E*/
1324: typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;
1326: PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1327: PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*);
1328: PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1329: PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1330: PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1331: PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1332: PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1333: PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1334: PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1335: PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1336: PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1337: PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1338: PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1339: PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType);
1340: PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*);
1341: PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm);
1342: PETSC_EXTERN PetscErrorCode MatColoringTest(MatColoring,ISColoring);
1343: PETSC_DEPRECATED_FUNCTION("Use MatColoringTest() (since version 3.10)") PETSC_STATIC_INLINE PetscErrorCode MatColoringTestValid(MatColoring matcoloring,ISColoring iscoloring) {return MatColoringTest(matcoloring,iscoloring);}
1344: PETSC_EXTERN PetscErrorCode MatISColoringTest(Mat,ISColoring);
1346: /*S
1347: MatFDColoring - Object for computing a sparse Jacobian via finite differences
1348: and coloring
1350: Level: beginner
1352: Notes:
1353: This object is creating utilizing a coloring provided by the MatColoring object or DMCreateColoring()
1355: .seealso: MatFDColoringCreate(), MatColoring, DMCreateColoring()
1356: S*/
1357: typedef struct _p_MatFDColoring* MatFDColoring;
1359: PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1360: PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1361: PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1362: PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1363: PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1364: PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1365: PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1366: PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *);
1367: PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1368: PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,const PetscInt*[]);
1369: PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1370: PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1371: PETSC_EXTERN PetscErrorCode MatFDColoringSetValues(Mat,MatFDColoring,const PetscScalar*);
1373: /*S
1374: MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1376: Level: beginner
1378: .seealso: MatTransposeColoringCreate()
1379: S*/
1380: typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1382: PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1383: PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1384: PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1385: PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1387: /*
1388: These routines are for partitioning matrices: currently used only
1389: for adjacency matrix, MatCreateMPIAdj().
1390: */
1392: /*S
1393: MatPartitioning - Object for managing the partitioning of a matrix or graph
1395: Level: beginner
1397: Notes:
1398: There is also a PetscPartitioner object that provides the same functionality. It can utilize the MatPartitioning operations
1399: via PetscPartitionerSetType(p,PETSCPARTITIONERMATPARTITIONING)
1401: Developers Note:
1402: It is an extra maintainance and documentation cost to have two objects with the same functionality.
1404: .seealso: MatPartitioningCreate(), MatPartitioningType, MatColoring, MatGetOrdering()
1405: S*/
1406: typedef struct _p_MatPartitioning* MatPartitioning;
1408: /*J
1409: MatPartitioningType - String with the name of a PETSc matrix partitioning
1411: Level: beginner
1412: dm
1413: .seealso: MatPartitioningCreate(), MatPartitioning
1414: J*/
1415: typedef const char* MatPartitioningType;
1416: #define MATPARTITIONINGCURRENT "current"
1417: #define MATPARTITIONINGAVERAGE "average"
1418: #define MATPARTITIONINGSQUARE "square"
1419: #define MATPARTITIONINGPARMETIS "parmetis"
1420: #define MATPARTITIONINGCHACO "chaco"
1421: #define MATPARTITIONINGPARTY "party"
1422: #define MATPARTITIONINGPTSCOTCH "ptscotch"
1423: #define MATPARTITIONINGHIERARCH "hierarch"
1425: PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1426: PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1427: PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1428: PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1429: PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1430: PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1431: PETSC_EXTERN PetscErrorCode MatPartitioningSetUseEdgeWeights(MatPartitioning,PetscBool);
1432: PETSC_EXTERN PetscErrorCode MatPartitioningGetUseEdgeWeights(MatPartitioning,PetscBool*);
1433: PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1434: PETSC_EXTERN PetscErrorCode MatPartitioningImprove(MatPartitioning,IS*);
1435: PETSC_EXTERN PetscErrorCode MatPartitioningViewImbalance(MatPartitioning,IS);
1436: PETSC_EXTERN PetscErrorCode MatPartitioningApplyND(MatPartitioning,IS*);
1437: PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
1438: PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
1439: PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1440: PETSC_EXTERN PetscErrorCode MatPartitioningViewFromOptions(MatPartitioning,PetscObject,const char[]);
1441: PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
1442: PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1444: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1445: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1446: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
1448: typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
1449: PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1450: typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
1451: PETSC_EXTERN const char *const MPChacoLocalTypes[];
1452: typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
1453: PETSC_EXTERN const char *const MPChacoEigenTypes[];
1455: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1456: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1457: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1458: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1459: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1460: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1461: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1462: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1463: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1464: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1465: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
1467: #define MP_PARTY_OPT "opt"
1468: #define MP_PARTY_LIN "lin"
1469: #define MP_PARTY_SCA "sca"
1470: #define MP_PARTY_RAN "ran"
1471: #define MP_PARTY_GBF "gbf"
1472: #define MP_PARTY_GCF "gcf"
1473: #define MP_PARTY_BUB "bub"
1474: #define MP_PARTY_DEF "def"
1475: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
1476: #define MP_PARTY_HELPFUL_SETS "hs"
1477: #define MP_PARTY_KERNIGHAN_LIN "kl"
1478: #define MP_PARTY_NONE "no"
1479: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1480: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1481: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1482: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
1484: typedef enum { MP_PTSCOTCH_DEFAULT,MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
1485: PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1487: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1488: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1489: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1490: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
1492: /*
1493: * hierarchical partitioning
1494: */
1495: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*);
1496: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*);
1497: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt);
1498: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);
1500: PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1501: PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1503: /*
1504: If you add entries here you must also add them to include/petscmat.h
1505: and src/mat/f90-mod/petscmat.h
1506: */
1507: typedef enum { MATOP_SET_VALUES=0,
1508: MATOP_GET_ROW=1,
1509: MATOP_RESTORE_ROW=2,
1510: MATOP_MULT=3,
1511: MATOP_MULT_ADD=4,
1512: MATOP_MULT_TRANSPOSE=5,
1513: MATOP_MULT_TRANSPOSE_ADD=6,
1514: MATOP_SOLVE=7,
1515: MATOP_SOLVE_ADD=8,
1516: MATOP_SOLVE_TRANSPOSE=9,
1517: MATOP_SOLVE_TRANSPOSE_ADD=10,
1518: MATOP_LUFACTOR=11,
1519: MATOP_CHOLESKYFACTOR=12,
1520: MATOP_SOR=13,
1521: MATOP_TRANSPOSE=14,
1522: MATOP_GETINFO=15,
1523: MATOP_EQUAL=16,
1524: MATOP_GET_DIAGONAL=17,
1525: MATOP_DIAGONAL_SCALE=18,
1526: MATOP_NORM=19,
1527: MATOP_ASSEMBLY_BEGIN=20,
1528: MATOP_ASSEMBLY_END=21,
1529: MATOP_SET_OPTION=22,
1530: MATOP_ZERO_ENTRIES=23,
1531: MATOP_ZERO_ROWS=24,
1532: MATOP_LUFACTOR_SYMBOLIC=25,
1533: MATOP_LUFACTOR_NUMERIC=26,
1534: MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1535: MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1536: MATOP_SETUP_PREALLOCATION=29,
1537: MATOP_ILUFACTOR_SYMBOLIC=30,
1538: MATOP_ICCFACTOR_SYMBOLIC=31,
1539: MATOP_GET_DIAGONAL_BLOCK=32,
1540: MATOP_FREE_INTER_STRUCT=33,
1541: MATOP_DUPLICATE=34,
1542: MATOP_FORWARD_SOLVE=35,
1543: MATOP_BACKWARD_SOLVE=36,
1544: MATOP_ILUFACTOR=37,
1545: MATOP_ICCFACTOR=38,
1546: MATOP_AXPY=39,
1547: MATOP_CREATE_SUBMATRICES=40,
1548: MATOP_INCREASE_OVERLAP=41,
1549: MATOP_GET_VALUES=42,
1550: MATOP_COPY=43,
1551: MATOP_GET_ROW_MAX=44,
1552: MATOP_SCALE=45,
1553: MATOP_SHIFT=46,
1554: MATOP_DIAGONAL_SET=47,
1555: MATOP_ZERO_ROWS_COLUMNS=48,
1556: MATOP_SET_RANDOM=49,
1557: MATOP_GET_ROW_IJ=50,
1558: MATOP_RESTORE_ROW_IJ=51,
1559: MATOP_GET_COLUMN_IJ=52,
1560: MATOP_RESTORE_COLUMN_IJ=53,
1561: MATOP_FDCOLORING_CREATE=54,
1562: MATOP_COLORING_PATCH=55,
1563: MATOP_SET_UNFACTORED=56,
1564: MATOP_PERMUTE=57,
1565: MATOP_SET_VALUES_BLOCKED=58,
1566: MATOP_CREATE_SUBMATRIX=59,
1567: MATOP_DESTROY=60,
1568: MATOP_VIEW=61,
1569: MATOP_CONVERT_FROM=62,
1570: MATOP_MATMAT_MULT=63,
1571: MATOP_MATMAT_MULT_SYMBOLIC=64,
1572: MATOP_MATMAT_MULT_NUMERIC=65,
1573: MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1574: MATOP_SET_VALUES_LOCAL=67,
1575: MATOP_ZERO_ROWS_LOCAL=68,
1576: MATOP_GET_ROW_MAX_ABS=69,
1577: MATOP_GET_ROW_MIN_ABS=70,
1578: MATOP_CONVERT=71,
1579: MATOP_SET_COLORING=72,
1580: /* MATOP_PLACEHOLDER_73=73, */
1581: MATOP_SET_VALUES_ADIFOR=74,
1582: MATOP_FD_COLORING_APPLY=75,
1583: MATOP_SET_FROM_OPTIONS=76,
1584: MATOP_MULT_CONSTRAINED=77,
1585: MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
1586: MATOP_FIND_ZERO_DIAGONALS=79,
1587: MATOP_MULT_MULTIPLE=80,
1588: MATOP_SOLVE_MULTIPLE=81,
1589: MATOP_GET_INERTIA=82,
1590: MATOP_LOAD=83,
1591: MATOP_IS_SYMMETRIC=84,
1592: MATOP_IS_HERMITIAN=85,
1593: MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1594: MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1595: MATOP_CREATE_VECS=88,
1596: MATOP_MAT_MULT=89,
1597: MATOP_MAT_MULT_SYMBOLIC=90,
1598: MATOP_MAT_MULT_NUMERIC=91,
1599: MATOP_PTAP=92,
1600: MATOP_PTAP_SYMBOLIC=93,
1601: MATOP_PTAP_NUMERIC=94,
1602: MATOP_MAT_TRANSPOSE_MULT=95,
1603: MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1604: MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
1605: /* MATOP_PLACEHOLDER_98=98, */
1606: MATOP_PRODUCTSETFROMOPTIONS=99,
1607: MATOP_PRODUCTSYMBOLIC=100,
1608: MATOP_PRODUCTNUMERIC=101,
1609: MATOP_CONJUGATE=102,
1610: /* MATOP_PLACEHOLDER_103=103, */
1611: MATOP_SET_VALUES_ROW=104,
1612: MATOP_REAL_PART=105,
1613: MATOP_IMAGINARY_PART=106,
1614: MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1615: MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1616: MATOP_MAT_SOLVE=109,
1617: MATOP_MAT_SOLVE_TRANSPOSE=110,
1618: MATOP_GET_ROW_MIN=111,
1619: MATOP_GET_COLUMN_VECTOR=112,
1620: MATOP_MISSING_DIAGONAL=113,
1621: MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
1622: MATOP_CREATE=115,
1623: MATOP_GET_GHOSTS=116,
1624: MATOP_GET_LOCAL_SUB_MATRIX=117,
1625: MATOP_RESTORE_LOCALSUB_MATRIX=118,
1626: MATOP_MULT_DIAGONAL_BLOCK=119,
1627: MATOP_HERMITIAN_TRANSPOSE=120,
1628: MATOP_MULT_HERMITIAN_TRANSPOSE=121,
1629: MATOP_MULT_HERMITIAN_TRANS_ADD=122,
1630: MATOP_GET_MULTI_PROC_BLOCK=123,
1631: MATOP_FIND_NONZERO_ROWS=124,
1632: MATOP_GET_COLUMN_NORMS=125,
1633: MATOP_INVERT_BLOCK_DIAGONAL=126,
1634: /* MATOP_PLACEHOLDER_127=127, */
1635: MATOP_CREATE_SUB_MATRICES_MPI=128,
1636: MATOP_SET_VALUES_BATCH=129,
1637: MATOP_TRANSPOSE_MAT_MULT=130,
1638: MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
1639: MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
1640: MATOP_TRANSPOSE_COLORING_CREAT=133,
1641: MATOP_TRANS_COLORING_APPLY_SPT=134,
1642: MATOP_TRANS_COLORING_APPLY_DEN=135,
1643: MATOP_RART=136,
1644: MATOP_RART_SYMBOLIC=137,
1645: MATOP_RART_NUMERIC=138,
1646: MATOP_SET_BLOCK_SIZES=139,
1647: MATOP_AYPX=140,
1648: MATOP_RESIDUAL=141,
1649: MATOP_FDCOLORING_SETUP=142,
1650: MATOP_MPICONCATENATESEQ=144,
1651: MATOP_DESTROYSUBMATRICES=145,
1652: MATOP_TRANSPOSE_SOLVE=146,
1653: MATOP_GET_VALUES_LOCAL=147
1654: } MatOperation;
1655: PETSC_EXTERN PetscErrorCode MatSetOperation(Mat,MatOperation,void(*)(void));
1656: PETSC_EXTERN PetscErrorCode MatGetOperation(Mat,MatOperation,void(**)(void));
1657: PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool*);
1658: PETSC_EXTERN PetscErrorCode MatHasCongruentLayouts(Mat,PetscBool*);
1659: PETSC_DEPRECATED_FUNCTION("Use MatProductClear() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode MatFreeIntermediateDataStructures(Mat A) {return MatProductClear(A);}
1660: PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1661: PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1662: PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1663: PETSC_EXTERN PetscErrorCode MatShellSetVecType(Mat,VecType);
1664: PETSC_EXTERN PetscErrorCode MatShellTestMult(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*);
1665: PETSC_EXTERN PetscErrorCode MatShellTestMultTranspose(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*);
1666: PETSC_EXTERN PetscErrorCode MatShellSetManageScalingShifts(Mat);
1667: PETSC_EXTERN PetscErrorCode MatShellSetMatProductOperation(Mat,MatProductType,PetscErrorCode (*)(Mat,Mat,Mat,void**),PetscErrorCode (*)(Mat,Mat,Mat,void*),PetscErrorCode (*)(void*),MatType,MatType);
1668: PETSC_EXTERN PetscErrorCode MatIsShell(Mat,PetscBool*);
1670: /*
1671: Codes for matrices stored on disk. By default they are
1672: stored in a universal format. By changing the format with
1673: PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
1674: be stored in a way natural for the matrix, for example dense matrices
1675: would be stored as dense. Matrices stored this way may only be
1676: read into matrices of the same type.
1677: */
1678: #define MATRIX_BINARY_FORMAT_DENSE -1
1680: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1682: PETSC_EXTERN PetscErrorCode MatISSetLocalMatType(Mat,MatType);
1683: PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1684: PETSC_EXTERN PetscErrorCode MatISStoreL2L(Mat,PetscBool);
1685: PETSC_EXTERN PetscErrorCode MatISFixLocalEmpty(Mat,PetscBool);
1686: PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1687: PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat,Mat*);
1688: PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
1689: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use the MatConvert() interface (since version 3.10)") PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*);
1691: /*S
1692: MatNullSpace - Object that removes a null space from a vector, i.e.
1693: orthogonalizes the vector to a subsapce
1695: Level: advanced
1697: Users manual sections:
1698: . sec_singular
1700: .seealso: MatNullSpaceCreate()
1701: S*/
1702: typedef struct _p_MatNullSpace* MatNullSpace;
1704: PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1705: PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1706: PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1707: PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1708: PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1709: PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *);
1710: PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace);
1711: PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1712: PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1713: PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1714: PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *);
1715: PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1716: PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1717: PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
1719: PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1720: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1721: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1722: PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
1724: PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1725: PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1726: PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1728: PETSC_EXTERN PetscErrorCode MatComputeOperator(Mat,MatType,Mat*);
1729: PETSC_EXTERN PetscErrorCode MatComputeOperatorTranspose(Mat,MatType,Mat*);
1731: PETSC_DEPRECATED_FUNCTION("Use MatComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode MatComputeExplicitOperator(Mat A,Mat* B) { return MatComputeOperator(A,NULL,B); }
1732: PETSC_DEPRECATED_FUNCTION("Use MatComputeOperatorTranspose() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode MatComputeExplicitOperatorTranspose(Mat A,Mat* B) { return MatComputeOperatorTranspose(A,NULL,B); }
1734: PETSC_EXTERN PetscErrorCode MatCreateKAIJ(Mat,PetscInt,PetscInt,const PetscScalar[],const PetscScalar[],Mat*);
1735: PETSC_EXTERN PetscErrorCode MatKAIJGetAIJ(Mat,Mat*);
1736: PETSC_EXTERN PetscErrorCode MatKAIJGetS(Mat,PetscInt*,PetscInt*,PetscScalar**);
1737: PETSC_EXTERN PetscErrorCode MatKAIJGetSRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1738: PETSC_EXTERN PetscErrorCode MatKAIJRestoreS(Mat,PetscScalar**);
1739: PETSC_EXTERN PetscErrorCode MatKAIJRestoreSRead(Mat,const PetscScalar**);
1740: PETSC_EXTERN PetscErrorCode MatKAIJGetT(Mat,PetscInt*,PetscInt*,PetscScalar**);
1741: PETSC_EXTERN PetscErrorCode MatKAIJGetTRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1742: PETSC_EXTERN PetscErrorCode MatKAIJRestoreT(Mat,PetscScalar**);
1743: PETSC_EXTERN PetscErrorCode MatKAIJRestoreTRead(Mat,const PetscScalar**);
1744: PETSC_EXTERN PetscErrorCode MatKAIJSetAIJ(Mat,Mat);
1745: PETSC_EXTERN PetscErrorCode MatKAIJSetS(Mat,PetscInt,PetscInt,const PetscScalar[]);
1746: PETSC_EXTERN PetscErrorCode MatKAIJSetT(Mat,PetscInt,PetscInt,const PetscScalar[]);
1747: PETSC_EXTERN PetscErrorCode MatKAIJGetScaledIdentity(Mat,PetscBool*);
1749: PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
1751: PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1752: PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1753: PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1754: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1755: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1756: PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1757: PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1758: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1759: PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1760: PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1761: PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1762: PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1763: PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1765: /*S
1766: MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
1767: Jacobian vector products
1769: Notes:
1770: MATMFFD is a specific MatType which uses the MatMFFD data structure
1772: MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
1774: Level: developer
1776: .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
1777: S*/
1778: typedef struct _p_MatMFFD* MatMFFD;
1780: /*J
1781: MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1783: Level: beginner
1785: .seealso: MatMFFDSetType(), MatMFFDRegister()
1786: J*/
1787: typedef const char* MatMFFDType;
1788: #define MATMFFD_DS "ds"
1789: #define MATMFFD_WP "wp"
1791: PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1792: PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));
1794: PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1795: PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool);
1797: PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType);
1799: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1800: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
1802: #ifdef PETSC_HAVE_HARA
1803: PETSC_EXTERN_TYPEDEF typedef PetscScalar (*MatHaraKernel)(PetscInt,PetscReal[],PetscReal[],void*);
1804: PETSC_EXTERN PetscErrorCode MatCreateHaraFromKernel(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscReal[],MatHaraKernel,void*,PetscReal,PetscInt,PetscInt,Mat*);
1805: PETSC_EXTERN PetscErrorCode MatCreateHaraFromMat(Mat,PetscInt,const PetscReal[],PetscReal,PetscInt,PetscInt,PetscInt,PetscReal,Mat*);
1806: PETSC_EXTERN PetscErrorCode MatHaraSetSamplingMat(Mat,Mat,PetscInt,PetscReal);
1807: #endif
1809: /*
1810: PETSc interface to MUMPS
1811: */
1812: #ifdef PETSC_HAVE_MUMPS
1813: PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1814: PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*);
1815: PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
1816: PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*);
1818: PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1819: PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1820: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1821: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
1822: PETSC_EXTERN PetscErrorCode MatMumpsGetInverse(Mat,Mat);
1823: PETSC_EXTERN PetscErrorCode MatMumpsGetInverseTranspose(Mat,Mat);
1824: #endif
1826: /*
1827: PETSc interface to Mkl_Pardiso
1828: */
1829: #ifdef PETSC_HAVE_MKL_PARDISO
1830: PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1831: #endif
1833: /*
1834: PETSc interface to Mkl_CPardiso
1835: */
1836: #ifdef PETSC_HAVE_MKL_CPARDISO
1837: PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt);
1838: #endif
1840: /*
1841: PETSc interface to SUPERLU
1842: */
1843: #ifdef PETSC_HAVE_SUPERLU
1844: PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1845: #endif
1847: /*
1848: PETSc interface to SUPERLU_DIST
1849: */
1850: #ifdef PETSC_HAVE_SUPERLU_DIST
1851: PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*);
1852: #endif
1854: /*
1855: PETSc interface to STRUMPACK
1856: */
1857: #ifdef PETSC_HAVE_STRUMPACK
1858: /*E
1859: MatSTRUMPACKReordering - sparsity reducing ordering to be used in STRUMPACK
1861: Level: intermediate
1862: E*/
1863: typedef enum {MAT_STRUMPACK_NATURAL=0,
1864: MAT_STRUMPACK_METIS=1,
1865: MAT_STRUMPACK_PARMETIS=2,
1866: MAT_STRUMPACK_SCOTCH=3,
1867: MAT_STRUMPACK_PTSCOTCH=4,
1868: MAT_STRUMPACK_RCM=5} MatSTRUMPACKReordering;
1869: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetReordering(Mat,MatSTRUMPACKReordering);
1870: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool);
1871: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat,PetscReal);
1872: PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSRelTol() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat mat,PetscReal rtol) {return MatSTRUMPACKSetHSSRelTol(mat,rtol);}
1873: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat,PetscReal);
1874: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat,PetscInt);
1875: PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSMinSepSize() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat mat,PetscInt hssminsize) {return MatSTRUMPACKSetHSSMinSepSize(mat,hssminsize);}
1876: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat,PetscInt);
1877: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat,PetscInt);
1878: #endif
1881: PETSC_EXTERN PetscErrorCode MatBindToCPU(Mat,PetscBool);
1882: PETSC_DEPRECATED_FUNCTION("Use MatBindToCPU (since v3.13)") PETSC_STATIC_INLINE PetscErrorCode MatPinToCPU(Mat A,PetscBool flg) {return MatBindToCPU(A,flg);}
1884: #ifdef PETSC_HAVE_CUDA
1885: /*E
1886: MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
1887: matrices.
1889: Not Collective
1891: + MAT_CUSPARSE_CSR - Compressed Sparse Row
1892: . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later).
1893: - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).
1895: Level: intermediate
1897: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1899: .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1900: E*/
1902: typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;
1904: /* these will be strings associated with enumerated type defined above */
1905: PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
1907: /*E
1908: MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
1909: matrices whose operation should use a particular storage format.
1911: Not Collective
1913: + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
1914: . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
1915: . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
1916: - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices
1918: Level: intermediate
1920: .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1921: E*/
1922: typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;
1924: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1925: PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1926: PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
1927: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSESetGenerateTranspose(Mat,PetscBool);
1929: typedef struct _p_SplitCSRMat PetscSplitCSRDataStructure;
1930: PETSC_EXTERN PetscErrorCode MatCUSPARSEGetDeviceMatWrite(Mat,PetscSplitCSRDataStructure**);
1932: PETSC_EXTERN PetscErrorCode MatCreateDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
1933: PETSC_EXTERN PetscErrorCode MatCreateSeqDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
1934: PETSC_EXTERN PetscErrorCode MatMPIDenseCUDASetPreallocation(Mat,PetscScalar[]);
1935: PETSC_EXTERN PetscErrorCode MatSeqDenseCUDASetPreallocation(Mat,PetscScalar[]);
1936: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayWrite(Mat,PetscScalar**);
1937: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayRead(Mat,const PetscScalar**);
1938: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArray(Mat,PetscScalar**);
1939: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayWrite(Mat,PetscScalar**);
1940: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayRead(Mat,const PetscScalar**);
1941: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArray(Mat,PetscScalar**);
1942: PETSC_EXTERN PetscErrorCode MatDenseCUDAPlaceArray(Mat,const PetscScalar*);
1943: PETSC_EXTERN PetscErrorCode MatDenseCUDAReplaceArray(Mat,const PetscScalar*);
1944: PETSC_EXTERN PetscErrorCode MatDenseCUDAResetArray(Mat);
1946: #endif
1948: #if defined(PETSC_HAVE_VIENNACL)
1949: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1950: PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1951: #endif
1953: /*
1954: PETSc interface to FFTW
1955: */
1956: #if defined(PETSC_HAVE_FFTW)
1957: PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1958: PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
1959: PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*);
1960: #endif
1962: #if defined(PETSC_HAVE_SCALAPACK)
1963: PETSC_EXTERN PetscErrorCode MatCreateScaLAPACK(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1964: PETSC_EXTERN PetscErrorCode MatScaLAPACKSetBlockSizes(Mat,PetscInt,PetscInt);
1965: PETSC_EXTERN PetscErrorCode MatScaLAPACKGetBlockSizes(Mat,PetscInt*,PetscInt*);
1966: #endif
1968: PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1969: PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1970: PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1971: PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1972: PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1973: PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
1974: PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1975: PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1976: PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1978: PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
1979: PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);
1981: PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**);
1983: PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat);
1985: PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*);
1986: PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*);
1987: #endif