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