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