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