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