xref: /petsc/include/petscmat.h (revision d7cc930e14e615e9907267aaa472dd0ccceeab82)
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 MatSetLayouts(Mat,PetscLayout,PetscLayout);
685 PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
686 PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
687 PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
688 PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
689 PETSC_EXTERN PetscErrorCode MatGetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
690 PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
691 PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
692 
693 PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
694 PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
695 
696 PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
697 PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
698 PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
699 PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*);
700 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);}
701 PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
702 PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
703 PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
704 PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
705 PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
706 
707 /*MC
708    MatSetValue - Set a single entry into a matrix.
709 
710    Not collective
711 
712    Synopsis:
713      #include <petscmat.h>
714      PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode)
715 
716    Input Parameters:
717 +  m - the matrix
718 .  row - the row location of the entry
719 .  col - the column location of the entry
720 .  value - the value to insert
721 -  mode - either INSERT_VALUES or ADD_VALUES
722 
723    Notes:
724    For efficiency one should use MatSetValues() and set several or many
725    values simultaneously if possible.
726 
727    Level: beginner
728 
729 .seealso: MatSetValues(), MatSetValueLocal()
730 M*/
731 PETSC_STATIC_INLINE PetscErrorCode MatSetValue(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValues(v,1,&i,1,&j,&va,mode);}
732 
733 PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
734 
735 PETSC_STATIC_INLINE PetscErrorCode MatSetValueLocal(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValuesLocal(v,1,&i,1,&j,&va,mode);}
736 
737 /*MC
738    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
739        row in a matrix providing the data that one can use to correctly preallocate the matrix.
740 
741    Synopsis:
742    #include <petscmat.h>
743    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
744 
745    Collective
746 
747    Input Parameters:
748 +  comm - the communicator that will share the eventually allocated matrix
749 .  nrows - the number of LOCAL rows in the matrix
750 -  ncols - the number of LOCAL columns in the matrix
751 
752    Output Parameters:
753 +  dnz - the array that will be passed to the matrix preallocation routines
754 -  onz - the other array passed to the matrix preallocation routines
755 
756    Level: intermediate
757 
758    Notes:
759     See Users-Manual: ch_performance for more details.
760 
761    Do not malloc or free dnz and onz, that is handled internally by these routines
762 
763    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
764 
765 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
766           MatPreallocateSymmetricSetLocalBlock()
767 M*/
768 #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
769 do { \
770   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ncols = (ncols),__rstart,__start,__end; \
771   _4_ierr = PetscCalloc2((size_t)__nrows,&dnz,(size_t)__nrows,&onz);CHKERRQ(_4_ierr); \
772   __start = 0; __end = __start;                                         \
773   _4_ierr = MPI_Scan(&__ncols,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ncols;\
774   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; \
775   do { } while (0)
776 
777 /*MC
778    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
779        inserted using a local number of the rows and columns
780 
781    Synopsis:
782    #include <petscmat.h>
783    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
784 
785    Not Collective
786 
787    Input Parameters:
788 +  map - the row mapping from local numbering to global numbering
789 .  nrows - the number of rows indicated
790 .  rows - the indices of the rows
791 .  cmap - the column mapping from local to global numbering
792 .  ncols - the number of columns in the matrix
793 .  cols - the columns indicated
794 .  dnz - the array that will be passed to the matrix preallocation routines
795 -  onz - the other array passed to the matrix preallocation routines
796 
797    Level: intermediate
798 
799    Notes:
800     See Users-Manual: ch_performance for more details.
801 
802    Do not malloc or free dnz and onz, that is handled internally by these routines
803 
804 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
805           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups()
806 M*/
807 #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
808 do {\
809   PetscInt __l;\
810   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
811   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
812   for (__l=0;__l<nrows;__l++) {\
813     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
814   }\
815  } while (0)
816 
817 /*MC
818    MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be
819        inserted using a local number of the rows and columns. This version removes any duplicate columns in cols
820 
821    Synopsis:
822    #include <petscmat.h>
823    PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
824 
825    Not Collective
826 
827    Input Parameters:
828 +  map - the row mapping from local numbering to global numbering
829 .  nrows - the number of rows indicated
830 .  rows - the indices of the rows (these values are mapped to the global values)
831 .  cmap - the column mapping from local to global numbering
832 .  ncols - the number of columns in the matrix   (this value will be changed if duplicate columns are found)
833 .  cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed)
834 .  dnz - the array that will be passed to the matrix preallocation routines
835 -  onz - the other array passed to the matrix preallocation routines
836 
837    Level: intermediate
838 
839    Notes:
840     See Users-Manual: ch_performance for more details.
841 
842    Do not malloc or free dnz and onz, that is handled internally by these routines
843 
844 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
845           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
846 M*/
847 #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
848 do {\
849   PetscInt __l;\
850   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
851   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
852   _4_ierr = PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\
853   for (__l=0;__l<nrows;__l++) {\
854     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
855   }\
856  } while (0)
857 
858 /*MC
859    MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
860        inserted using a local number of the rows and columns
861 
862    Synopsis:
863    #include <petscmat.h>
864    PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
865 
866    Not Collective
867 
868    Input Parameters:
869 +  map - the row mapping from local numbering to global numbering
870 .  nrows - the number of rows indicated
871 .  rows - the indices of the rows
872 .  cmap - the column mapping from local to global numbering
873 .  ncols - the number of columns in the matrix
874 .  cols - the columns indicated
875 .  dnz - the array that will be passed to the matrix preallocation routines
876 -  onz - the other array passed to the matrix preallocation routines
877 
878    Level: intermediate
879 
880    Notes:
881     See Users-Manual: ch_performance for more details.
882 
883    Do not malloc or free dnz and onz, that is handled internally by these routines
884 
885 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
886           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
887 M*/
888 #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
889 do {\
890   PetscInt __l;\
891   _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
892   _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
893   for (__l=0;__l<nrows;__l++) {\
894     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
895   }\
896 } while (0)
897 
898 /*MC
899    MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
900        inserted using a local number of the rows and columns
901 
902    Synopsis:
903    #include <petscmat.h>
904    PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
905 
906    Not Collective
907 
908    Input Parameters:
909 +  map - the mapping between local numbering and global numbering
910 .  nrows - the number of rows indicated
911 .  rows - the indices of the rows
912 .  ncols - the number of columns in the matrix
913 .  cols - the columns indicated
914 .  dnz - the array that will be passed to the matrix preallocation routines
915 -  onz - the other array passed to the matrix preallocation routines
916 
917    Level: intermediate
918 
919    Notes:
920     See Users-Manual: ch_performance for more details.
921 
922    Do not malloc or free dnz and onz that is handled internally by these routines
923 
924 .seealso: MatPreallocateFinalize(), MatPreallocateSet()
925           MatPreallocateInitialize(),  MatPreallocateSetLocal()
926 M*/
927 #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\
928 do {\
929   PetscInt __l;\
930   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
931   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
932   for (__l=0;__l<nrows;__l++) {\
933     _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
934   }\
935 } while (0)
936 
937 /*MC
938    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
939        inserted using a local number of the rows and columns
940 
941    Synopsis:
942    #include <petscmat.h>
943    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
944 
945    Not Collective
946 
947    Input Parameters:
948 +  row - the row
949 .  ncols - the number of columns in the matrix
950 -  cols - the columns indicated
951 
952    Output Parameters:
953 +  dnz - the array that will be passed to the matrix preallocation routines
954 -  onz - the other array passed to the matrix preallocation routines
955 
956    Level: intermediate
957 
958    Notes:
959     See Users-Manual: ch_performance for more details.
960 
961    Do not malloc or free dnz and onz that is handled internally by these routines
962 
963    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
964 
965 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
966           MatPreallocateInitialize(), MatPreallocateSetLocal()
967 M*/
968 #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
969 do { PetscInt __i; \
970   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);\
971   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);\
972   for (__i=0; __i<nc; __i++) {\
973     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
974     else if (dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\
975   }\
976 } while (0)
977 
978 /*MC
979    MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
980        inserted using a local number of the rows and columns
981 
982    Synopsis:
983    #include <petscmat.h>
984    PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
985 
986    Not Collective
987 
988    Input Parameters:
989 +  nrows - the number of rows indicated
990 .  rows - the indices of the rows
991 .  ncols - the number of columns in the matrix
992 .  cols - the columns indicated
993 .  dnz - the array that will be passed to the matrix preallocation routines
994 -  onz - the other array passed to the matrix preallocation routines
995 
996    Level: intermediate
997 
998    Notes:
999     See Users-Manual: ch_performance for more details.
1000 
1001    Do not malloc or free dnz and onz that is handled internally by these routines
1002 
1003    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
1004 
1005 .seealso: MatPreallocateFinalize(), MatPreallocateSet(),  MatPreallocateInitialize(),
1006           MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
1007 M*/
1008 #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\
1009 do { PetscInt __i; \
1010   for (__i=0; __i<nc; __i++) {\
1011     if (cols[__i] >= __end) onz[row - __rstart]++; \
1012     else if (cols[__i] >= row && dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\
1013   }\
1014 } while (0)
1015 
1016 /*MC
1017    MatPreallocateLocation -  An alternative to MatPreallocateSet() that puts the nonzero locations into the matrix if it exists
1018 
1019    Synopsis:
1020    #include <petscmat.h>
1021    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
1022 
1023    Not Collective
1024 
1025    Input Parameters:
1026 +  A - matrix
1027 .  row - row where values exist (must be local to this process)
1028 .  ncols - number of columns
1029 .  cols - columns with nonzeros
1030 .  dnz - the array that will be passed to the matrix preallocation routines
1031 -  onz - the other array passed to the matrix preallocation routines
1032 
1033    Level: intermediate
1034 
1035    Notes:
1036     See Users-Manual: ch_performance for more details.
1037 
1038    Do not malloc or free dnz and onz that is handled internally by these routines
1039 
1040    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
1041 
1042 .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
1043           MatPreallocateSymmetricSetLocalBlock()
1044 M*/
1045 #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)
1046 
1047 
1048 /*MC
1049    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
1050        row in a matrix providing the data that one can use to correctly preallocate the matrix.
1051 
1052    Synopsis:
1053    #include <petscmat.h>
1054    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
1055 
1056    Collective
1057 
1058    Input Parameters:
1059 +  dnz - the array that was be passed to the matrix preallocation routines
1060 -  onz - the other array passed to the matrix preallocation routines
1061 
1062    Level: intermediate
1063 
1064    Notes:
1065     See Users-Manual: ch_performance for more details.
1066 
1067    Do not malloc or free dnz and onz that is handled internally by these routines
1068 
1069    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
1070 
1071 .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
1072           MatPreallocateSymmetricSetLocalBlock()
1073 M*/
1074 #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} while (0)
1075 
1076 /* Routines unique to particular data structures */
1077 PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
1078 
1079 PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
1080 PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
1081 
1082 PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
1083 PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
1084 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1085 PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1086 PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1087 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
1088 
1089 #define MAT_SKIP_ALLOCATION -4
1090 
1091 PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1092 PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1093 PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
1094 PETSC_EXTERN PetscErrorCode MatSeqAIJSetTotalPreallocation(Mat,PetscInt);
1095 
1096 PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1097 PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1098 PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1099 PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
1100 PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1101 PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
1102 PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1103 PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
1104 PETSC_EXTERN PetscErrorCode MatMPIAdjToSeq(Mat,Mat*);
1105 PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
1106 PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
1107 PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1108 PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1109 PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
1110 
1111 PETSC_EXTERN PetscErrorCode MatDenseGetLDA(Mat,PetscInt*);
1112 PETSC_EXTERN PetscErrorCode MatDenseSetLDA(Mat,PetscInt);
1113 PETSC_DEPRECATED_FUNCTION("Use MatDenseSetLDA() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode MatSeqDenseSetLDA(Mat A,PetscInt lda) {return MatDenseSetLDA(A,lda);}
1114 PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
1115 
1116 PETSC_EXTERN PetscErrorCode MatBlockMatSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1117 
1118 PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
1119 PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
1120 
1121 PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
1122 
1123 PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
1124 PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*);
1125 /*
1126   These routines are not usually accessed directly, rather solving is
1127   done through the KSP and PC interfaces.
1128 */
1129 
1130 /*J
1131     MatOrderingType - String with the name of a PETSc matrix ordering
1132 
1133    Level: beginner
1134 
1135    Notes:
1136       If MATORDERINGEXTERNAL is used then PETSc does not compute an ordering and utilizes one built into the factorization package
1137 
1138 .seealso: MatGetOrdering()
1139 J*/
1140 typedef const char* MatOrderingType;
1141 #define MATORDERINGNATURAL        "natural"
1142 #define MATORDERINGND             "nd"
1143 #define MATORDERING1WD            "1wd"
1144 #define MATORDERINGRCM            "rcm"
1145 #define MATORDERINGQMD            "qmd"
1146 #define MATORDERINGROWLENGTH      "rowlength"
1147 #define MATORDERINGWBM            "wbm"
1148 #define MATORDERINGSPECTRAL       "spectral"
1149 #define MATORDERINGAMD            "amd"            /* only works if UMFPACK is installed with PETSc */
1150 #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 */
1151 #define MATORDERINGEXTERNAL       "external"       /* uses an ordering type internal to the factorization package */
1152 
1153 PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
1154 PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
1155 PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
1156 PETSC_EXTERN PetscFunctionList MatOrderingList;
1157 
1158 PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1159 PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);
1160 
1161 /*S
1162     MatFactorShiftType - Numeric Shift for factorizations
1163 
1164    Level: beginner
1165 
1166 .seealso: MatGetFactor()
1167 S*/
1168 typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1169 PETSC_EXTERN const char *const MatFactorShiftTypes[];
1170 PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];
1171 
1172 /*S
1173     MatFactorError - indicates what type of error was generated in a matrix factorization
1174 
1175     Level: beginner
1176 
1177     Developer Notes:
1178     Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1179 
1180 .seealso: MatGetFactor()
1181 S*/
1182 typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError;
1183 
1184 PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*);
1185 PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat);
1186 PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*);
1187 
1188 /*S
1189    MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization
1190 
1191    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
1192 $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1193 
1194    Notes:
1195     These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1196 
1197       You can use MatFactorInfoInitialize() to set default values.
1198 
1199    Level: developer
1200 
1201 .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1202           MatFactorInfoInitialize()
1203 
1204 S*/
1205 typedef struct {
1206   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
1207   PetscReal     usedt;
1208   PetscReal     dt;             /* drop tolerance */
1209   PetscReal     dtcol;          /* tolerance for pivoting */
1210   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
1211   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1212   PetscReal     levels;         /* ICC/ILU(levels) */
1213   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1214                                    factorization may be faster if do not pivot */
1215   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1216   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1217   PetscReal     shiftamount;     /* how large the shift is */
1218 } MatFactorInfo;
1219 
1220 PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
1221 PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
1222 PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1223 PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
1224 PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
1225 PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
1226 PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1227 PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1228 PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1229 PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
1230 PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1231 PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1232 PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1233 PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1234 PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1235 PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1236 PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1237 PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1238 PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
1239 PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
1240 
1241 typedef enum {MAT_FACTOR_SCHUR_UNFACTORED, MAT_FACTOR_SCHUR_FACTORED, MAT_FACTOR_SCHUR_INVERTED} MatFactorSchurStatus;
1242 PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS);
1243 PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1244 PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*,MatFactorSchurStatus);
1245 PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat);
1246 PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1247 PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec);
1248 PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec);
1249 PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat);
1250 
1251 /*E
1252     MatSORType - What type of (S)SOR to perform
1253 
1254     Level: beginner
1255 
1256    May be bitwise ORd together
1257 
1258    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1259 
1260    MatSORType may be bitwise ORd together, so do not change the numbers
1261 
1262 .seealso: MatSOR()
1263 E*/
1264 typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1265               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1266               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
1267               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1268 PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1269 
1270 /*
1271     These routines are for efficiently computing Jacobians via finite differences.
1272 */
1273 
1274 /*S
1275      MatColoring - Object for managing the coloring of matrices.
1276 
1277    Level: beginner
1278 
1279     Notes:
1280        Coloring of matrices can be computed directly from the sparse matrix nonzero structure via the MatColoring object or from the mesh from which the
1281        matrix comes from via DMCreateColoring(). In general using the mesh produces a more optimal coloring (fewer colors).
1282 
1283        Once a coloring is available MatFDColoringCreate() creates an object that can be used to efficiently compute Jacobians using that coloring. This
1284        same object can also be used to efficiently convert data created by Automatic Differentation tools to PETSc sparse matrices.
1285 
1286 .seealso:  MatFDColoringCreate(), MatColoringWeightType, ISColoring, MatFDColoring, DMCreateColoring(), MatColoringCreate(), MatOrdering, MatPartitioning
1287 S*/
1288 typedef struct _p_MatColoring* MatColoring;
1289 
1290 /*J
1291     MatColoringType - String with the name of a PETSc matrix coloring
1292 
1293    Level: beginner
1294 
1295 .seealso: MatColoringSetType(), MatColoring
1296 J*/
1297 typedef const  char*           MatColoringType;
1298 #define MATCOLORINGJP      "jp"
1299 #define MATCOLORINGPOWER   "power"
1300 #define MATCOLORINGNATURAL "natural"
1301 #define MATCOLORINGSL      "sl"
1302 #define MATCOLORINGLF      "lf"
1303 #define MATCOLORINGID      "id"
1304 #define MATCOLORINGGREEDY  "greedy"
1305 
1306 /*E
1307    MatColoringWeightType - Type of weight scheme
1308 
1309     Not Collective
1310 
1311 +   MAT_COLORING_RANDOM  - Random weights
1312 .   MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
1313 -   MAT_COLORING_LF      - Last-first weighting.
1314 
1315     Level: intermediate
1316 
1317    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1318 
1319 .seealso: MatColoring, MatColoringCreate()
1320 E*/
1321 typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;
1322 
1323 PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1324 PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*);
1325 PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1326 PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1327 PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1328 PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1329 PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1330 PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1331 PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1332 PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1333 PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1334 PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1335 PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1336 PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType);
1337 PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*);
1338 PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm);
1339 PETSC_EXTERN PetscErrorCode MatColoringTest(MatColoring,ISColoring);
1340 PETSC_DEPRECATED_FUNCTION("Use MatColoringTest() (since version 3.10)") PETSC_STATIC_INLINE PetscErrorCode MatColoringTestValid(MatColoring matcoloring,ISColoring iscoloring) {return MatColoringTest(matcoloring,iscoloring);}
1341 PETSC_EXTERN PetscErrorCode MatISColoringTest(Mat,ISColoring);
1342 
1343 /*S
1344      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1345         and coloring
1346 
1347    Level: beginner
1348 
1349    Notes:
1350       This object is creating utilizing a coloring provided by the MatColoring object or DMCreateColoring()
1351 
1352 .seealso:  MatFDColoringCreate(), MatColoring, DMCreateColoring()
1353 S*/
1354 typedef struct _p_MatFDColoring* MatFDColoring;
1355 
1356 PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1357 PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1358 PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1359 PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1360 PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1361 PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1362 PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1363 PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *);
1364 PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1365 PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,const PetscInt*[]);
1366 PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1367 PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1368 PETSC_EXTERN PetscErrorCode MatFDColoringSetValues(Mat,MatFDColoring,const PetscScalar*);
1369 
1370 /*S
1371      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1372 
1373    Level: beginner
1374 
1375 .seealso:  MatTransposeColoringCreate()
1376 S*/
1377 typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1378 
1379 PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1380 PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1381 PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1382 PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1383 
1384 /*
1385     These routines are for partitioning matrices: currently used only
1386   for adjacency matrix, MatCreateMPIAdj().
1387 */
1388 
1389 /*S
1390      MatPartitioning - Object for managing the partitioning of a matrix or graph
1391 
1392    Level: beginner
1393 
1394    Notes:
1395      There is also a PetscPartitioner object that provides the same functionality. It can utilize the MatPartitioning operations
1396      via PetscPartitionerSetType(p,PETSCPARTITIONERMATPARTITIONING)
1397 
1398    Developers Note:
1399      It is an extra maintainance and documentation cost to have two objects with the same functionality.
1400 
1401 .seealso:  MatPartitioningCreate(), MatPartitioningType, MatColoring, MatGetOrdering()
1402 S*/
1403 typedef struct _p_MatPartitioning* MatPartitioning;
1404 
1405 /*J
1406     MatPartitioningType - String with the name of a PETSc matrix partitioning
1407 
1408    Level: beginner
1409 dm
1410 .seealso: MatPartitioningCreate(), MatPartitioning
1411 J*/
1412 typedef const char* MatPartitioningType;
1413 #define MATPARTITIONINGCURRENT  "current"
1414 #define MATPARTITIONINGAVERAGE   "average"
1415 #define MATPARTITIONINGSQUARE   "square"
1416 #define MATPARTITIONINGPARMETIS "parmetis"
1417 #define MATPARTITIONINGCHACO    "chaco"
1418 #define MATPARTITIONINGPARTY    "party"
1419 #define MATPARTITIONINGPTSCOTCH "ptscotch"
1420 #define MATPARTITIONINGHIERARCH  "hierarch"
1421 
1422 PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1423 PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1424 PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1425 PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1426 PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1427 PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1428 PETSC_EXTERN PetscErrorCode MatPartitioningSetUseEdgeWeights(MatPartitioning,PetscBool);
1429 PETSC_EXTERN PetscErrorCode MatPartitioningGetUseEdgeWeights(MatPartitioning,PetscBool*);
1430 PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1431 PETSC_EXTERN PetscErrorCode MatPartitioningImprove(MatPartitioning,IS*);
1432 PETSC_EXTERN PetscErrorCode MatPartitioningViewImbalance(MatPartitioning,IS);
1433 PETSC_EXTERN PetscErrorCode MatPartitioningApplyND(MatPartitioning,IS*);
1434 PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
1435 PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
1436 PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1437 PETSC_EXTERN PetscErrorCode MatPartitioningViewFromOptions(MatPartitioning,PetscObject,const char[]);
1438 PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
1439 PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1440 
1441 PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1442 PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1443 PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
1444 
1445 typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
1446 PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1447 typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
1448 PETSC_EXTERN const char *const MPChacoLocalTypes[];
1449 typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
1450 PETSC_EXTERN const char *const MPChacoEigenTypes[];
1451 
1452 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1453 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1454 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1455 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1456 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1457 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1458 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1459 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1460 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1461 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1462 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
1463 
1464 #define MP_PARTY_OPT "opt"
1465 #define MP_PARTY_LIN "lin"
1466 #define MP_PARTY_SCA "sca"
1467 #define MP_PARTY_RAN "ran"
1468 #define MP_PARTY_GBF "gbf"
1469 #define MP_PARTY_GCF "gcf"
1470 #define MP_PARTY_BUB "bub"
1471 #define MP_PARTY_DEF "def"
1472 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
1473 #define MP_PARTY_HELPFUL_SETS "hs"
1474 #define MP_PARTY_KERNIGHAN_LIN "kl"
1475 #define MP_PARTY_NONE "no"
1476 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1477 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1478 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1479 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
1480 
1481 typedef enum { MP_PTSCOTCH_DEFAULT,MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
1482 PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1483 
1484 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1485 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1486 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1487 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
1488 
1489 /*
1490  * hierarchical partitioning
1491  */
1492 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*);
1493 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*);
1494 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt);
1495 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);
1496 
1497 PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1498 PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1499 
1500 /*
1501     If you add entries here you must also add them to include/petscmat.h
1502     and src/mat/f90-mod/petscmat.h
1503 */
1504 typedef enum { MATOP_SET_VALUES=0,
1505                MATOP_GET_ROW=1,
1506                MATOP_RESTORE_ROW=2,
1507                MATOP_MULT=3,
1508                MATOP_MULT_ADD=4,
1509                MATOP_MULT_TRANSPOSE=5,
1510                MATOP_MULT_TRANSPOSE_ADD=6,
1511                MATOP_SOLVE=7,
1512                MATOP_SOLVE_ADD=8,
1513                MATOP_SOLVE_TRANSPOSE=9,
1514                MATOP_SOLVE_TRANSPOSE_ADD=10,
1515                MATOP_LUFACTOR=11,
1516                MATOP_CHOLESKYFACTOR=12,
1517                MATOP_SOR=13,
1518                MATOP_TRANSPOSE=14,
1519                MATOP_GETINFO=15,
1520                MATOP_EQUAL=16,
1521                MATOP_GET_DIAGONAL=17,
1522                MATOP_DIAGONAL_SCALE=18,
1523                MATOP_NORM=19,
1524                MATOP_ASSEMBLY_BEGIN=20,
1525                MATOP_ASSEMBLY_END=21,
1526                MATOP_SET_OPTION=22,
1527                MATOP_ZERO_ENTRIES=23,
1528                MATOP_ZERO_ROWS=24,
1529                MATOP_LUFACTOR_SYMBOLIC=25,
1530                MATOP_LUFACTOR_NUMERIC=26,
1531                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1532                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1533                MATOP_SETUP_PREALLOCATION=29,
1534                MATOP_ILUFACTOR_SYMBOLIC=30,
1535                MATOP_ICCFACTOR_SYMBOLIC=31,
1536                MATOP_GET_DIAGONAL_BLOCK=32,
1537                MATOP_FREE_INTER_STRUCT=33,
1538                MATOP_DUPLICATE=34,
1539                MATOP_FORWARD_SOLVE=35,
1540                MATOP_BACKWARD_SOLVE=36,
1541                MATOP_ILUFACTOR=37,
1542                MATOP_ICCFACTOR=38,
1543                MATOP_AXPY=39,
1544                MATOP_CREATE_SUBMATRICES=40,
1545                MATOP_INCREASE_OVERLAP=41,
1546                MATOP_GET_VALUES=42,
1547                MATOP_COPY=43,
1548                MATOP_GET_ROW_MAX=44,
1549                MATOP_SCALE=45,
1550                MATOP_SHIFT=46,
1551                MATOP_DIAGONAL_SET=47,
1552                MATOP_ZERO_ROWS_COLUMNS=48,
1553                MATOP_SET_RANDOM=49,
1554                MATOP_GET_ROW_IJ=50,
1555                MATOP_RESTORE_ROW_IJ=51,
1556                MATOP_GET_COLUMN_IJ=52,
1557                MATOP_RESTORE_COLUMN_IJ=53,
1558                MATOP_FDCOLORING_CREATE=54,
1559                MATOP_COLORING_PATCH=55,
1560                MATOP_SET_UNFACTORED=56,
1561                MATOP_PERMUTE=57,
1562                MATOP_SET_VALUES_BLOCKED=58,
1563                MATOP_CREATE_SUBMATRIX=59,
1564                MATOP_DESTROY=60,
1565                MATOP_VIEW=61,
1566                MATOP_CONVERT_FROM=62,
1567                MATOP_MATMAT_MULT=63,
1568                MATOP_MATMAT_MULT_SYMBOLIC=64,
1569                MATOP_MATMAT_MULT_NUMERIC=65,
1570                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1571                MATOP_SET_VALUES_LOCAL=67,
1572                MATOP_ZERO_ROWS_LOCAL=68,
1573                MATOP_GET_ROW_MAX_ABS=69,
1574                MATOP_GET_ROW_MIN_ABS=70,
1575                MATOP_CONVERT=71,
1576                MATOP_SET_COLORING=72,
1577                /* MATOP_PLACEHOLDER_73=73, */
1578                MATOP_SET_VALUES_ADIFOR=74,
1579                MATOP_FD_COLORING_APPLY=75,
1580                MATOP_SET_FROM_OPTIONS=76,
1581                MATOP_MULT_CONSTRAINED=77,
1582                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
1583                MATOP_FIND_ZERO_DIAGONALS=79,
1584                MATOP_MULT_MULTIPLE=80,
1585                MATOP_SOLVE_MULTIPLE=81,
1586                MATOP_GET_INERTIA=82,
1587                MATOP_LOAD=83,
1588                MATOP_IS_SYMMETRIC=84,
1589                MATOP_IS_HERMITIAN=85,
1590                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1591                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1592                MATOP_CREATE_VECS=88,
1593                MATOP_MAT_MULT=89,
1594                MATOP_MAT_MULT_SYMBOLIC=90,
1595                MATOP_MAT_MULT_NUMERIC=91,
1596                MATOP_PTAP=92,
1597                MATOP_PTAP_SYMBOLIC=93,
1598                MATOP_PTAP_NUMERIC=94,
1599                MATOP_MAT_TRANSPOSE_MULT=95,
1600                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1601                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
1602                /* MATOP_PLACEHOLDER_98=98, */
1603                MATOP_PRODUCTSETFROMOPTIONS=99,
1604                MATOP_PRODUCTSYMBOLIC=100,
1605                MATOP_PRODUCTNUMERIC=101,
1606                MATOP_CONJUGATE=102,
1607                /* MATOP_PLACEHOLDER_103=103, */
1608                MATOP_SET_VALUES_ROW=104,
1609                MATOP_REAL_PART=105,
1610                MATOP_IMAGINARY_PART=106,
1611                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1612                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1613                MATOP_MAT_SOLVE=109,
1614                MATOP_MAT_SOLVE_TRANSPOSE=110,
1615                MATOP_GET_ROW_MIN=111,
1616                MATOP_GET_COLUMN_VECTOR=112,
1617                MATOP_MISSING_DIAGONAL=113,
1618                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
1619                MATOP_CREATE=115,
1620                MATOP_GET_GHOSTS=116,
1621                MATOP_GET_LOCAL_SUB_MATRIX=117,
1622                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1623                MATOP_MULT_DIAGONAL_BLOCK=119,
1624                MATOP_HERMITIAN_TRANSPOSE=120,
1625                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
1626                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
1627                MATOP_GET_MULTI_PROC_BLOCK=123,
1628                MATOP_FIND_NONZERO_ROWS=124,
1629                MATOP_GET_COLUMN_NORMS=125,
1630                MATOP_INVERT_BLOCK_DIAGONAL=126,
1631                /* MATOP_PLACEHOLDER_127=127, */
1632                MATOP_CREATE_SUB_MATRICES_MPI=128,
1633                MATOP_SET_VALUES_BATCH=129,
1634                MATOP_TRANSPOSE_MAT_MULT=130,
1635                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
1636                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
1637                MATOP_TRANSPOSE_COLORING_CREAT=133,
1638                MATOP_TRANS_COLORING_APPLY_SPT=134,
1639                MATOP_TRANS_COLORING_APPLY_DEN=135,
1640                MATOP_RART=136,
1641                MATOP_RART_SYMBOLIC=137,
1642                MATOP_RART_NUMERIC=138,
1643                MATOP_SET_BLOCK_SIZES=139,
1644                MATOP_AYPX=140,
1645                MATOP_RESIDUAL=141,
1646                MATOP_FDCOLORING_SETUP=142,
1647                MATOP_MPICONCATENATESEQ=144,
1648                MATOP_DESTROYSUBMATRICES=145,
1649                MATOP_TRANSPOSE_SOLVE=146,
1650                MATOP_GET_VALUES_LOCAL=147
1651              } MatOperation;
1652 PETSC_EXTERN PetscErrorCode MatSetOperation(Mat,MatOperation,void(*)(void));
1653 PETSC_EXTERN PetscErrorCode MatGetOperation(Mat,MatOperation,void(**)(void));
1654 PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool*);
1655 PETSC_EXTERN PetscErrorCode MatHasCongruentLayouts(Mat,PetscBool*);
1656 PETSC_DEPRECATED_FUNCTION("Use MatProductClear() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode MatFreeIntermediateDataStructures(Mat A) {return MatProductClear(A);}
1657 PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1658 PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1659 PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1660 PETSC_EXTERN PetscErrorCode MatShellSetVecType(Mat,VecType);
1661 PETSC_EXTERN PetscErrorCode MatShellTestMult(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*);
1662 PETSC_EXTERN PetscErrorCode MatShellTestMultTranspose(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*);
1663 PETSC_EXTERN PetscErrorCode MatShellSetManageScalingShifts(Mat);
1664 PETSC_EXTERN PetscErrorCode MatShellSetMatProductOperation(Mat,MatProductType,PetscErrorCode (*)(Mat,Mat,Mat,void**),PetscErrorCode (*)(Mat,Mat,Mat,void*),PetscErrorCode (*)(void*),MatType,MatType);
1665 PETSC_EXTERN PetscErrorCode MatIsShell(Mat,PetscBool*);
1666 
1667 /*
1668    Codes for matrices stored on disk. By default they are
1669    stored in a universal format. By changing the format with
1670    PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
1671    be stored in a way natural for the matrix, for example dense matrices
1672    would be stored as dense. Matrices stored this way may only be
1673    read into matrices of the same type.
1674 */
1675 #define MATRIX_BINARY_FORMAT_DENSE -1
1676 
1677 PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1678 
1679 PETSC_EXTERN PetscErrorCode MatISSetLocalMatType(Mat,MatType);
1680 PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1681 PETSC_EXTERN PetscErrorCode MatISStoreL2L(Mat,PetscBool);
1682 PETSC_EXTERN PetscErrorCode MatISFixLocalEmpty(Mat,PetscBool);
1683 PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1684 PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat,Mat*);
1685 PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
1686 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use the MatConvert() interface (since version 3.10)") PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*);
1687 
1688 /*S
1689      MatNullSpace - Object that removes a null space from a vector, i.e.
1690          orthogonalizes the vector to a subsapce
1691 
1692    Level: advanced
1693 
1694   Users manual sections:
1695 .   sec_singular
1696 
1697 .seealso:  MatNullSpaceCreate()
1698 S*/
1699 typedef struct _p_MatNullSpace* MatNullSpace;
1700 
1701 PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1702 PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1703 PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1704 PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1705 PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1706 PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *);
1707 PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace);
1708 PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1709 PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1710 PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1711 PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1712 PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1713 PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1714 PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
1715 
1716 PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1717 PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1718 PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1719 PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
1720 
1721 PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1722 PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1723 PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1724 
1725 PETSC_EXTERN PetscErrorCode MatComputeOperator(Mat,MatType,Mat*);
1726 PETSC_EXTERN PetscErrorCode MatComputeOperatorTranspose(Mat,MatType,Mat*);
1727 
1728 PETSC_DEPRECATED_FUNCTION("Use MatComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode MatComputeExplicitOperator(Mat A,Mat* B) { return MatComputeOperator(A,NULL,B); }
1729 PETSC_DEPRECATED_FUNCTION("Use MatComputeOperatorTranspose() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode MatComputeExplicitOperatorTranspose(Mat A,Mat* B) { return MatComputeOperatorTranspose(A,NULL,B); }
1730 
1731 PETSC_EXTERN PetscErrorCode MatCreateKAIJ(Mat,PetscInt,PetscInt,const PetscScalar[],const PetscScalar[],Mat*);
1732 PETSC_EXTERN PetscErrorCode MatKAIJGetAIJ(Mat,Mat*);
1733 PETSC_EXTERN PetscErrorCode MatKAIJGetS(Mat,PetscInt*,PetscInt*,PetscScalar**);
1734 PETSC_EXTERN PetscErrorCode MatKAIJGetSRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1735 PETSC_EXTERN PetscErrorCode MatKAIJRestoreS(Mat,PetscScalar**);
1736 PETSC_EXTERN PetscErrorCode MatKAIJRestoreSRead(Mat,const PetscScalar**);
1737 PETSC_EXTERN PetscErrorCode MatKAIJGetT(Mat,PetscInt*,PetscInt*,PetscScalar**);
1738 PETSC_EXTERN PetscErrorCode MatKAIJGetTRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1739 PETSC_EXTERN PetscErrorCode MatKAIJRestoreT(Mat,PetscScalar**);
1740 PETSC_EXTERN PetscErrorCode MatKAIJRestoreTRead(Mat,const PetscScalar**);
1741 PETSC_EXTERN PetscErrorCode MatKAIJSetAIJ(Mat,Mat);
1742 PETSC_EXTERN PetscErrorCode MatKAIJSetS(Mat,PetscInt,PetscInt,const PetscScalar[]);
1743 PETSC_EXTERN PetscErrorCode MatKAIJSetT(Mat,PetscInt,PetscInt,const PetscScalar[]);
1744 PETSC_EXTERN PetscErrorCode MatKAIJGetScaledIdentity(Mat,PetscBool*);
1745 
1746 PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
1747 
1748 PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1749 PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1750 PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1751 PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1752 PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1753 PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1754 PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1755 PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1756 PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1757 PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1758 PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1759 PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1760 PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1761 
1762 /*S
1763     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
1764               Jacobian vector products
1765 
1766     Notes:
1767     MATMFFD is a specific MatType which uses the MatMFFD data structure
1768 
1769            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
1770 
1771     Level: developer
1772 
1773 .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
1774 S*/
1775 typedef struct _p_MatMFFD* MatMFFD;
1776 
1777 /*J
1778     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1779 
1780    Level: beginner
1781 
1782 .seealso: MatMFFDSetType(), MatMFFDRegister()
1783 J*/
1784 typedef const char* MatMFFDType;
1785 #define MATMFFD_DS  "ds"
1786 #define MATMFFD_WP  "wp"
1787 
1788 PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1789 PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));
1790 
1791 PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1792 PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool);
1793 
1794 PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType);
1795 
1796 PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1797 PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
1798 
1799 #ifdef PETSC_HAVE_HARA
1800 PETSC_EXTERN_TYPEDEF typedef PetscScalar (*MatHaraKernel)(PetscInt,PetscReal[],PetscReal[],void*);
1801 PETSC_EXTERN PetscErrorCode MatCreateHaraFromKernel(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscReal[],MatHaraKernel,void*,PetscReal,PetscInt,PetscInt,Mat*);
1802 PETSC_EXTERN PetscErrorCode MatCreateHaraFromMat(Mat,PetscInt,const PetscReal[],PetscReal,PetscInt,PetscInt,PetscInt,PetscReal,Mat*);
1803 PETSC_EXTERN PetscErrorCode MatHaraSetSamplingMat(Mat,Mat,PetscInt,PetscReal);
1804 #endif
1805 
1806 /*
1807    PETSc interface to MUMPS
1808 */
1809 #ifdef PETSC_HAVE_MUMPS
1810 PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1811 PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*);
1812 PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
1813 PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*);
1814 
1815 PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1816 PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1817 PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1818 PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
1819 PETSC_EXTERN PetscErrorCode MatMumpsGetInverse(Mat,Mat);
1820 PETSC_EXTERN PetscErrorCode MatMumpsGetInverseTranspose(Mat,Mat);
1821 #endif
1822 
1823 /*
1824    PETSc interface to Mkl_Pardiso
1825 */
1826 #ifdef PETSC_HAVE_MKL_PARDISO
1827 PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1828 #endif
1829 
1830 /*
1831    PETSc interface to Mkl_CPardiso
1832 */
1833 #ifdef PETSC_HAVE_MKL_CPARDISO
1834 PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt);
1835 #endif
1836 
1837 /*
1838    PETSc interface to SUPERLU
1839 */
1840 #ifdef PETSC_HAVE_SUPERLU
1841 PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1842 #endif
1843 
1844 /*
1845    PETSc interface to SUPERLU_DIST
1846 */
1847 #ifdef PETSC_HAVE_SUPERLU_DIST
1848 PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*);
1849 #endif
1850 
1851 /*
1852    PETSc interface to STRUMPACK
1853 */
1854 #ifdef PETSC_HAVE_STRUMPACK
1855 /*E
1856     MatSTRUMPACKReordering - sparsity reducing ordering to be used in STRUMPACK
1857 
1858     Level: intermediate
1859 E*/
1860 typedef enum {MAT_STRUMPACK_NATURAL=0,
1861               MAT_STRUMPACK_METIS=1,
1862               MAT_STRUMPACK_PARMETIS=2,
1863               MAT_STRUMPACK_SCOTCH=3,
1864               MAT_STRUMPACK_PTSCOTCH=4,
1865               MAT_STRUMPACK_RCM=5} MatSTRUMPACKReordering;
1866 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetReordering(Mat,MatSTRUMPACKReordering);
1867 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool);
1868 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat,PetscReal);
1869 PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSRelTol() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat mat,PetscReal rtol) {return MatSTRUMPACKSetHSSRelTol(mat,rtol);}
1870 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat,PetscReal);
1871 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat,PetscInt);
1872 PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSMinSepSize() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat mat,PetscInt hssminsize) {return MatSTRUMPACKSetHSSMinSepSize(mat,hssminsize);}
1873 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat,PetscInt);
1874 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat,PetscInt);
1875 #endif
1876 
1877 
1878 PETSC_EXTERN PetscErrorCode MatBindToCPU(Mat,PetscBool);
1879 PETSC_DEPRECATED_FUNCTION("Use MatBindToCPU (since v3.13)") PETSC_STATIC_INLINE PetscErrorCode MatPinToCPU(Mat A,PetscBool flg) {return MatBindToCPU(A,flg);}
1880 
1881 #ifdef PETSC_HAVE_CUDA
1882 /*E
1883     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
1884     matrices.
1885 
1886     Not Collective
1887 
1888 +   MAT_CUSPARSE_CSR - Compressed Sparse Row
1889 .   MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later).
1890 -   MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).
1891 
1892     Level: intermediate
1893 
1894    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1895 
1896 .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1897 E*/
1898 
1899 typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;
1900 
1901 /* these will be strings associated with enumerated type defined above */
1902 PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
1903 
1904 /*E
1905     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
1906     matrices whose operation should use a particular storage format.
1907 
1908     Not Collective
1909 
1910 +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
1911 .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
1912 .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
1913 -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices
1914 
1915     Level: intermediate
1916 
1917 .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1918 E*/
1919 typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;
1920 
1921 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1922 PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1923 PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
1924 PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSESetGenerateTranspose(Mat,PetscBool);
1925 
1926 PETSC_EXTERN PetscErrorCode MatCreateDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
1927 PETSC_EXTERN PetscErrorCode MatCreateSeqDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
1928 PETSC_EXTERN PetscErrorCode MatMPIDenseCUDASetPreallocation(Mat,PetscScalar[]);
1929 PETSC_EXTERN PetscErrorCode MatSeqDenseCUDASetPreallocation(Mat,PetscScalar[]);
1930 PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayWrite(Mat,PetscScalar**);
1931 PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayRead(Mat,const PetscScalar**);
1932 PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArray(Mat,PetscScalar**);
1933 PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayWrite(Mat,PetscScalar**);
1934 PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayRead(Mat,const PetscScalar**);
1935 PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArray(Mat,PetscScalar**);
1936 PETSC_EXTERN PetscErrorCode MatDenseCUDAPlaceArray(Mat,const PetscScalar*);
1937 PETSC_EXTERN PetscErrorCode MatDenseCUDAReplaceArray(Mat,const PetscScalar*);
1938 PETSC_EXTERN PetscErrorCode MatDenseCUDAResetArray(Mat);
1939 #endif
1940 
1941 #if defined(PETSC_HAVE_VIENNACL)
1942 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1943 PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1944 #endif
1945 
1946 /*
1947    PETSc interface to FFTW
1948 */
1949 #if defined(PETSC_HAVE_FFTW)
1950 PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1951 PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
1952 PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*);
1953 #endif
1954 
1955 #if defined(PETSC_HAVE_SCALAPACK)
1956 PETSC_EXTERN PetscErrorCode MatCreateScaLAPACK(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1957 PETSC_EXTERN PetscErrorCode MatScaLAPACKSetBlockSizes(Mat,PetscInt,PetscInt);
1958 PETSC_EXTERN PetscErrorCode MatScaLAPACKGetBlockSizes(Mat,PetscInt*,PetscInt*);
1959 #endif
1960 
1961 PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1962 PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1963 PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1964 PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1965 PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1966 PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
1967 PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1968 PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1969 PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1970 
1971 PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
1972 PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);
1973 
1974 PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**);
1975 
1976 PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat);
1977 
1978 PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*);
1979 PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*);
1980 
1981 #endif
1982