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