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