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