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