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