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