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