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