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