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