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