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