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