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