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