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