xref: /petsc/include/petscmat.h (revision 5bd1e5768a3b141382a3ebb2e780cd2d3b3bfbf0)
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_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
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 MatGetType(Mat,MatType*);
415 
416 PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
417 PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
418 PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
419 PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
420 PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
421 PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
422 PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
423 PETSC_EXTERN PetscErrorCode MatSeqAIJGetArrayRead(Mat,const PetscScalar *[]);
424 PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
425 PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArrayRead(Mat,const PetscScalar *[]);
426 PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*);
427 PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
428 PETSC_EXTERN PetscErrorCode MatSeqAIJSetType(Mat,MatType);
429 PETSC_EXTERN PetscErrorCode MatSeqAIJRegister(const char[],PetscErrorCode (*)(Mat,MatType,MatReuse,Mat *));
430 PETSC_EXTERN PetscFunctionList MatSeqAIJList;
431 PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat,PetscScalar *[]);
432 PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat,PetscScalar *[]);
433 PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
434 PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
435 PETSC_EXTERN PetscErrorCode MatDensePlaceArray(Mat,const PetscScalar[]);
436 PETSC_EXTERN PetscErrorCode MatDenseResetArray(Mat);
437 PETSC_EXTERN PetscErrorCode MatDenseGetArrayRead(Mat,const PetscScalar *[]);
438 PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayRead(Mat,const PetscScalar *[]);
439 PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
440 PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
441 PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
442 PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
443 PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat);
444 PETSC_EXTERN PetscErrorCode MatSetVariableBlockSizes(Mat,PetscInt,PetscInt*);
445 PETSC_EXTERN PetscErrorCode MatGetVariableBlockSizes(Mat,PetscInt*,const PetscInt**);
446 
447 PETSC_EXTERN PetscErrorCode MatDenseGetColumn(Mat,PetscInt,PetscScalar *[]);
448 PETSC_EXTERN PetscErrorCode MatDenseRestoreColumn(Mat,PetscScalar *[]);
449 
450 PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
451 PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
452 PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
453 PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
454 PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
455 PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
456 PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
457 PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
458 PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
459 PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
460 PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
461 PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
462 PETSC_EXTERN PetscErrorCode MatMatSolveTranspose(Mat,Mat,Mat);
463 PETSC_EXTERN PetscErrorCode MatMatTransposeSolve(Mat,Mat,Mat);
464 PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);
465 
466 /*E
467     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
468   its numerical values copied over or just its nonzero structure.
469 
470     Level: beginner
471 
472    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
473 
474 $   MAT_DO_NOT_COPY_VALUES    - Create a matrix using the same nonzero pattern as the original matrix,
475 $                               with zeros for the numerical values.
476 $   MAT_COPY_VALUES           - Create a matrix with the same nonzero pattern as the original matrix
477 $                               and with the same numerical values.
478 $   MAT_SHARE_NONZERO_PATTERN - Create a matrix that shares the nonzero structure with the previous matrix
479 $                               and does not copy it, using zeros for the numerical values. The parent and
480 $                               child matrices will share their index (i and j) arrays, and you cannot
481 $                               insert new nonzero entries into either matrix.
482 
483 Notes:
484     Many matrix types (including SeqAIJ) do not support the MAT_SHARE_NONZERO_PATTERN optimization; in
485 this case the behavior is as if MAT_DO_NOT_COPY_VALUES has been specified.
486 
487 .seealso: MatDuplicate()
488 E*/
489 typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
490 
491 PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
492 PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);
493 
494 
495 PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
496 PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
497 PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
498 PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
499 PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
500 PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
501 PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
502 PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
503 PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
504 
505 PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
506 PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
507 PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
508 PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
509 
510 /*S
511      MatInfo - Context of matrix information, used with MatGetInfo()
512 
513    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
514 
515    Level: intermediate
516 
517 .seealso:  MatGetInfo(), MatInfoType
518 S*/
519 typedef struct {
520   PetscLogDouble block_size;                         /* block size */
521   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
522   PetscLogDouble memory;                             /* memory allocated */
523   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
524   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
525   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
526   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
527 } MatInfo;
528 
529 /*E
530     MatInfoType - Indicates if you want information about the local part of the matrix,
531      the entire parallel matrix or the maximum over all the local parts.
532 
533     Level: beginner
534 
535    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
536 
537 .seealso: MatGetInfo(), MatInfo
538 E*/
539 typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
540 PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
541 PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
542 PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
543 PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
544 PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
545 PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
546 PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
547 PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
548 PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
549 PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*);
550 PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
551 PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);
552 
553 PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*);
554 PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*);
555 PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*);
556 PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*);
557 PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*);
558 PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
559 PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
560 PETSC_EXTERN PetscErrorCode MatMatTransposeMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
561 PETSC_EXTERN PetscErrorCode MatIsLinear(Mat,PetscInt,PetscBool*);
562 
563 PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*);
564 PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*);
565 PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
566 PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
567 PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
568 PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
569 PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
570 PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
571 PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
572 
573 PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
574 PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
575 PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
576 PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
577 PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
578 PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
579 PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);
580 
581 PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
582 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);}
583 PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
584 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);}
585 PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
586 PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]);
587 PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *);
588 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);}
589 PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
590 PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
591 PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
592 PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);
593 
594 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
595 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
596 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
597 PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
598 PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
599 PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
600 PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
601 
602 PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
603 PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov);
604 PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool);
605 
606 PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
607 PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
608 PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat);
609 
610 PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
611 PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
612 
613 PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
614 PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
615 PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat);
616 PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
617 PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
618 PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat);
619 
620 PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
621 PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
622 PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat);
623 PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
624 PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
625 PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat);
626 
627 PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
628 PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);
629 
630 PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
631 PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);
632 
633 PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
634 PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
635 PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*);
636 PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
637 PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
638 PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
639 PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
640 PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
641 PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
642 
643 PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
644 PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
645 
646 PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
647 PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
648 PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
649 PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*);
650 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);}
651 PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
652 PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
653 PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
654 PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
655 PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
656 
657 /*MC
658    MatSetValue - Set a single entry into a matrix.
659 
660    Not collective
661 
662    Synopsis:
663      #include <petscmat.h>
664      PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode)
665 
666    Input Parameters:
667 +  m - the matrix
668 .  row - the row location of the entry
669 .  col - the column location of the entry
670 .  value - the value to insert
671 -  mode - either INSERT_VALUES or ADD_VALUES
672 
673    Notes:
674    For efficiency one should use MatSetValues() and set several or many
675    values simultaneously if possible.
676 
677    Level: beginner
678 
679 .seealso: MatSetValues(), MatSetValueLocal()
680 M*/
681 PETSC_STATIC_INLINE PetscErrorCode MatSetValue(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValues(v,1,&i,1,&j,&va,mode);}
682 
683 PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
684 
685 PETSC_STATIC_INLINE PetscErrorCode MatSetValueLocal(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValuesLocal(v,1,&i,1,&j,&va,mode);}
686 
687 /*MC
688    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
689        row in a matrix providing the data that one can use to correctly preallocate the matrix.
690 
691    Synopsis:
692    #include <petscmat.h>
693    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
694 
695    Collective
696 
697    Input Parameters:
698 +  comm - the communicator that will share the eventually allocated matrix
699 .  nrows - the number of LOCAL rows in the matrix
700 -  ncols - the number of LOCAL columns in the matrix
701 
702    Output Parameters:
703 +  dnz - the array that will be passed to the matrix preallocation routines
704 -  onz - the other array passed to the matrix preallocation routines
705 
706    Level: intermediate
707 
708    Notes:
709     See Users-Manual: ch_performance for more details.
710 
711    Do not malloc or free dnz and onz, that is handled internally by these routines
712 
713    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
714 
715 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
716           MatPreallocateSymmetricSetLocalBlock()
717 M*/
718 #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
719 do { \
720   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ncols = (ncols),__rstart,__start,__end; \
721   _4_ierr = PetscCalloc2((size_t)__nrows,&dnz,(size_t)__nrows,&onz);CHKERRQ(_4_ierr); \
722   __start = 0; __end = __start;                                         \
723   _4_ierr = MPI_Scan(&__ncols,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ncols;\
724   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; \
725   do { } while(0)
726 
727 /*MC
728    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
729        inserted using a local number of the rows and columns
730 
731    Synopsis:
732    #include <petscmat.h>
733    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
734 
735    Not Collective
736 
737    Input Parameters:
738 +  map - the row mapping from local numbering to global numbering
739 .  nrows - the number of rows indicated
740 .  rows - the indices of the rows
741 .  cmap - the column mapping from local to global numbering
742 .  ncols - the number of columns in the matrix
743 .  cols - the columns indicated
744 .  dnz - the array that will be passed to the matrix preallocation routines
745 -  onz - the other array passed to the matrix preallocation routines
746 
747    Level: intermediate
748 
749    Notes:
750     See Users-Manual: ch_performance for more details.
751 
752    Do not malloc or free dnz and onz, that is handled internally by these routines
753 
754 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
755           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups()
756 M*/
757 #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
758 do {\
759   PetscInt __l;\
760   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
761   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
762   for (__l=0;__l<nrows;__l++) {\
763     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
764   }\
765  } while (0)
766 
767 /*MC
768    MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be
769        inserted using a local number of the rows and columns. This version removes any duplicate columns in cols
770 
771    Synopsis:
772    #include <petscmat.h>
773    PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
774 
775    Not Collective
776 
777    Input Parameters:
778 +  map - the row mapping from local numbering to global numbering
779 .  nrows - the number of rows indicated
780 .  rows - the indices of the rows (these values are mapped to the global values)
781 .  cmap - the column mapping from local to global numbering
782 .  ncols - the number of columns in the matrix   (this value will be changed if duplicate columns are found)
783 .  cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed)
784 .  dnz - the array that will be passed to the matrix preallocation routines
785 -  onz - the other array passed to the matrix preallocation routines
786 
787    Level: intermediate
788 
789    Notes:
790     See Users-Manual: ch_performance for more details.
791 
792    Do not malloc or free dnz and onz, that is handled internally by these routines
793 
794 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
795           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
796 M*/
797 #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
798 do {\
799   PetscInt __l;\
800   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
801   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
802   _4_ierr = PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\
803   for (__l=0;__l<nrows;__l++) {\
804     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
805   }\
806  } while (0)
807 
808 /*MC
809    MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
810        inserted using a local number of the rows and columns
811 
812    Synopsis:
813    #include <petscmat.h>
814    PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
815 
816    Not Collective
817 
818    Input Parameters:
819 +  map - the row mapping from local numbering to global numbering
820 .  nrows - the number of rows indicated
821 .  rows - the indices of the rows
822 .  cmap - the column mapping from local to global numbering
823 .  ncols - the number of columns in the matrix
824 .  cols - the columns indicated
825 .  dnz - the array that will be passed to the matrix preallocation routines
826 -  onz - the other array passed to the matrix preallocation routines
827 
828    Level: intermediate
829 
830    Notes:
831     See Users-Manual: ch_performance for more details.
832 
833    Do not malloc or free dnz and onz, that is handled internally by these routines
834 
835 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
836           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
837 M*/
838 #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
839 do {\
840   PetscInt __l;\
841   _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
842   _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
843   for (__l=0;__l<nrows;__l++) {\
844     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
845   }\
846 } while (0)
847 
848 /*MC
849    MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
850        inserted using a local number of the rows and columns
851 
852    Synopsis:
853    #include <petscmat.h>
854    PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
855 
856    Not Collective
857 
858    Input Parameters:
859 +  map - the mapping between local numbering and global numbering
860 .  nrows - the number of rows indicated
861 .  rows - the indices of the rows
862 .  ncols - the number of columns in the matrix
863 .  cols - the columns indicated
864 .  dnz - the array that will be passed to the matrix preallocation routines
865 -  onz - the other array passed to the matrix preallocation routines
866 
867    Level: intermediate
868 
869    Notes:
870     See Users-Manual: ch_performance for more details.
871 
872    Do not malloc or free dnz and onz that is handled internally by these routines
873 
874 .seealso: MatPreallocateFinalize(), MatPreallocateSet()
875           MatPreallocateInitialize(),  MatPreallocateSetLocal()
876 M*/
877 #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\
878 do {\
879   PetscInt __l;\
880   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
881   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
882   for (__l=0;__l<nrows;__l++) {\
883     _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
884   }\
885 } while (0)
886 
887 /*MC
888    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
889        inserted using a local number of the rows and columns
890 
891    Synopsis:
892    #include <petscmat.h>
893    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
894 
895    Not Collective
896 
897    Input Parameters:
898 +  row - the row
899 .  ncols - the number of columns in the matrix
900 -  cols - the columns indicated
901 
902    Output Parameters:
903 +  dnz - the array that will be passed to the matrix preallocation routines
904 -  onz - the other array passed to the matrix preallocation routines
905 
906    Level: intermediate
907 
908    Notes:
909     See Users-Manual: ch_performance for more details.
910 
911    Do not malloc or free dnz and onz that is handled internally by these routines
912 
913    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
914 
915 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
916           MatPreallocateInitialize(), MatPreallocateSetLocal()
917 M*/
918 #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
919 do { PetscInt __i; \
920   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);\
921   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);\
922   for (__i=0; __i<nc; __i++) {\
923     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
924     else if (dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\
925   }\
926 } while (0)
927 
928 /*MC
929    MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
930        inserted using a local number of the rows and columns
931 
932    Synopsis:
933    #include <petscmat.h>
934    PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
935 
936    Not Collective
937 
938    Input Parameters:
939 +  nrows - the number of rows indicated
940 .  rows - the indices of the rows
941 .  ncols - the number of columns in the matrix
942 .  cols - the columns indicated
943 .  dnz - the array that will be passed to the matrix preallocation routines
944 -  onz - the other array passed to the matrix preallocation routines
945 
946    Level: intermediate
947 
948    Notes:
949     See Users-Manual: ch_performance for more details.
950 
951    Do not malloc or free dnz and onz that is handled internally by these routines
952 
953    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
954 
955 .seealso: MatPreallocateFinalize(), MatPreallocateSet(),  MatPreallocateInitialize(),
956           MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
957 M*/
958 #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\
959 do { PetscInt __i; \
960   for (__i=0; __i<nc; __i++) {\
961     if (cols[__i] >= __end) onz[row - __rstart]++; \
962     else if (cols[__i] >= row && dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\
963   }\
964 } while (0)
965 
966 /*MC
967    MatPreallocateLocation -  An alternative to MatPreallocateSet() that puts the nonzero locations into the matrix if it exists
968 
969    Synopsis:
970    #include <petscmat.h>
971    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
972 
973    Not Collective
974 
975    Input Parameters:
976 +  A - matrix
977 .  row - row where values exist (must be local to this process)
978 .  ncols - number of columns
979 .  cols - columns with nonzeros
980 .  dnz - the array that will be passed to the matrix preallocation routines
981 -  onz - the other array passed to the matrix preallocation routines
982 
983    Level: intermediate
984 
985    Notes:
986     See Users-Manual: ch_performance for more details.
987 
988    Do not malloc or free dnz and onz that is handled internally by these routines
989 
990    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
991 
992 .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
993           MatPreallocateSymmetricSetLocalBlock()
994 M*/
995 #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)
996 
997 
998 /*MC
999    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
1000        row in a matrix providing the data that one can use to correctly preallocate the matrix.
1001 
1002    Synopsis:
1003    #include <petscmat.h>
1004    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
1005 
1006    Collective
1007 
1008    Input Parameters:
1009 +  dnz - the array that was be passed to the matrix preallocation routines
1010 -  onz - the other array passed to the matrix preallocation routines
1011 
1012    Level: intermediate
1013 
1014    Notes:
1015     See Users-Manual: ch_performance for more details.
1016 
1017    Do not malloc or free dnz and onz that is handled internally by these routines
1018 
1019    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
1020 
1021 .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
1022           MatPreallocateSymmetricSetLocalBlock()
1023 M*/
1024 #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} while(0)
1025 
1026 /* Routines unique to particular data structures */
1027 PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
1028 
1029 PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
1030 PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
1031 
1032 PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
1033 PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
1034 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1035 PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1036 PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1037 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
1038 
1039 #define MAT_SKIP_ALLOCATION -4
1040 
1041 PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1042 PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1043 PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
1044 
1045 PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1046 PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1047 PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1048 PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
1049 PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1050 PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
1051 PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1052 PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
1053 PETSC_EXTERN PetscErrorCode MatMPIAdjToSeq(Mat,Mat*);
1054 PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
1055 PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
1056 PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1057 PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1058 PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
1059 
1060 PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt);
1061 PETSC_EXTERN PetscErrorCode MatDenseGetLDA(Mat,PetscInt*);
1062 PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
1063 
1064 PETSC_EXTERN PetscErrorCode MatBlockMatSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1065 
1066 PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
1067 PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
1068 
1069 PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
1070 
1071 PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
1072 PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*);
1073 /*
1074   These routines are not usually accessed directly, rather solving is
1075   done through the KSP and PC interfaces.
1076 */
1077 
1078 /*J
1079     MatOrderingType - String with the name of a PETSc matrix ordering
1080 
1081    Level: beginner
1082 
1083 .seealso: MatGetOrdering()
1084 J*/
1085 typedef const char* MatOrderingType;
1086 #define MATORDERINGNATURAL     "natural"
1087 #define MATORDERINGND          "nd"
1088 #define MATORDERING1WD         "1wd"
1089 #define MATORDERINGRCM         "rcm"
1090 #define MATORDERINGQMD         "qmd"
1091 #define MATORDERINGROWLENGTH   "rowlength"
1092 #define MATORDERINGWBM         "wbm"
1093 #define MATORDERINGSPECTRAL    "spectral"
1094 #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
1095 
1096 PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
1097 PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
1098 PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
1099 PETSC_EXTERN PetscFunctionList MatOrderingList;
1100 
1101 PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1102 PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);
1103 
1104 /*S
1105     MatFactorShiftType - Numeric Shift.
1106 
1107    Level: beginner
1108 
1109 S*/
1110 typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1111 PETSC_EXTERN const char *const MatFactorShiftTypes[];
1112 PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];
1113 
1114 /*S
1115     MatFactorError - indicates what type of error in matrix factor
1116 
1117     Level: beginner
1118 
1119     Developer Notes:
1120     Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1121 
1122 .seealso: MatGetFactor()
1123 S*/
1124 typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError;
1125 
1126 PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*);
1127 PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat);
1128 PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*);
1129 
1130 /*S
1131    MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization
1132 
1133    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
1134 $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1135 
1136    Notes:
1137     These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1138 
1139       You can use MatFactorInfoInitialize() to set default values.
1140 
1141    Level: developer
1142 
1143 .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1144           MatFactorInfoInitialize()
1145 
1146 S*/
1147 typedef struct {
1148   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
1149   PetscReal     usedt;
1150   PetscReal     dt;             /* drop tolerance */
1151   PetscReal     dtcol;          /* tolerance for pivoting */
1152   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
1153   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1154   PetscReal     levels;         /* ICC/ILU(levels) */
1155   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1156                                    factorization may be faster if do not pivot */
1157   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1158   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1159   PetscReal     shiftamount;     /* how large the shift is */
1160 } MatFactorInfo;
1161 
1162 PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
1163 PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
1164 PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1165 PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
1166 PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
1167 PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
1168 PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1169 PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1170 PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1171 PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
1172 PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1173 PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1174 PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1175 PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1176 PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1177 PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1178 PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1179 PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1180 PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
1181 PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
1182 
1183 typedef enum {MAT_FACTOR_SCHUR_UNFACTORED, MAT_FACTOR_SCHUR_FACTORED, MAT_FACTOR_SCHUR_INVERTED} MatFactorSchurStatus;
1184 PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS);
1185 PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1186 PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*,MatFactorSchurStatus);
1187 PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat);
1188 PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1189 PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec);
1190 PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec);
1191 PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat);
1192 
1193 /*E
1194     MatSORType - What type of (S)SOR to perform
1195 
1196     Level: beginner
1197 
1198    May be bitwise ORd together
1199 
1200    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1201 
1202    MatSORType may be bitwise ORd together, so do not change the numbers
1203 
1204 .seealso: MatSOR()
1205 E*/
1206 typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1207               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1208               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
1209               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1210 PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1211 
1212 /*
1213     These routines are for efficiently computing Jacobians via finite differences.
1214 */
1215 
1216 /*S
1217      MatColoring - Object for managing the coloring of matrices.
1218 
1219    Level: beginner
1220 
1221 .seealso:  MatFDColoringCreate() ISColoring MatFDColoring
1222 S*/
1223 typedef struct _p_MatColoring* MatColoring;
1224 /*J
1225     MatColoringType - String with the name of a PETSc matrix coloring
1226 
1227    Level: beginner
1228 
1229 .seealso: MatColoringSetType(), MatColoring
1230 J*/
1231 
1232 typedef const  char*           MatColoringType;
1233 #define MATCOLORINGJP      "jp"
1234 #define MATCOLORINGPOWER   "power"
1235 #define MATCOLORINGNATURAL "natural"
1236 #define MATCOLORINGSL      "sl"
1237 #define MATCOLORINGLF      "lf"
1238 #define MATCOLORINGID      "id"
1239 #define MATCOLORINGGREEDY  "greedy"
1240 
1241 /*E
1242    MatColoringWeightType - Type of weight scheme
1243 
1244     Not Collective
1245 
1246 +   MAT_COLORING_RANDOM  - Random weights
1247 .   MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
1248 -   MAT_COLORING_LF      - Last-first weighting.
1249 
1250     Level: intermediate
1251 
1252    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1253 E*/
1254 typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;
1255 
1256 PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1257 PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*);
1258 PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1259 PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1260 PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1261 PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1262 PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1263 PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1264 PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1265 PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1266 PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1267 PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1268 PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1269 PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType);
1270 PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*);
1271 PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm);
1272 PETSC_EXTERN PetscErrorCode MatColoringTest(MatColoring,ISColoring);
1273 PETSC_DEPRECATED_FUNCTION("Use MatColoringTest() (since version 3.10)") PETSC_STATIC_INLINE PetscErrorCode MatColoringTestValid(MatColoring matcoloring,ISColoring iscoloring) {return MatColoringTest(matcoloring,iscoloring);}
1274 PETSC_EXTERN PetscErrorCode MatISColoringTest(Mat,ISColoring);
1275 
1276 /*S
1277      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1278         and coloring
1279 
1280    Level: beginner
1281 
1282 .seealso:  MatFDColoringCreate()
1283 S*/
1284 typedef struct _p_MatFDColoring* MatFDColoring;
1285 
1286 PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1287 PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1288 PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1289 PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1290 PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1291 PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1292 PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1293 PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *);
1294 PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1295 PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,const PetscInt*[]);
1296 PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1297 PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1298 PETSC_EXTERN PetscErrorCode MatFDColoringSetValues(Mat,MatFDColoring,const PetscScalar*);
1299 
1300 /*S
1301      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1302 
1303    Level: beginner
1304 
1305 .seealso:  MatTransposeColoringCreate()
1306 S*/
1307 typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1308 
1309 PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1310 PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1311 PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1312 PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1313 
1314 /*
1315     These routines are for partitioning matrices: currently used only
1316   for adjacency matrix, MatCreateMPIAdj().
1317 */
1318 
1319 /*S
1320      MatPartitioning - Object for managing the partitioning of a matrix or graph
1321 
1322    Level: beginner
1323 
1324 .seealso:  MatPartitioningCreate(), MatPartitioningType
1325 S*/
1326 typedef struct _p_MatPartitioning* MatPartitioning;
1327 
1328 /*J
1329     MatPartitioningType - String with the name of a PETSc matrix partitioning
1330 
1331    Level: beginner
1332 dm
1333 .seealso: MatPartitioningCreate(), MatPartitioning
1334 J*/
1335 typedef const char* MatPartitioningType;
1336 #define MATPARTITIONINGCURRENT  "current"
1337 #define MATPARTITIONINGAVERAGE   "average"
1338 #define MATPARTITIONINGSQUARE   "square"
1339 #define MATPARTITIONINGPARMETIS "parmetis"
1340 #define MATPARTITIONINGCHACO    "chaco"
1341 #define MATPARTITIONINGPARTY    "party"
1342 #define MATPARTITIONINGPTSCOTCH "ptscotch"
1343 #define MATPARTITIONINGHIERARCH  "hierarch"
1344 
1345 
1346 PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1347 PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1348 PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1349 PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1350 PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1351 PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1352 PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1353 PETSC_EXTERN PetscErrorCode MatPartitioningImprove(MatPartitioning,IS*);
1354 PETSC_EXTERN PetscErrorCode MatPartitioningViewImbalance(MatPartitioning,IS);
1355 PETSC_EXTERN PetscErrorCode MatPartitioningApplyND(MatPartitioning,IS*);
1356 PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
1357 
1358 PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
1359 
1360 
1361 
1362 PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1363 PETSC_STATIC_INLINE PetscErrorCode MatPartitioningViewFromOptions(MatPartitioning A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);}
1364 
1365 PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
1366 PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1367 
1368 PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1369 PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1370 PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
1371 
1372 typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
1373 PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1374 typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
1375 PETSC_EXTERN const char *const MPChacoLocalTypes[];
1376 typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
1377 PETSC_EXTERN const char *const MPChacoEigenTypes[];
1378 
1379 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1380 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1381 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1382 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1383 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1384 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1385 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1386 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1387 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1388 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1389 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
1390 
1391 #define MP_PARTY_OPT "opt"
1392 #define MP_PARTY_LIN "lin"
1393 #define MP_PARTY_SCA "sca"
1394 #define MP_PARTY_RAN "ran"
1395 #define MP_PARTY_GBF "gbf"
1396 #define MP_PARTY_GCF "gcf"
1397 #define MP_PARTY_BUB "bub"
1398 #define MP_PARTY_DEF "def"
1399 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
1400 #define MP_PARTY_HELPFUL_SETS "hs"
1401 #define MP_PARTY_KERNIGHAN_LIN "kl"
1402 #define MP_PARTY_NONE "no"
1403 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1404 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1405 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1406 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
1407 
1408 typedef enum { MP_PTSCOTCH_DEFAULT,MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
1409 PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1410 
1411 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1412 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1413 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1414 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
1415 
1416 /*
1417  * hierarchical partitioning
1418  */
1419 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*);
1420 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*);
1421 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt);
1422 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);
1423 
1424 PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1425 PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1426 
1427 /*
1428     If you add entries here you must also add them to petsc/finclude/petscmat.h
1429 */
1430 typedef enum { MATOP_SET_VALUES=0,
1431                MATOP_GET_ROW=1,
1432                MATOP_RESTORE_ROW=2,
1433                MATOP_MULT=3,
1434                MATOP_MULT_ADD=4,
1435                MATOP_MULT_TRANSPOSE=5,
1436                MATOP_MULT_TRANSPOSE_ADD=6,
1437                MATOP_SOLVE=7,
1438                MATOP_SOLVE_ADD=8,
1439                MATOP_SOLVE_TRANSPOSE=9,
1440                MATOP_SOLVE_TRANSPOSE_ADD=10,
1441                MATOP_LUFACTOR=11,
1442                MATOP_CHOLESKYFACTOR=12,
1443                MATOP_SOR=13,
1444                MATOP_TRANSPOSE=14,
1445                MATOP_GETINFO=15,
1446                MATOP_EQUAL=16,
1447                MATOP_GET_DIAGONAL=17,
1448                MATOP_DIAGONAL_SCALE=18,
1449                MATOP_NORM=19,
1450                MATOP_ASSEMBLY_BEGIN=20,
1451                MATOP_ASSEMBLY_END=21,
1452                MATOP_SET_OPTION=22,
1453                MATOP_ZERO_ENTRIES=23,
1454                MATOP_ZERO_ROWS=24,
1455                MATOP_LUFACTOR_SYMBOLIC=25,
1456                MATOP_LUFACTOR_NUMERIC=26,
1457                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1458                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1459                MATOP_SETUP_PREALLOCATION=29,
1460                MATOP_ILUFACTOR_SYMBOLIC=30,
1461                MATOP_ICCFACTOR_SYMBOLIC=31,
1462                MATOP_GET_DIAGONAL_BLOCK=32,
1463                MATOP_FREE_INTER_STRUCT=33,
1464                MATOP_DUPLICATE=34,
1465                MATOP_FORWARD_SOLVE=35,
1466                MATOP_BACKWARD_SOLVE=36,
1467                MATOP_ILUFACTOR=37,
1468                MATOP_ICCFACTOR=38,
1469                MATOP_AXPY=39,
1470                MATOP_CREATE_SUBMATRICES=40,
1471                MATOP_INCREASE_OVERLAP=41,
1472                MATOP_GET_VALUES=42,
1473                MATOP_COPY=43,
1474                MATOP_GET_ROW_MAX=44,
1475                MATOP_SCALE=45,
1476                MATOP_SHIFT=46,
1477                MATOP_DIAGONAL_SET=47,
1478                MATOP_ZERO_ROWS_COLUMNS=48,
1479                MATOP_SET_RANDOM=49,
1480                MATOP_GET_ROW_IJ=50,
1481                MATOP_RESTORE_ROW_IJ=51,
1482                MATOP_GET_COLUMN_IJ=52,
1483                MATOP_RESTORE_COLUMN_IJ=53,
1484                MATOP_FDCOLORING_CREATE=54,
1485                MATOP_COLORING_PATCH=55,
1486                MATOP_SET_UNFACTORED=56,
1487                MATOP_PERMUTE=57,
1488                MATOP_SET_VALUES_BLOCKED=58,
1489                MATOP_CREATE_SUBMATRIX=59,
1490                MATOP_DESTROY=60,
1491                MATOP_VIEW=61,
1492                MATOP_CONVERT_FROM=62,
1493                MATOP_MATMAT_MULT=63,
1494                MATOP_MATMAT_MULT_SYMBOLIC=64,
1495                MATOP_MATMAT_MULT_NUMERIC=65,
1496                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1497                MATOP_SET_VALUES_LOCAL=67,
1498                MATOP_ZERO_ROWS_LOCAL=68,
1499                MATOP_GET_ROW_MAX_ABS=69,
1500                MATOP_GET_ROW_MIN_ABS=70,
1501                MATOP_CONVERT=71,
1502                MATOP_SET_COLORING=72,
1503                /* MATOP_PLACEHOLDER_73=73, */
1504                MATOP_SET_VALUES_ADIFOR=74,
1505                MATOP_FD_COLORING_APPLY=75,
1506                MATOP_SET_FROM_OPTIONS=76,
1507                MATOP_MULT_CONSTRAINED=77,
1508                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
1509                MATOP_FIND_ZERO_DIAGONALS=79,
1510                MATOP_MULT_MULTIPLE=80,
1511                MATOP_SOLVE_MULTIPLE=81,
1512                MATOP_GET_INERTIA=82,
1513                MATOP_LOAD=83,
1514                MATOP_IS_SYMMETRIC=84,
1515                MATOP_IS_HERMITIAN=85,
1516                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1517                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1518                MATOP_CREATE_VECS=88,
1519                MATOP_MAT_MULT=89,
1520                MATOP_MAT_MULT_SYMBOLIC=90,
1521                MATOP_MAT_MULT_NUMERIC=91,
1522                MATOP_PTAP=92,
1523                MATOP_PTAP_SYMBOLIC=93,
1524                MATOP_PTAP_NUMERIC=94,
1525                MATOP_MAT_TRANSPOSE_MULT=95,
1526                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1527                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
1528                /* MATOP_PLACEHOLDER_98=98, */
1529                /* MATOP_PLACEHOLDER_99=99, */
1530                /* MATOP_PLACEHOLDER_100=100, */
1531                /* MATOP_PLACEHOLDER_101=101, */
1532                MATOP_CONJUGATE=102,
1533                /* MATOP_PLACEHOLDER_103=103, */
1534                MATOP_SET_VALUES_ROW=104,
1535                MATOP_REAL_PART=105,
1536                MATOP_IMAGINARY_PART=106,
1537                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1538                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1539                MATOP_MAT_SOLVE=109,
1540                MATOP_MAT_SOLVE_TRANSPOSE=110,
1541                MATOP_GET_ROW_MIN=111,
1542                MATOP_GET_COLUMN_VECTOR=112,
1543                MATOP_MISSING_DIAGONAL=113,
1544                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
1545                MATOP_CREATE=115,
1546                MATOP_GET_GHOSTS=116,
1547                MATOP_GET_LOCAL_SUB_MATRIX=117,
1548                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1549                MATOP_MULT_DIAGONAL_BLOCK=119,
1550                MATOP_HERMITIAN_TRANSPOSE=120,
1551                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
1552                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
1553                MATOP_GET_MULTI_PROC_BLOCK=123,
1554                MATOP_FIND_NONZERO_ROWS=124,
1555                MATOP_GET_COLUMN_NORMS=125,
1556                MATOP_INVERT_BLOCK_DIAGONAL=126,
1557                /* MATOP_PLACEHOLDER_127=127, */
1558                MATOP_CREATE_SUB_MATRICES_MPI=128,
1559                MATOP_SET_VALUES_BATCH=129,
1560                MATOP_TRANSPOSE_MAT_MULT=130,
1561                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
1562                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
1563                MATOP_TRANSPOSE_COLORING_CREAT=133,
1564                MATOP_TRANS_COLORING_APPLY_SPT=134,
1565                MATOP_TRANS_COLORING_APPLY_DEN=135,
1566                MATOP_RART=136,
1567                MATOP_RART_SYMBOLIC=137,
1568                MATOP_RART_NUMERIC=138,
1569                MATOP_SET_BLOCK_SIZES=139,
1570                MATOP_AYPX=140,
1571                MATOP_RESIDUAL=141,
1572                MATOP_FDCOLORING_SETUP=142,
1573                MATOP_MPICONCATENATESEQ=144,
1574                MATOP_DESTROYSUBMATRICES=145
1575              } MatOperation;
1576 PETSC_EXTERN PetscErrorCode MatSetOperation(Mat,MatOperation,void(*)(void));
1577 PETSC_EXTERN PetscErrorCode MatGetOperation(Mat,MatOperation,void(**)(void));
1578 PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *);
1579 PETSC_EXTERN PetscErrorCode MatHasCongruentLayouts(Mat,PetscBool *);
1580 PETSC_EXTERN PetscErrorCode MatFreeIntermediateDataStructures(Mat);
1581 
1582 PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1583 PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1584 PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1585 PETSC_EXTERN PetscErrorCode MatShellTestMult(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*);
1586 PETSC_EXTERN PetscErrorCode MatShellTestMultTranspose(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*);
1587 PETSC_EXTERN PetscErrorCode MatShellSetManageScalingShifts(Mat);
1588 
1589 /*
1590    Codes for matrices stored on disk. By default they are
1591    stored in a universal format. By changing the format with
1592    PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
1593    be stored in a way natural for the matrix, for example dense matrices
1594    would be stored as dense. Matrices stored this way may only be
1595    read into matrices of the same type.
1596 */
1597 #define MATRIX_BINARY_FORMAT_DENSE -1
1598 
1599 PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1600 
1601 PETSC_EXTERN PetscErrorCode MatISSetLocalMatType(Mat,MatType);
1602 PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1603 PETSC_EXTERN PetscErrorCode MatISStoreL2L(Mat,PetscBool);
1604 PETSC_EXTERN PetscErrorCode MatISFixLocalEmpty(Mat,PetscBool);
1605 PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1606 PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat,Mat*);
1607 PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
1608 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use the MatConvert() interface (since version 3.10)") PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*);
1609 
1610 /*S
1611      MatNullSpace - Object that removes a null space from a vector, i.e.
1612          orthogonalizes the vector to a subsapce
1613 
1614    Level: advanced
1615 
1616   Users manual sections:
1617 .   sec_singular
1618 
1619 .seealso:  MatNullSpaceCreate()
1620 S*/
1621 typedef struct _p_MatNullSpace* MatNullSpace;
1622 
1623 PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1624 PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1625 PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1626 PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1627 PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1628 PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *);
1629 PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace);
1630 PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1631 PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1632 PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1633 PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1634 PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1635 PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1636 PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
1637 
1638 PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1639 PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1640 PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1641 PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
1642 
1643 PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1644 PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1645 PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1646 
1647 PETSC_EXTERN PetscErrorCode MatComputeOperator(Mat,MatType,Mat*);
1648 PETSC_EXTERN PetscErrorCode MatComputeOperatorTranspose(Mat,MatType,Mat*);
1649 
1650 PETSC_DEPRECATED_FUNCTION("Use MatComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode MatComputeExplicitOperator(Mat A,Mat* B) { return MatComputeOperator(A,NULL,B); }
1651 PETSC_DEPRECATED_FUNCTION("Use MatComputeOperatorTranspose() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode MatComputeExplicitOperatorTranspose(Mat A,Mat* B) { return MatComputeOperatorTranspose(A,NULL,B); }
1652 
1653 PETSC_EXTERN PetscErrorCode MatCreateKAIJ(Mat,PetscInt,PetscInt,const PetscScalar[],const PetscScalar[],Mat*);
1654 PETSC_EXTERN PetscErrorCode MatKAIJGetAIJ(Mat,Mat*);
1655 PETSC_EXTERN PetscErrorCode MatKAIJGetS(Mat,PetscInt*,PetscInt*,PetscScalar**);
1656 PETSC_EXTERN PetscErrorCode MatKAIJGetSRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1657 PETSC_EXTERN PetscErrorCode MatKAIJRestoreS(Mat,PetscScalar**);
1658 PETSC_EXTERN PetscErrorCode MatKAIJRestoreSRead(Mat,const PetscScalar**);
1659 PETSC_EXTERN PetscErrorCode MatKAIJGetT(Mat,PetscInt*,PetscInt*,PetscScalar**);
1660 PETSC_EXTERN PetscErrorCode MatKAIJGetTRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1661 PETSC_EXTERN PetscErrorCode MatKAIJRestoreT(Mat,PetscScalar**);
1662 PETSC_EXTERN PetscErrorCode MatKAIJRestoreTRead(Mat,const PetscScalar**);
1663 PETSC_EXTERN PetscErrorCode MatKAIJSetAIJ(Mat,Mat);
1664 PETSC_EXTERN PetscErrorCode MatKAIJSetS(Mat,PetscInt,PetscInt,const PetscScalar []);
1665 PETSC_EXTERN PetscErrorCode MatKAIJSetT(Mat,PetscInt,PetscInt,const PetscScalar []);
1666 
1667 PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
1668 
1669 PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1670 PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1671 PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1672 PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1673 PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1674 PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1675 PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1676 PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1677 PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1678 PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1679 PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1680 PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1681 PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1682 
1683 /*S
1684     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
1685               Jacobian vector products
1686 
1687     Notes:
1688     MATMFFD is a specific MatType which uses the MatMFFD data structure
1689 
1690            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
1691 
1692     Level: developer
1693 
1694 .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
1695 S*/
1696 typedef struct _p_MatMFFD* MatMFFD;
1697 
1698 /*J
1699     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1700 
1701    Level: beginner
1702 
1703 .seealso: MatMFFDSetType(), MatMFFDRegister()
1704 J*/
1705 typedef const char* MatMFFDType;
1706 #define MATMFFD_DS  "ds"
1707 #define MATMFFD_WP  "wp"
1708 
1709 PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1710 PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));
1711 
1712 PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1713 PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool );
1714 
1715 PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType);
1716 
1717 PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1718 PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
1719 
1720 /*
1721    PETSc interface to MUMPS
1722 */
1723 #ifdef PETSC_HAVE_MUMPS
1724 PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1725 PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*);
1726 PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
1727 PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*);
1728 
1729 PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1730 PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1731 PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1732 PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
1733 PETSC_EXTERN PetscErrorCode MatMumpsGetInverse(Mat,Mat);
1734 PETSC_EXTERN PetscErrorCode MatMumpsGetInverseTranspose(Mat,Mat);
1735 #endif
1736 
1737 /*
1738    PETSc interface to Mkl_Pardiso
1739 */
1740 #ifdef PETSC_HAVE_MKL_PARDISO
1741 PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1742 #endif
1743 
1744 /*
1745    PETSc interface to Mkl_CPardiso
1746 */
1747 #ifdef PETSC_HAVE_MKL_CPARDISO
1748 PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt);
1749 #endif
1750 
1751 /*
1752    PETSc interface to SUPERLU
1753 */
1754 #ifdef PETSC_HAVE_SUPERLU
1755 PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1756 #endif
1757 
1758 /*
1759    PETSc interface to SUPERLU_DIST
1760 */
1761 #ifdef PETSC_HAVE_SUPERLU_DIST
1762 PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*);
1763 #endif
1764 
1765 /*
1766    PETSc interface to STRUMPACK
1767 */
1768 #ifdef PETSC_HAVE_STRUMPACK
1769 /*E
1770     MatSTRUMPACKReordering - sparsity reducing ordering to be used in STRUMPACK
1771 
1772     Level: intermediate
1773 E*/
1774 typedef enum {MAT_STRUMPACK_NATURAL=0,
1775               MAT_STRUMPACK_METIS=1,
1776               MAT_STRUMPACK_PARMETIS=2,
1777               MAT_STRUMPACK_SCOTCH=3,
1778               MAT_STRUMPACK_PTSCOTCH=4,
1779               MAT_STRUMPACK_RCM=5} MatSTRUMPACKReordering;
1780 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetReordering(Mat,MatSTRUMPACKReordering);
1781 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool);
1782 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat,PetscReal);
1783 PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSRelTol() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat mat,PetscReal rtol) {return MatSTRUMPACKSetHSSRelTol(mat,rtol);}
1784 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat,PetscReal);
1785 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat,PetscInt);
1786 PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSMinSepSize() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat mat,PetscInt hssminsize) {return MatSTRUMPACKSetHSSMinSepSize(mat,hssminsize);}
1787 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat,PetscInt);
1788 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat,PetscInt);
1789 #endif
1790 
1791 
1792 PETSC_EXTERN PetscErrorCode MatPinToCPU(Mat,PetscBool);
1793 
1794 #ifdef PETSC_HAVE_CUDA
1795 /*E
1796     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
1797     matrices.
1798 
1799     Not Collective
1800 
1801 +   MAT_CUSPARSE_CSR - Compressed Sparse Row
1802 .   MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later).
1803 -   MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).
1804 
1805     Level: intermediate
1806 
1807    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1808 
1809 .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1810 E*/
1811 
1812 typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;
1813 
1814 /* these will be strings associated with enumerated type defined above */
1815 PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
1816 
1817 /*E
1818     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
1819     matrices whose operation should use a particular storage format.
1820 
1821     Not Collective
1822 
1823 +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
1824 .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
1825 .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
1826 -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices
1827 
1828     Level: intermediate
1829 
1830 .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1831 E*/
1832 typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;
1833 
1834 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1835 PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1836 PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
1837 #endif
1838 
1839 #if defined(PETSC_HAVE_VIENNACL)
1840 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1841 PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1842 #endif
1843 
1844 /*
1845    PETSc interface to FFTW
1846 */
1847 #if defined(PETSC_HAVE_FFTW)
1848 PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1849 PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
1850 PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*);
1851 #endif
1852 
1853 PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1854 PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1855 PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1856 PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1857 PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1858 PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
1859 PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1860 PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1861 PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1862 
1863 PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
1864 PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);
1865 
1866 PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**);
1867 
1868 PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat);
1869 
1870 PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*);
1871 PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*);
1872 
1873 #endif
1874