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