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