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