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