xref: /petsc/include/petscmat.h (revision 939b8a20d39cb19df6ece01ab2377b8403c87da5)
1 /*
2      Include file for the matrix component of PETSc
3 */
4 #ifndef __PETSCMAT_H
5 #define __PETSCMAT_H
6 #include "petscvec.h"
7 PETSC_EXTERN_CXX_BEGIN
8 
9 /*S
10      Mat - Abstract PETSc matrix object
11 
12    Level: beginner
13 
14   Concepts: matrix; linear operator
15 
16 .seealso:  MatCreate(), MatType, MatSetType()
17 S*/
18 typedef struct _p_Mat*           Mat;
19 
20 /*J
21     MatType - String with the name of a PETSc matrix or the creation function
22        with an optional dynamic library name, for example
23        http://www.mcs.anl.gov/petsc/lib.a:mymatcreate()
24 
25    Level: beginner
26 
27 .seealso: MatSetType(), Mat, MatSolverPackage
28 J*/
29 #define MatType char*
30 #define MATSAME            "same"
31 #define MATMAIJ            "maij"
32 #define MATSEQMAIJ         "seqmaij"
33 #define MATMPIMAIJ         "mpimaij"
34 #define MATIS              "is"
35 #define MATAIJ             "aij"
36 #define MATSEQAIJ          "seqaij"
37 #define MATSEQAIJPTHREAD   "seqaijpthread"
38 #define MATAIJPTHREAD      "aijpthread"
39 #define MATMPIAIJ          "mpiaij"
40 #define MATAIJCRL          "aijcrl"
41 #define MATSEQAIJCRL       "seqaijcrl"
42 #define MATMPIAIJCRL       "mpiaijcrl"
43 #define MATAIJCUSP         "aijcusp"
44 #define MATSEQAIJCUSP      "seqaijcusp"
45 #define MATMPIAIJCUSP      "mpiaijcusp"
46 #define MATAIJPERM         "aijperm"
47 #define MATSEQAIJPERM      "seqaijperm"
48 #define MATMPIAIJPERM      "mpiaijperm"
49 #define MATSHELL           "shell"
50 #define MATDENSE           "dense"
51 #define MATSEQDENSE        "seqdense"
52 #define MATMPIDENSE        "mpidense"
53 #define MATBAIJ            "baij"
54 #define MATSEQBAIJ         "seqbaij"
55 #define MATMPIBAIJ         "mpibaij"
56 #define MATMPIADJ          "mpiadj"
57 #define MATSBAIJ           "sbaij"
58 #define MATSEQSBAIJ        "seqsbaij"
59 #define MATMPISBAIJ        "mpisbaij"
60 #define MATSEQBSTRM        "seqbstrm"
61 #define MATMPIBSTRM        "mpibstrm"
62 #define MATBSTRM           "bstrm"
63 #define MATSEQSBSTRM       "seqsbstrm"
64 #define MATMPISBSTRM       "mpisbstrm"
65 #define MATSBSTRM          "sbstrm"
66 #define MATDAAD            "daad"
67 #define MATMFFD            "mffd"
68 #define MATNORMAL          "normal"
69 #define MATLRC             "lrc"
70 #define MATSCATTER         "scatter"
71 #define MATBLOCKMAT        "blockmat"
72 #define MATCOMPOSITE       "composite"
73 #define MATFFT             "fft"
74 #define MATFFTW            "fftw"
75 #define MATSEQCUFFT        "seqcufft"
76 #define MATTRANSPOSEMAT    "transpose"
77 #define MATSCHURCOMPLEMENT "schurcomplement"
78 #define MATPYTHON          "python"
79 #define MATHYPRESTRUCT     "hyprestruct"
80 #define MATHYPRESSTRUCT    "hypresstruct"
81 #define MATSUBMATRIX       "submatrix"
82 #define MATLOCALREF        "localref"
83 #define MATNEST            "nest"
84 #define MATIJ              "ij"
85 
86 /*J
87     MatSolverPackage - String with the name of a PETSc matrix solver type.
88 
89     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
90        SuperLU or SuperLU_Dist etc.
91 
92 
93    Level: beginner
94 
95 .seealso: MatGetFactor(), Mat, MatSetType(), MatType
96 J*/
97 #define MatSolverPackage char*
98 #define MATSOLVERSPOOLES      "spooles"
99 #define MATSOLVERSUPERLU      "superlu"
100 #define MATSOLVERSUPERLU_DIST "superlu_dist"
101 #define MATSOLVERUMFPACK      "umfpack"
102 #define MATSOLVERCHOLMOD      "cholmod"
103 #define MATSOLVERESSL         "essl"
104 #define MATSOLVERLUSOL        "lusol"
105 #define MATSOLVERMUMPS        "mumps"
106 #define MATSOLVERPASTIX       "pastix"
107 #define MATSOLVERMATLAB       "matlab"
108 #define MATSOLVERPETSC        "petsc"
109 #define MATSOLVERPLAPACK      "plapack"
110 #define MATSOLVERBAS          "bas"
111 
112 #define MATSOLVERBSTRM        "bstrm"
113 #define MATSOLVERSBSTRM       "sbstrm"
114 
115 /*E
116     MatFactorType - indicates what type of factorization is requested
117 
118     Level: beginner
119 
120    Any additions/changes here MUST also be made in include/finclude/petscmat.h
121 
122 .seealso: MatSolverPackage, MatGetFactor()
123 E*/
124 typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
125 extern const char *const MatFactorTypes[];
126 
127 extern PetscErrorCode  MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
128 extern PetscErrorCode  MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
129 extern PetscErrorCode  MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
130 extern PetscErrorCode  MatGetFactorType(Mat,MatFactorType*);
131 
132 /* Logging support */
133 #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
134 extern PetscClassId  MAT_CLASSID;
135 extern PetscClassId  MAT_FDCOLORING_CLASSID;
136 extern PetscClassId  MAT_TRANSPOSECOLORING_CLASSID;
137 extern PetscClassId  MAT_PARTITIONING_CLASSID;
138 extern PetscClassId  MAT_COARSEN_CLASSID;
139 extern PetscClassId  MAT_NULLSPACE_CLASSID;
140 extern PetscClassId  MATMFFD_CLASSID;
141 
142 /*E
143     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
144      or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate
145      that the input matrix is to be replaced with the converted matrix.
146 
147     Level: beginner
148 
149    Any additions/changes here MUST also be made in include/finclude/petscmat.h
150 
151 .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
152 E*/
153 typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse;
154 
155 /*E
156     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
157      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
158 
159     Level: beginner
160 
161 .seealso: MatGetSeqNonzerostructure()
162 E*/
163 typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
164 
165 extern PetscErrorCode  MatInitializePackage(const char[]);
166 
167 extern PetscErrorCode  MatCreate(MPI_Comm,Mat*);
168 PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
169 PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
170 extern PetscErrorCode  MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
171 extern PetscErrorCode  MatSetType(Mat,const MatType);
172 extern PetscErrorCode  MatSetFromOptions(Mat);
173 extern PetscErrorCode  MatRegisterAll(const char[]);
174 extern PetscErrorCode  MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
175 extern PetscErrorCode  MatRegisterBaseName(const char[],const char[],const char[]);
176 extern PetscErrorCode  MatSetOptionsPrefix(Mat,const char[]);
177 extern PetscErrorCode  MatAppendOptionsPrefix(Mat,const char[]);
178 extern PetscErrorCode  MatGetOptionsPrefix(Mat,const char*[]);
179 
180 /*MC
181    MatRegisterDynamic - Adds a new matrix type
182 
183    Synopsis:
184    PetscErrorCode MatRegisterDynamic(const char *name,const char *path,const char *name_create,PetscErrorCode (*routine_create)(Mat))
185 
186    Not Collective
187 
188    Input Parameters:
189 +  name - name of a new user-defined matrix type
190 .  path - path (either absolute or relative) the library containing this solver
191 .  name_create - name of routine to create method context
192 -  routine_create - routine to create method context
193 
194    Notes:
195    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
196 
197    If dynamic libraries are used, then the fourth input argument (routine_create)
198    is ignored.
199 
200    Sample usage:
201 .vb
202    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
203                "MyMatCreate",MyMatCreate);
204 .ve
205 
206    Then, your solver can be chosen with the procedural interface via
207 $     MatSetType(Mat,"my_mat")
208    or at runtime via the option
209 $     -mat_type my_mat
210 
211    Level: advanced
212 
213    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
214          If your function is not being put into a shared library then use VecRegister() instead
215 
216 .keywords: Mat, register
217 
218 .seealso: MatRegisterAll(), MatRegisterDestroy()
219 
220 M*/
221 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
222 #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
223 #else
224 #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
225 #endif
226 
227 extern PetscBool  MatRegisterAllCalled;
228 extern PetscFList MatList;
229 extern PetscFList MatColoringList;
230 extern PetscFList MatPartitioningList;
231 extern PetscFList MatCoarsenList;
232 
233 /*E
234     MatStructure - Indicates if the matrix has the same nonzero structure
235 
236     Level: beginner
237 
238    Any additions/changes here MUST also be made in include/finclude/petscmat.h
239 
240 .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
241 E*/
242 typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN,SAME_PRECONDITIONER} MatStructure;
243 
244 extern PetscErrorCode  MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
245 extern PetscErrorCode  MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
246 extern PetscErrorCode  MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
247 PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
248 PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
249 PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
250 PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
251 PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
252 PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
253 PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
254 extern PetscErrorCode  MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
255 PetscPolymorphicFunction(MatCreateAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(comm,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
256 PetscPolymorphicFunction(MatCreateAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(comm,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
257 PetscPolymorphicFunction(MatCreateAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(comm,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
258 PetscPolymorphicFunction(MatCreateAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(comm,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
259 PetscPolymorphicSubroutine(MatCreateAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(comm,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
260 PetscPolymorphicSubroutine(MatCreateAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(comm,m,n,M,N,0,nnz,0,onz,A))
261 PetscPolymorphicSubroutine(MatCreateAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(comm,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
262 PetscPolymorphicFunction(MatCreateAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(PETSC_COMM_WORLD,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
263 PetscPolymorphicFunction(MatCreateAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(PETSC_COMM_WORLD,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
264 PetscPolymorphicFunction(MatCreateAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(PETSC_COMM_WORLD,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
265 PetscPolymorphicFunction(MatCreateAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N),(PETSC_COMM_WORLD,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
266 PetscPolymorphicSubroutine(MatCreateAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(PETSC_COMM_WORLD,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
267 PetscPolymorphicSubroutine(MatCreateAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(PETSC_COMM_WORLD,m,n,M,N,0,nnz,0,onz,A))
268 PetscPolymorphicSubroutine(MatCreateAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(PETSC_COMM_WORLD,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
269 extern PetscErrorCode  MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
270 extern PetscErrorCode  MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
271 
272 extern PetscErrorCode  MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
273 PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,nz,nnz,&A),Mat,A)
274 PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
275 PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
276 PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
277 PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
278 PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
279 PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
280 extern PetscErrorCode  MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
281 PetscPolymorphicFunction(MatCreateBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(comm,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
282 PetscPolymorphicFunction(MatCreateBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
283 PetscPolymorphicFunction(MatCreateBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(comm,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
284 PetscPolymorphicFunction(MatCreateBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
285 PetscPolymorphicSubroutine(MatCreateBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
286 PetscPolymorphicSubroutine(MatCreateBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(comm,bs,m,n,M,N,0,nnz,0,onz,A))
287 PetscPolymorphicSubroutine(MatCreateBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
288 PetscPolymorphicFunction(MatCreateBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
289 PetscPolymorphicFunction(MatCreateBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
290 PetscPolymorphicFunction(MatCreateBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
291 PetscPolymorphicFunction(MatCreateBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
292 PetscPolymorphicSubroutine(MatCreateBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
293 PetscPolymorphicSubroutine(MatCreateBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,A))
294 PetscPolymorphicSubroutine(MatCreateBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
295 extern PetscErrorCode  MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
296 
297 extern PetscErrorCode  MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
298 extern PetscErrorCode  MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
299 PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,nz,nnz,&A),Mat,A)
300 PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
301 PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
302 PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
303 PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
304 PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
305 PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
306 
307 extern PetscErrorCode  MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
308 PetscPolymorphicFunction(MatCreateSBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(comm,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
309 PetscPolymorphicFunction(MatCreateSBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
310 PetscPolymorphicFunction(MatCreateSBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(comm,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
311 PetscPolymorphicFunction(MatCreateSBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
312 PetscPolymorphicSubroutine(MatCreateSBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
313 PetscPolymorphicSubroutine(MatCreateSBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(comm,bs,m,n,M,N,0,nnz,0,onz,A))
314 PetscPolymorphicSubroutine(MatCreateSBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
315 PetscPolymorphicFunction(MatCreateSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
316 PetscPolymorphicFunction(MatCreateSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
317 PetscPolymorphicFunction(MatCreateSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
318 PetscPolymorphicFunction(MatCreateSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
319 PetscPolymorphicSubroutine(MatCreateSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
320 PetscPolymorphicSubroutine(MatCreateSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,A))
321 PetscPolymorphicSubroutine(MatCreateSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
322 extern PetscErrorCode  MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
323 extern PetscErrorCode  MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
324 extern PetscErrorCode  MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*);
325 
326 extern PetscErrorCode  MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
327 PetscPolymorphicFunction(MatCreateShell,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(comm,m,n,M,N,ctx,&A),Mat,A)
328 PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
329 extern PetscErrorCode  MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
330 extern PetscErrorCode  MatCreateNormal(Mat,Mat*);
331 PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
332 extern PetscErrorCode  MatCreateLRC(Mat,Mat,Mat,Mat*);
333 extern PetscErrorCode  MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
334 extern PetscErrorCode  MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
335 extern PetscErrorCode  MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
336 
337 extern PetscErrorCode  MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
338 extern PetscErrorCode  MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
339 extern PetscErrorCode  MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
340 extern PetscErrorCode  MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
341 
342 extern PetscErrorCode  MatCreateScatter(MPI_Comm,VecScatter,Mat*);
343 extern PetscErrorCode  MatScatterSetVecScatter(Mat,VecScatter);
344 extern PetscErrorCode  MatScatterGetVecScatter(Mat,VecScatter*);
345 extern PetscErrorCode  MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
346 extern PetscErrorCode  MatCompositeAddMat(Mat,Mat);
347 extern PetscErrorCode  MatCompositeMerge(Mat);
348 extern PetscErrorCode  MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
349 typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
350 extern PetscErrorCode  MatCompositeSetType(Mat,MatCompositeType);
351 
352 extern PetscErrorCode  MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],const MatType,Mat*);
353 extern PetscErrorCode  MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
354 
355 extern PetscErrorCode  MatCreateTranspose(Mat,Mat*);
356 extern PetscErrorCode  MatCreateSubMatrix(Mat,IS,IS,Mat*);
357 extern PetscErrorCode  MatSubMatrixUpdate(Mat,Mat,IS,IS);
358 extern PetscErrorCode  MatCreateLocalRef(Mat,IS,IS,Mat*);
359 
360 extern PetscErrorCode  MatPythonSetType(Mat,const char[]);
361 
362 extern PetscErrorCode  MatSetUp(Mat);
363 extern PetscErrorCode  MatDestroy(Mat*);
364 
365 extern PetscErrorCode  MatConjugate(Mat);
366 extern PetscErrorCode  MatRealPart(Mat);
367 extern PetscErrorCode  MatImaginaryPart(Mat);
368 extern PetscErrorCode  MatGetDiagonalBlock(Mat,Mat*);
369 extern PetscErrorCode  MatGetTrace(Mat,PetscScalar*);
370 extern PetscErrorCode  MatInvertBlockDiagonal(Mat,PetscScalar **);
371 
372 /* ------------------------------------------------------------*/
373 extern PetscErrorCode  MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
374 extern PetscErrorCode  MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
375 extern PetscErrorCode  MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
376 extern PetscErrorCode  MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
377 extern PetscErrorCode  MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
378 
379 /*S
380      MatStencil - Data structure (C struct) for storing information about a single row or
381         column of a matrix as indexed on an associated grid.
382 
383    Fortran usage is different, see MatSetValuesStencil() for details.
384 
385    Level: beginner
386 
387   Concepts: matrix; linear operator
388 
389 .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil()
390 S*/
391 typedef struct {
392   PetscInt k,j,i,c;
393 } MatStencil;
394 
395 extern PetscErrorCode  MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
396 extern PetscErrorCode  MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
397 extern PetscErrorCode  MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
398 
399 extern PetscErrorCode  MatSetColoring(Mat,ISColoring);
400 extern PetscErrorCode  MatSetValuesAdic(Mat,void*);
401 extern PetscErrorCode  MatSetValuesAdifor(Mat,PetscInt,void*);
402 
403 /*E
404     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
405      to continue to add values to it
406 
407     Level: beginner
408 
409 .seealso: MatAssemblyBegin(), MatAssemblyEnd()
410 E*/
411 typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
412 extern PetscErrorCode  MatAssemblyBegin(Mat,MatAssemblyType);
413 extern PetscErrorCode  MatAssemblyEnd(Mat,MatAssemblyType);
414 extern PetscErrorCode  MatAssembled(Mat,PetscBool *);
415 
416 
417 
418 /*E
419     MatOption - Options that may be set for a matrix and its behavior or storage
420 
421     Level: beginner
422 
423    Any additions/changes here MUST also be made in include/finclude/petscmat.h
424 
425 .seealso: MatSetOption()
426 E*/
427 typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
428               MAT_SYMMETRIC,
429               MAT_STRUCTURALLY_SYMMETRIC,
430               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
431               MAT_NEW_NONZERO_LOCATION_ERR,
432               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
433               MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES,
434               MAT_USE_INODES,
435               MAT_HERMITIAN,
436               MAT_SYMMETRY_ETERNAL,
437               MAT_CHECK_COMPRESSED_ROW,
438               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
439               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
440               MAT_SPD,MAT_NO_OFF_PROC_ENTRIES,MAT_NO_OFF_PROC_ZERO_ROWS,
441               NUM_MAT_OPTIONS} MatOption;
442 extern const char *MatOptions[];
443 extern PetscErrorCode  MatSetOption(Mat,MatOption,PetscBool );
444 extern PetscErrorCode  MatGetType(Mat,const MatType*);
445 PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
446 
447 extern PetscErrorCode  MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
448 extern PetscErrorCode  MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
449 extern PetscErrorCode  MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
450 extern PetscErrorCode  MatGetRowUpperTriangular(Mat);
451 extern PetscErrorCode  MatRestoreRowUpperTriangular(Mat);
452 extern PetscErrorCode  MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
453 extern PetscErrorCode  MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
454 extern PetscErrorCode  MatGetColumnVector(Mat,Vec,PetscInt);
455 extern PetscErrorCode  MatGetArray(Mat,PetscScalar *[]);
456 PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
457 extern PetscErrorCode  MatRestoreArray(Mat,PetscScalar *[]);
458 extern PetscErrorCode  MatGetBlockSize(Mat,PetscInt *);
459 PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
460 extern PetscErrorCode  MatSetBlockSize(Mat,PetscInt);
461 extern PetscErrorCode MatSetNThreads(Mat,PetscInt);
462 extern PetscErrorCode MatGetNThreads(Mat,PetscInt*);
463 
464 extern PetscErrorCode  MatMult(Mat,Vec,Vec);
465 extern PetscErrorCode  MatMultDiagonalBlock(Mat,Vec,Vec);
466 extern PetscErrorCode  MatMultAdd(Mat,Vec,Vec,Vec);
467 PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
468 extern PetscErrorCode  MatMultTranspose(Mat,Vec,Vec);
469 extern PetscErrorCode  MatMultHermitianTranspose(Mat,Vec,Vec);
470 extern PetscErrorCode  MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
471 PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscBool ,t)
472 PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscBool ,t)
473 extern PetscErrorCode  MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
474 extern PetscErrorCode  MatMultTransposeAdd(Mat,Vec,Vec,Vec);
475 extern PetscErrorCode  MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
476 PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
477 extern PetscErrorCode  MatMultConstrained(Mat,Vec,Vec);
478 extern PetscErrorCode  MatMultTransposeConstrained(Mat,Vec,Vec);
479 extern PetscErrorCode  MatMatSolve(Mat,Mat,Mat);
480 
481 /*E
482     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
483   its numerical values copied over or just its nonzero structure.
484 
485     Level: beginner
486 
487    Any additions/changes here MUST also be made in include/finclude/petscmat.h
488 
489 $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
490 $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
491 $                               have several matrices with the same nonzero pattern.
492 
493 .seealso: MatDuplicate()
494 E*/
495 typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
496 
497 extern PetscErrorCode  MatConvert(Mat,const MatType,MatReuse,Mat*);
498 PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
499 extern PetscErrorCode  MatDuplicate(Mat,MatDuplicateOption,Mat*);
500 PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
501 PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
502 
503 
504 extern PetscErrorCode  MatCopy(Mat,Mat,MatStructure);
505 extern PetscErrorCode  MatView(Mat,PetscViewer);
506 extern PetscErrorCode  MatIsSymmetric(Mat,PetscReal,PetscBool *);
507 PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscBool ,t)
508 PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscBool ,t)
509 extern PetscErrorCode  MatIsStructurallySymmetric(Mat,PetscBool *);
510 PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscBool ,t)
511 extern PetscErrorCode  MatIsHermitian(Mat,PetscReal,PetscBool *);
512 PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscBool ,t)
513 extern PetscErrorCode  MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
514 extern PetscErrorCode  MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
515 extern PetscErrorCode  MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
516 extern PetscErrorCode  MatLoad(Mat, PetscViewer);
517 
518 extern PetscErrorCode  MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
519 extern PetscErrorCode  MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
520 extern PetscErrorCode  MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
521 extern PetscErrorCode  MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
522 
523 /*S
524      MatInfo - Context of matrix information, used with MatGetInfo()
525 
526    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
527 
528    Level: intermediate
529 
530   Concepts: matrix^nonzero information
531 
532 .seealso:  MatGetInfo(), MatInfoType
533 S*/
534 typedef struct {
535   PetscLogDouble block_size;                         /* block size */
536   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
537   PetscLogDouble memory;                             /* memory allocated */
538   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
539   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
540   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
541   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
542 } MatInfo;
543 
544 /*E
545     MatInfoType - Indicates if you want information about the local part of the matrix,
546      the entire parallel matrix or the maximum over all the local parts.
547 
548     Level: beginner
549 
550    Any additions/changes here MUST also be made in include/finclude/petscmat.h
551 
552 .seealso: MatGetInfo(), MatInfo
553 E*/
554 typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
555 extern PetscErrorCode  MatGetInfo(Mat,MatInfoType,MatInfo*);
556 extern PetscErrorCode  MatGetDiagonal(Mat,Vec);
557 extern PetscErrorCode  MatGetRowMax(Mat,Vec,PetscInt[]);
558 extern PetscErrorCode  MatGetRowMin(Mat,Vec,PetscInt[]);
559 extern PetscErrorCode  MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
560 extern PetscErrorCode  MatGetRowMinAbs(Mat,Vec,PetscInt[]);
561 extern PetscErrorCode  MatGetRowSum(Mat,Vec);
562 extern PetscErrorCode  MatTranspose(Mat,MatReuse,Mat*);
563 PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
564 extern PetscErrorCode  MatHermitianTranspose(Mat,MatReuse,Mat*);
565 extern PetscErrorCode  MatPermute(Mat,IS,IS,Mat *);
566 PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
567 extern PetscErrorCode  MatDiagonalScale(Mat,Vec,Vec);
568 extern PetscErrorCode  MatDiagonalSet(Mat,Vec,InsertMode);
569 extern PetscErrorCode  MatEqual(Mat,Mat,PetscBool *);
570 PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscBool ,t)
571 extern PetscErrorCode  MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
572 extern PetscErrorCode  MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
573 extern PetscErrorCode  MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
574 extern PetscErrorCode  MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
575 
576 extern PetscErrorCode  MatNorm(Mat,NormType,PetscReal *);
577 PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
578 extern PetscErrorCode  MatGetColumnNorms(Mat,NormType,PetscReal *);
579 extern PetscErrorCode  MatZeroEntries(Mat);
580 extern PetscErrorCode  MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
581 extern PetscErrorCode  MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
582 extern PetscErrorCode  MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
583 extern PetscErrorCode  MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
584 extern PetscErrorCode  MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
585 extern PetscErrorCode  MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
586 
587 extern PetscErrorCode  MatUseScaledForm(Mat,PetscBool );
588 extern PetscErrorCode  MatScaleSystem(Mat,Vec,Vec);
589 extern PetscErrorCode  MatUnScaleSystem(Mat,Vec,Vec);
590 
591 extern PetscErrorCode  MatGetSize(Mat,PetscInt*,PetscInt*);
592 extern PetscErrorCode  MatGetLocalSize(Mat,PetscInt*,PetscInt*);
593 extern PetscErrorCode  MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
594 extern PetscErrorCode  MatGetOwnershipRanges(Mat,const PetscInt**);
595 extern PetscErrorCode  MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
596 extern PetscErrorCode  MatGetOwnershipRangesColumn(Mat,const PetscInt**);
597 
598 extern PetscErrorCode  MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
599 extern PetscErrorCode  MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
600 extern PetscErrorCode  MatDestroyMatrices(PetscInt,Mat *[]);
601 extern PetscErrorCode  MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
602 extern PetscErrorCode  MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
603 extern PetscErrorCode  MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
604 extern PetscErrorCode  MatGetSeqNonzeroStructure(Mat,Mat*);
605 extern PetscErrorCode  MatDestroySeqNonzeroStructure(Mat*);
606 
607 extern PetscErrorCode  MatCreateMPIAIJConcatenateSeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
608 extern PetscErrorCode  MatCreateMPIAIJConcatenateSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,Mat*);
609 extern PetscErrorCode  MatCreateMPIAIJConcatenateSeqAIJNumeric(MPI_Comm,Mat,PetscInt,Mat);
610 extern PetscErrorCode  MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
611 extern PetscErrorCode  MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
612 extern PetscErrorCode  MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
613 extern PetscErrorCode  MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
614 extern PetscErrorCode  MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
615 extern PetscErrorCode  MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
616 #if defined (PETSC_USE_CTABLE)
617 extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
618 #else
619 extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
620 #endif
621 extern PetscErrorCode  MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
622 
623 extern PetscErrorCode  MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
624 
625 extern PetscErrorCode  MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
626 extern PetscErrorCode  MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
627 extern PetscErrorCode  MatMatMultNumeric(Mat,Mat,Mat);
628 
629 extern PetscErrorCode  MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
630 extern PetscErrorCode  MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
631 extern PetscErrorCode  MatPtAPNumeric(Mat,Mat,Mat);
632 extern PetscErrorCode  MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
633 extern PetscErrorCode  MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
634 extern PetscErrorCode  MatRARtNumeric(Mat,Mat,Mat);
635 
636 extern PetscErrorCode  MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
637 extern PetscErrorCode  MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
638 extern PetscErrorCode  MatTransposetMatMultNumeric(Mat,Mat,Mat);
639 extern PetscErrorCode  MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
640 extern PetscErrorCode  MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
641 extern PetscErrorCode  MatMatTransposeMultNumeric(Mat,Mat,Mat);
642 
643 extern PetscErrorCode  MatAXPY(Mat,PetscScalar,Mat,MatStructure);
644 extern PetscErrorCode  MatAYPX(Mat,PetscScalar,Mat,MatStructure);
645 
646 extern PetscErrorCode  MatScale(Mat,PetscScalar);
647 extern PetscErrorCode  MatShift(Mat,PetscScalar);
648 
649 extern PetscErrorCode  MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
650 extern PetscErrorCode  MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
651 extern PetscErrorCode  MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
652 extern PetscErrorCode  MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
653 extern PetscErrorCode  MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
654 extern PetscErrorCode  MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
655 extern PetscErrorCode  MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
656 extern PetscErrorCode  MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
657 extern PetscErrorCode  MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
658 extern PetscErrorCode  MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
659 
660 extern PetscErrorCode  MatStashSetInitialSize(Mat,PetscInt,PetscInt);
661 extern PetscErrorCode  MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
662 
663 extern PetscErrorCode  MatInterpolate(Mat,Vec,Vec);
664 extern PetscErrorCode  MatInterpolateAdd(Mat,Vec,Vec,Vec);
665 PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
666 extern PetscErrorCode  MatRestrict(Mat,Vec,Vec);
667 extern PetscErrorCode  MatGetVecs(Mat,Vec*,Vec*);
668 extern PetscErrorCode  MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
669 extern PetscErrorCode  MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
670 extern PetscErrorCode  MatFindZeroDiagonals(Mat,IS*);
671 
672 /*MC
673    MatSetValue - Set a single entry into a matrix.
674 
675    Not collective
676 
677    Input Parameters:
678 +  m - the matrix
679 .  row - the row location of the entry
680 .  col - the column location of the entry
681 .  value - the value to insert
682 -  mode - either INSERT_VALUES or ADD_VALUES
683 
684    Notes:
685    For efficiency one should use MatSetValues() and set several or many
686    values simultaneously if possible.
687 
688    Level: beginner
689 
690 .seealso: MatSetValues(), MatSetValueLocal()
691 M*/
692 PETSC_STATIC_INLINE PetscErrorCode MatSetValue(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValues(v,1,&i,1,&j,&va,mode);}
693 
694 PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
695 
696 PETSC_STATIC_INLINE PetscErrorCode MatSetValueLocal(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValuesLocal(v,1,&i,1,&j,&va,mode);}
697 
698 /*MC
699    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
700        row in a matrix providing the data that one can use to correctly preallocate the matrix.
701 
702    Synopsis:
703    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
704 
705    Collective on MPI_Comm
706 
707    Input Parameters:
708 +  comm - the communicator that will share the eventually allocated matrix
709 .  nrows - the number of LOCAL rows in the matrix
710 -  ncols - the number of LOCAL columns in the matrix
711 
712    Output Parameters:
713 +  dnz - the array that will be passed to the matrix preallocation routines
714 -  ozn - the other array passed to the matrix preallocation routines
715 
716 
717    Level: intermediate
718 
719    Notes:
720     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
721 
722    Do not malloc or free dnz and onz, that is handled internally by these routines
723 
724    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
725 
726    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
727 
728   Concepts: preallocation^Matrix
729 
730 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
731           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
732 M*/
733 #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
734 { \
735   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
736   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
737   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
738   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr); __start = 0; __end = __start; \
739   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
740   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
741 
742 /*MC
743    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
744        inserted using a local number of the rows and columns
745 
746    Synopsis:
747    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
748 
749    Not Collective
750 
751    Input Parameters:
752 +  map - the row mapping from local numbering to global numbering
753 .  nrows - the number of rows indicated
754 .  rows - the indices of the rows
755 .  cmap - the column mapping from local to global numbering
756 .  ncols - the number of columns in the matrix
757 .  cols - the columns indicated
758 .  dnz - the array that will be passed to the matrix preallocation routines
759 -  ozn - the other array passed to the matrix preallocation routines
760 
761 
762    Level: intermediate
763 
764    Notes:
765     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
766 
767    Do not malloc or free dnz and onz, that is handled internally by these routines
768 
769   Concepts: preallocation^Matrix
770 
771 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
772           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
773 M*/
774 #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
775 {\
776   PetscInt __l;\
777   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
778   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
779   for (__l=0;__l<nrows;__l++) {\
780     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
781   }\
782 }
783 
784 /*MC
785    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
786        inserted using a local number of the rows and columns
787 
788    Synopsis:
789    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
790 
791    Not Collective
792 
793    Input Parameters:
794 +  map - the mapping between local numbering and global numbering
795 .  nrows - the number of rows indicated
796 .  rows - the indices of the rows
797 .  ncols - the number of columns in the matrix
798 .  cols - the columns indicated
799 .  dnz - the array that will be passed to the matrix preallocation routines
800 -  ozn - the other array passed to the matrix preallocation routines
801 
802 
803    Level: intermediate
804 
805    Notes:
806     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
807 
808    Do not malloc or free dnz and onz that is handled internally by these routines
809 
810   Concepts: preallocation^Matrix
811 
812 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
813           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
814 M*/
815 #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
816 {\
817   PetscInt __l;\
818   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
819   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
820   for (__l=0;__l<nrows;__l++) {\
821     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
822   }\
823 }
824 
825 /*MC
826    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
827        inserted using a local number of the rows and columns
828 
829    Synopsis:
830    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
831 
832    Not Collective
833 
834    Input Parameters:
835 +  row - the row
836 .  ncols - the number of columns in the matrix
837 -  cols - the columns indicated
838 
839    Output Parameters:
840 +  dnz - the array that will be passed to the matrix preallocation routines
841 -  ozn - the other array passed to the matrix preallocation routines
842 
843 
844    Level: intermediate
845 
846    Notes:
847     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
848 
849    Do not malloc or free dnz and onz that is handled internally by these routines
850 
851    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
852 
853   Concepts: preallocation^Matrix
854 
855 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
856           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
857 M*/
858 #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
859 { PetscInt __i; \
860   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);\
861   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);\
862   for (__i=0; __i<nc; __i++) {\
863     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
864     else dnz[row - __rstart]++;\
865   }\
866 }
867 
868 /*MC
869    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
870        inserted using a local number of the rows and columns
871 
872    Synopsis:
873    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
874 
875    Not Collective
876 
877    Input Parameters:
878 +  nrows - the number of rows indicated
879 .  rows - the indices of the rows
880 .  ncols - the number of columns in the matrix
881 .  cols - the columns indicated
882 .  dnz - the array that will be passed to the matrix preallocation routines
883 -  ozn - the other array passed to the matrix preallocation routines
884 
885 
886    Level: intermediate
887 
888    Notes:
889     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
890 
891    Do not malloc or free dnz and onz that is handled internally by these routines
892 
893    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
894 
895   Concepts: preallocation^Matrix
896 
897 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
898           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
899 M*/
900 #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
901 { PetscInt __i; \
902   for (__i=0; __i<nc; __i++) {\
903     if (cols[__i] >= __end) onz[row - __rstart]++; \
904     else if (cols[__i] >= row) dnz[row - __rstart]++;\
905   }\
906 }
907 
908 /*MC
909    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
910 
911    Synopsis:
912    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
913 
914    Not Collective
915 
916    Input Parameters:
917 .  A - matrix
918 .  row - row where values exist (must be local to this process)
919 .  ncols - number of columns
920 .  cols - columns with nonzeros
921 .  dnz - the array that will be passed to the matrix preallocation routines
922 -  ozn - the other array passed to the matrix preallocation routines
923 
924 
925    Level: intermediate
926 
927    Notes:
928     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
929 
930    Do not malloc or free dnz and onz that is handled internally by these routines
931 
932    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
933 
934   Concepts: preallocation^Matrix
935 
936 .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
937           MatPreallocateSymmetricSetLocal()
938 M*/
939 #define MatPreallocateLocation(A,row,ncols,cols,dnz,onz) 0;if (A) {ierr = MatSetValues(A,1,&row,ncols,cols,PETSC_NULL,INSERT_VALUES);CHKERRQ(ierr);} else {ierr =  MatPreallocateSet(row,ncols,cols,dnz,onz);CHKERRQ(ierr);}
940 
941 
942 /*MC
943    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
944        row in a matrix providing the data that one can use to correctly preallocate the matrix.
945 
946    Synopsis:
947    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
948 
949    Collective on MPI_Comm
950 
951    Input Parameters:
952 +  dnz - the array that was be passed to the matrix preallocation routines
953 -  ozn - the other array passed to the matrix preallocation routines
954 
955 
956    Level: intermediate
957 
958    Notes:
959     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
960 
961    Do not malloc or free dnz and onz that is handled internally by these routines
962 
963    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
964 
965   Concepts: preallocation^Matrix
966 
967 .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
968           MatPreallocateSymmetricSetLocal()
969 M*/
970 #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
971 
972 
973 
974 /* Routines unique to particular data structures */
975 extern PetscErrorCode  MatShellGetContext(Mat,void *);
976 PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
977 
978 extern PetscErrorCode  MatInodeAdjustForInodes(Mat,IS*,IS*);
979 extern PetscErrorCode  MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
980 
981 extern PetscErrorCode  MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
982 extern PetscErrorCode  MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
983 extern PetscErrorCode  MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
984 extern PetscErrorCode  MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
985 extern PetscErrorCode  MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
986 extern PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
987 
988 #define MAT_SKIP_ALLOCATION -4
989 
990 extern PetscErrorCode  MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
991 PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
992 extern PetscErrorCode  MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
993 PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
994 extern PetscErrorCode  MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
995 PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
996 
997 extern PetscErrorCode  MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
998 PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
999 extern PetscErrorCode  MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1000 extern PetscErrorCode  MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1001 extern PetscErrorCode  MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
1002 extern PetscErrorCode  MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1003 extern PetscErrorCode  MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
1004 extern PetscErrorCode  MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1005 extern PetscErrorCode  MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
1006 extern PetscErrorCode  MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
1007 extern PetscErrorCode  MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
1008 extern PetscErrorCode  MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
1009 extern PetscErrorCode  MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
1010 extern PetscErrorCode  MatAdicSetLocalFunction(Mat,void (*)(void));
1011 extern PetscErrorCode  MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
1012 
1013 extern PetscErrorCode  MatSeqDenseSetLDA(Mat,PetscInt);
1014 extern PetscErrorCode  MatDenseGetLocalMatrix(Mat,Mat*);
1015 
1016 extern PetscErrorCode  MatStoreValues(Mat);
1017 extern PetscErrorCode  MatRetrieveValues(Mat);
1018 
1019 extern PetscErrorCode  MatDAADSetCtx(Mat,void*);
1020 
1021 extern PetscErrorCode  MatFindNonzeroRows(Mat,IS*);
1022 /*
1023   These routines are not usually accessed directly, rather solving is
1024   done through the KSP and PC interfaces.
1025 */
1026 
1027 /*J
1028     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
1029        with an optional dynamic library name, for example
1030        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
1031 
1032    Level: beginner
1033 
1034    Cannot use const because the PC objects manipulate the string
1035 
1036 .seealso: MatGetOrdering()
1037 J*/
1038 #define MatOrderingType char*
1039 #define MATORDERINGNATURAL     "natural"
1040 #define MATORDERINGND          "nd"
1041 #define MATORDERING1WD         "1wd"
1042 #define MATORDERINGRCM         "rcm"
1043 #define MATORDERINGQMD         "qmd"
1044 #define MATORDERINGROWLENGTH   "rowlength"
1045 #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
1046 
1047 extern PetscErrorCode  MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
1048 extern PetscErrorCode  MatGetOrderingList(PetscFList *list);
1049 extern PetscErrorCode  MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
1050 
1051 /*MC
1052    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package.
1053 
1054    Synopsis:
1055    PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
1056 
1057    Not Collective
1058 
1059    Input Parameters:
1060 +  sname - name of ordering (for example MATORDERINGND)
1061 .  path - location of library where creation routine is
1062 .  name - name of function that creates the ordering type,a string
1063 -  function - function pointer that creates the ordering
1064 
1065    Level: developer
1066 
1067    If dynamic libraries are used, then the fourth input argument (function)
1068    is ignored.
1069 
1070    Sample usage:
1071 .vb
1072    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
1073                "MyOrder",MyOrder);
1074 .ve
1075 
1076    Then, your partitioner can be chosen with the procedural interface via
1077 $     MatOrderingSetType(part,"my_order)
1078    or at runtime via the option
1079 $     -pc_factor_mat_ordering_type my_order
1080 
1081    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
1082 
1083 .keywords: matrix, ordering, register
1084 
1085 .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
1086 M*/
1087 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1088 #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1089 #else
1090 #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1091 #endif
1092 
1093 extern PetscErrorCode  MatOrderingRegisterDestroy(void);
1094 extern PetscErrorCode  MatOrderingRegisterAll(const char[]);
1095 extern PetscBool  MatOrderingRegisterAllCalled;
1096 extern PetscFList MatOrderingList;
1097 
1098 extern PetscErrorCode  MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1099 
1100 /*S
1101     MatFactorShiftType - Numeric Shift.
1102 
1103    Level: beginner
1104 
1105 S*/
1106 typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1107 extern const char *MatFactorShiftTypes[];
1108 
1109 /*S
1110    MatFactorInfo - Data passed into the matrix factorization routines
1111 
1112    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
1113 $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1114 
1115    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1116 
1117       You can use MatFactorInfoInitialize() to set default values.
1118 
1119    Level: developer
1120 
1121 .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1122           MatFactorInfoInitialize()
1123 
1124 S*/
1125 typedef struct {
1126   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
1127   PetscReal     usedt;
1128   PetscReal     dt;             /* drop tolerance */
1129   PetscReal     dtcol;          /* tolerance for pivoting */
1130   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
1131   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1132   PetscReal     levels;         /* ICC/ILU(levels) */
1133   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1134                                    factorization may be faster if do not pivot */
1135   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1136   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1137   PetscReal     shiftamount;     /* how large the shift is */
1138 } MatFactorInfo;
1139 
1140 extern PetscErrorCode  MatFactorInfoInitialize(MatFactorInfo*);
1141 extern PetscErrorCode  MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
1142 extern PetscErrorCode  MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1143 extern PetscErrorCode  MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
1144 extern PetscErrorCode  MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
1145 extern PetscErrorCode  MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
1146 extern PetscErrorCode  MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1147 extern PetscErrorCode  MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1148 extern PetscErrorCode  MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1149 extern PetscErrorCode  MatICCFactor(Mat,IS,const MatFactorInfo*);
1150 extern PetscErrorCode  MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1151 extern PetscErrorCode  MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1152 extern PetscErrorCode  MatSolve(Mat,Vec,Vec);
1153 extern PetscErrorCode  MatForwardSolve(Mat,Vec,Vec);
1154 extern PetscErrorCode  MatBackwardSolve(Mat,Vec,Vec);
1155 extern PetscErrorCode  MatSolveAdd(Mat,Vec,Vec,Vec);
1156 extern PetscErrorCode  MatSolveTranspose(Mat,Vec,Vec);
1157 extern PetscErrorCode  MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1158 extern PetscErrorCode  MatSolves(Mat,Vecs,Vecs);
1159 
1160 extern PetscErrorCode  MatSetUnfactored(Mat);
1161 
1162 /*E
1163     MatSORType - What type of (S)SOR to perform
1164 
1165     Level: beginner
1166 
1167    May be bitwise ORd together
1168 
1169    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1170 
1171    MatSORType may be bitwise ORd together, so do not change the numbers
1172 
1173 .seealso: MatSOR()
1174 E*/
1175 typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1176               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1177               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
1178               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1179 extern PetscErrorCode  MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1180 
1181 /*
1182     These routines are for efficiently computing Jacobians via finite differences.
1183 */
1184 
1185 /*J
1186     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1187        with an optional dynamic library name, for example
1188        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1189 
1190    Level: beginner
1191 
1192 .seealso: MatGetColoring()
1193 J*/
1194 #define MatColoringType char*
1195 #define MATCOLORINGNATURAL "natural"
1196 #define MATCOLORINGSL      "sl"
1197 #define MATCOLORINGLF      "lf"
1198 #define MATCOLORINGID      "id"
1199 
1200 extern PetscErrorCode  MatGetColoring(Mat,const MatColoringType,ISColoring*);
1201 extern PetscErrorCode  MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
1202 
1203 /*MC
1204    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
1205                                matrix package.
1206 
1207    Synopsis:
1208    PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring))
1209 
1210    Not Collective
1211 
1212    Input Parameters:
1213 +  sname - name of Coloring (for example MATCOLORINGSL)
1214 .  path - location of library where creation routine is
1215 .  name - name of function that creates the Coloring type, a string
1216 -  function - function pointer that creates the coloring
1217 
1218    Level: developer
1219 
1220    If dynamic libraries are used, then the fourth input argument (function)
1221    is ignored.
1222 
1223    Sample usage:
1224 .vb
1225    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
1226                "MyColor",MyColor);
1227 .ve
1228 
1229    Then, your partitioner can be chosen with the procedural interface via
1230 $     MatColoringSetType(part,"my_color")
1231    or at runtime via the option
1232 $     -mat_coloring_type my_color
1233 
1234    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
1235 
1236 .keywords: matrix, Coloring, register
1237 
1238 .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
1239 M*/
1240 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1241 #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1242 #else
1243 #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1244 #endif
1245 
1246 extern PetscBool  MatColoringRegisterAllCalled;
1247 
1248 extern PetscErrorCode  MatColoringRegisterAll(const char[]);
1249 extern PetscErrorCode  MatColoringRegisterDestroy(void);
1250 extern PetscErrorCode  MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1251 
1252 /*S
1253      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1254         and coloring
1255 
1256    Level: beginner
1257 
1258   Concepts: coloring, sparse Jacobian, finite differences
1259 
1260 .seealso:  MatFDColoringCreate()
1261 S*/
1262 typedef struct _p_MatFDColoring* MatFDColoring;
1263 
1264 extern PetscErrorCode  MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1265 extern PetscErrorCode  MatFDColoringDestroy(MatFDColoring*);
1266 extern PetscErrorCode  MatFDColoringView(MatFDColoring,PetscViewer);
1267 extern PetscErrorCode  MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1268 extern PetscErrorCode  MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1269 extern PetscErrorCode  MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1270 extern PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring);
1271 extern PetscErrorCode  MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1272 extern PetscErrorCode  MatFDColoringSetF(MatFDColoring,Vec);
1273 extern PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1274 
1275 /*S
1276      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1277 
1278    Level: beginner
1279 
1280   Concepts: coloring, sparse matrix product
1281 
1282 .seealso:  MatTransposeColoringCreate()
1283 S*/
1284 typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1285 
1286 extern PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1287 extern PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1288 extern PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1289 extern PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1290 
1291 /*
1292     These routines are for partitioning matrices: currently used only
1293   for adjacency matrix, MatCreateMPIAdj().
1294 */
1295 
1296 /*S
1297      MatPartitioning - Object for managing the partitioning of a matrix or graph
1298 
1299    Level: beginner
1300 
1301   Concepts: partitioning
1302 
1303 .seealso:  MatPartitioningCreate(), MatPartitioningType
1304 S*/
1305 typedef struct _p_MatPartitioning* MatPartitioning;
1306 
1307 /*J
1308     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1309        with an optional dynamic library name, for example
1310        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1311 
1312    Level: beginner
1313 dm
1314 .seealso: MatPartitioningCreate(), MatPartitioning
1315 J*/
1316 #define MatPartitioningType char*
1317 #define MATPARTITIONINGCURRENT  "current"
1318 #define MATPARTITIONINGSQUARE   "square"
1319 #define MATPARTITIONINGPARMETIS "parmetis"
1320 #define MATPARTITIONINGCHACO    "chaco"
1321 #define MATPARTITIONINGPARTY    "party"
1322 #define MATPARTITIONINGPTSCOTCH "ptscotch"
1323 
1324 
1325 extern PetscErrorCode  MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1326 extern PetscErrorCode  MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1327 extern PetscErrorCode  MatPartitioningSetNParts(MatPartitioning,PetscInt);
1328 extern PetscErrorCode  MatPartitioningSetAdjacency(MatPartitioning,Mat);
1329 extern PetscErrorCode  MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1330 extern PetscErrorCode  MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1331 extern PetscErrorCode  MatPartitioningApply(MatPartitioning,IS*);
1332 extern PetscErrorCode  MatPartitioningDestroy(MatPartitioning*);
1333 
1334 extern PetscErrorCode  MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
1335 
1336 /*MC
1337    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
1338    matrix package.
1339 
1340    Synopsis:
1341    PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
1342 
1343    Not Collective
1344 
1345    Input Parameters:
1346 +  sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis
1347 .  path - location of library where creation routine is
1348 .  name - name of function that creates the partitioning type, a string
1349 -  function - function pointer that creates the partitioning type
1350 
1351    Level: developer
1352 
1353    If dynamic libraries are used, then the fourth input argument (function)
1354    is ignored.
1355 
1356    Sample usage:
1357 .vb
1358    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
1359                "MyPartCreate",MyPartCreate);
1360 .ve
1361 
1362    Then, your partitioner can be chosen with the procedural interface via
1363 $     MatPartitioningSetType(part,"my_part")
1364    or at runtime via the option
1365 $     -mat_partitioning_type my_part
1366 
1367    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
1368 
1369 .keywords: matrix, partitioning, register
1370 
1371 .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
1372 M*/
1373 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1374 #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
1375 #else
1376 #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
1377 #endif
1378 
1379 extern PetscBool  MatPartitioningRegisterAllCalled;
1380 
1381 extern PetscErrorCode  MatPartitioningRegisterAll(const char[]);
1382 extern PetscErrorCode  MatPartitioningRegisterDestroy(void);
1383 
1384 extern PetscErrorCode  MatPartitioningView(MatPartitioning,PetscViewer);
1385 extern PetscErrorCode  MatPartitioningSetFromOptions(MatPartitioning);
1386 extern PetscErrorCode  MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1387 
1388 extern PetscErrorCode  MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1389 extern PetscErrorCode  MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
1390 
1391 typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
1392 extern const char *MPChacoGlobalTypes[];
1393 typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
1394 extern const char *MPChacoLocalTypes[];
1395 typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
1396 extern const char *MPChacoEigenTypes[];
1397 
1398 extern PetscErrorCode  MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1399 extern PetscErrorCode  MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1400 extern PetscErrorCode  MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1401 extern PetscErrorCode  MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1402 extern PetscErrorCode  MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1403 extern PetscErrorCode  MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1404 extern PetscErrorCode  MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1405 extern PetscErrorCode  MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1406 extern PetscErrorCode  MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1407 extern PetscErrorCode  MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1408 extern PetscErrorCode  MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
1409 
1410 #define MP_PARTY_OPT "opt"
1411 #define MP_PARTY_LIN "lin"
1412 #define MP_PARTY_SCA "sca"
1413 #define MP_PARTY_RAN "ran"
1414 #define MP_PARTY_GBF "gbf"
1415 #define MP_PARTY_GCF "gcf"
1416 #define MP_PARTY_BUB "bub"
1417 #define MP_PARTY_DEF "def"
1418 extern PetscErrorCode  MatPartitioningPartySetGlobal(MatPartitioning,const char*);
1419 #define MP_PARTY_HELPFUL_SETS "hs"
1420 #define MP_PARTY_KERNIGHAN_LIN "kl"
1421 #define MP_PARTY_NONE "no"
1422 extern PetscErrorCode  MatPartitioningPartySetLocal(MatPartitioning,const char*);
1423 extern PetscErrorCode  MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1424 extern PetscErrorCode  MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1425 extern PetscErrorCode  MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
1426 
1427 typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
1428 extern const char *MPPTScotchStrategyTypes[];
1429 
1430 extern PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1431 extern PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1432 extern PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1433 extern PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
1434 
1435 /*
1436     These routines are for coarsening matrices:
1437 */
1438 
1439 /*S
1440      MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix)
1441 
1442    Level: beginner
1443 
1444   Concepts: coarsen
1445 
1446 .seealso:  MatCoarsenCreate), MatCoarsenType
1447 S*/
1448 typedef struct _p_MatCoarsen* MatCoarsen;
1449 
1450 /*J
1451     MatCoarsenType - String with the name of a PETSc matrix coarsen or the creation function
1452        with an optional dynamic library name, for example
1453        http://www.mcs.anl.gov/petsc/lib.a:coarsencreate()
1454 
1455    Level: beginner
1456 dm
1457 .seealso: MatCoarsenCreate(), MatCoarsen
1458 J*/
1459 #define MatCoarsenType char*
1460 #define MATCOARSENMIS  "mis"
1461 #define MATCOARSENHEM  "hem"
1462 
1463 /* linked list for aggregates */
1464 typedef struct _PetscCDIntNd{
1465   struct _PetscCDIntNd *next;
1466   PetscInt   gid;
1467 }PetscCDIntNd;
1468 
1469 /* only used by node pool */
1470 typedef struct _PetscCDArrNd{
1471   struct _PetscCDArrNd *next;
1472   struct _PetscCDIntNd *array;
1473 }PetscCDArrNd;
1474 
1475 typedef struct _PetscCoarsenData{
1476   /* node pool */
1477   PetscCDArrNd  pool_list;
1478   PetscCDIntNd *new_node;
1479   PetscInt   new_left;
1480   PetscInt   chk_sz;
1481   PetscCDIntNd *extra_nodes;
1482   /* Array of lists */
1483   PetscCDIntNd **array;
1484   PetscInt size;
1485   /* cache a Mat for communication data */
1486   Mat mat;
1487   /* cache IS of removed equations */
1488   IS removedIS;
1489 }PetscCoarsenData;
1490 
1491 extern PetscErrorCode  MatCoarsenCreate(MPI_Comm,MatCoarsen*);
1492 extern PetscErrorCode  MatCoarsenSetType(MatCoarsen,const MatCoarsenType);
1493 extern PetscErrorCode  MatCoarsenSetAdjacency(MatCoarsen,Mat);
1494 extern PetscErrorCode  MatCoarsenSetGreedyOrdering(MatCoarsen,const IS);
1495 extern PetscErrorCode  MatCoarsenSetStrictAggs(MatCoarsen,PetscBool);
1496 extern PetscErrorCode  MatCoarsenSetVerbose(MatCoarsen,PetscInt);
1497 extern PetscErrorCode  MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** );
1498 extern PetscErrorCode  MatCoarsenApply(MatCoarsen);
1499 extern PetscErrorCode  MatCoarsenDestroy(MatCoarsen*);
1500 
1501 extern PetscErrorCode  MatCoarsenRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatCoarsen));
1502 
1503 /*MC
1504    MatCoarsenRegisterDynamic - Adds a new sparse matrix coarsen to the
1505    matrix package.
1506 
1507    Synopsis:
1508    PetscErrorCode MatCoarsenRegisterDynamic(const char *name_coarsen,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatCoarsen))
1509 
1510    Not Collective
1511 
1512    Input Parameters:
1513 +  sname - name of coarsen (for example MATCOARSENMIS)
1514 .  path - location of library where creation routine is
1515 .  name - name of function that creates the coarsen type, a string
1516 -  function - function pointer that creates the coarsen type
1517 
1518    Level: developer
1519 
1520    If dynamic libraries are used, then the fourth input argument (function)
1521    is ignored.
1522 
1523    Sample usage:
1524 .vb
1525    MatCoarsenRegisterDynamic("my_agg",/home/username/my_lib/lib/libO/solaris/mylib.a,
1526                "MyAggCreate",MyAggCreate);
1527 .ve
1528 
1529    Then, your aggregator can be chosen with the procedural interface via
1530 $     MatCoarsenSetType(agg,"my_agg")
1531    or at runtime via the option
1532 $     -mat_coarsen_type my_agg
1533 
1534    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
1535 
1536 .keywords: matrix, coarsen, register
1537 
1538 .seealso: MatCoarsenRegisterDestroy(), MatCoarsenRegisterAll()
1539 M*/
1540 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1541 #define MatCoarsenRegisterDynamic(a,b,c,d) MatCoarsenRegister(a,b,c,0)
1542 #else
1543 #define MatCoarsenRegisterDynamic(a,b,c,d) MatCoarsenRegister(a,b,c,d)
1544 #endif
1545 
1546 extern PetscBool  MatCoarsenRegisterAllCalled;
1547 
1548 extern PetscErrorCode  MatCoarsenRegisterAll(const char[]);
1549 extern PetscErrorCode  MatCoarsenRegisterDestroy(void);
1550 
1551 extern PetscErrorCode  MatCoarsenView(MatCoarsen,PetscViewer);
1552 extern PetscErrorCode  MatCoarsenSetFromOptions(MatCoarsen);
1553 extern PetscErrorCode  MatCoarsenGetType(MatCoarsen,const MatCoarsenType*);
1554 
1555 
1556 extern PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1557 extern PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1558 
1559 /*
1560     If you add entries here you must also add them to finclude/petscmat.h
1561 */
1562 typedef enum { MATOP_SET_VALUES=0,
1563                MATOP_GET_ROW=1,
1564                MATOP_RESTORE_ROW=2,
1565                MATOP_MULT=3,
1566                MATOP_MULT_ADD=4,
1567                MATOP_MULT_TRANSPOSE=5,
1568                MATOP_MULT_TRANSPOSE_ADD=6,
1569                MATOP_SOLVE=7,
1570                MATOP_SOLVE_ADD=8,
1571                MATOP_SOLVE_TRANSPOSE=9,
1572                MATOP_SOLVE_TRANSPOSE_ADD=10,
1573                MATOP_LUFACTOR=11,
1574                MATOP_CHOLESKYFACTOR=12,
1575                MATOP_SOR=13,
1576                MATOP_TRANSPOSE=14,
1577                MATOP_GETINFO=15,
1578                MATOP_EQUAL=16,
1579                MATOP_GET_DIAGONAL=17,
1580                MATOP_DIAGONAL_SCALE=18,
1581                MATOP_NORM=19,
1582                MATOP_ASSEMBLY_BEGIN=20,
1583                MATOP_ASSEMBLY_END=21,
1584                MATOP_SET_OPTION=22,
1585                MATOP_ZERO_ENTRIES=23,
1586                MATOP_ZERO_ROWS=24,
1587                MATOP_LUFACTOR_SYMBOLIC=25,
1588                MATOP_LUFACTOR_NUMERIC=26,
1589                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1590                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1591                MATOP_SETUP_PREALLOCATION=29,
1592                MATOP_ILUFACTOR_SYMBOLIC=30,
1593                MATOP_ICCFACTOR_SYMBOLIC=31,
1594                MATOP_GET_ARRAY=32,
1595                MATOP_RESTORE_ARRAY=33,
1596                MATOP_DUPLICATE=34,
1597                MATOP_FORWARD_SOLVE=35,
1598                MATOP_BACKWARD_SOLVE=36,
1599                MATOP_ILUFACTOR=37,
1600                MATOP_ICCFACTOR=38,
1601                MATOP_AXPY=39,
1602                MATOP_GET_SUBMATRICES=40,
1603                MATOP_INCREASE_OVERLAP=41,
1604                MATOP_GET_VALUES=42,
1605                MATOP_COPY=43,
1606                MATOP_GET_ROW_MAX=44,
1607                MATOP_SCALE=45,
1608                MATOP_SHIFT=46,
1609                MATOP_DIAGONAL_SET=47,
1610                MATOP_ILUDT_FACTOR=48,
1611                MATOP_SET_BLOCK_SIZE=49,
1612                MATOP_GET_ROW_IJ=50,
1613                MATOP_RESTORE_ROW_IJ=51,
1614                MATOP_GET_COLUMN_IJ=52,
1615                MATOP_RESTORE_COLUMN_IJ=53,
1616                MATOP_FDCOLORING_CREATE=54,
1617                MATOP_COLORING_PATCH=55,
1618                MATOP_SET_UNFACTORED=56,
1619                MATOP_PERMUTE=57,
1620                MATOP_SET_VALUES_BLOCKED=58,
1621                MATOP_GET_SUBMATRIX=59,
1622                MATOP_DESTROY=60,
1623                MATOP_VIEW=61,
1624                MATOP_CONVERT_FROM=62,
1625                MATOP_USE_SCALED_FORM=63,
1626                MATOP_SCALE_SYSTEM=64,
1627                MATOP_UNSCALE_SYSTEM=65,
1628                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1629                MATOP_SET_VALUES_LOCAL=67,
1630                MATOP_ZERO_ROWS_LOCAL=68,
1631                MATOP_GET_ROW_MAX_ABS=69,
1632                MATOP_GET_ROW_MIN_ABS=70,
1633                MATOP_CONVERT=71,
1634                MATOP_SET_COLORING=72,
1635                MATOP_SET_VALUES_ADIC=73,
1636                MATOP_SET_VALUES_ADIFOR=74,
1637                MATOP_FD_COLORING_APPLY=75,
1638                MATOP_SET_FROM_OPTIONS=76,
1639                MATOP_MULT_CON=77,
1640                MATOP_MULT_TRANSPOSE_CON=78,
1641                MATOP_PERMUTE_SPARSIFY=79,
1642                MATOP_MULT_MULTIPLE=80,
1643                MATOP_SOLVE_MULTIPLE=81,
1644                MATOP_GET_INERTIA=82,
1645                MATOP_LOAD=83,
1646                MATOP_IS_SYMMETRIC=84,
1647                MATOP_IS_HERMITIAN=85,
1648                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1649                MATOP_DUMMY=87,
1650                MATOP_GET_VECS=88,
1651                MATOP_MAT_MULT=89,
1652                MATOP_MAT_MULT_SYMBOLIC=90,
1653                MATOP_MAT_MULT_NUMERIC=91,
1654                MATOP_PTAP=92,
1655                MATOP_PTAP_SYMBOLIC=93,
1656                MATOP_PTAP_NUMERIC=94,
1657                MATOP_MAT_MULTTRANSPOSE=95,
1658                MATOP_MAT_MULTTRANSPOSE_SYM=96,
1659                MATOP_MAT_MULTTRANSPOSE_NUM=97,
1660                MATOP_PTAP_SYMBOLIC_SEQAIJ=98,
1661                MATOP_PTAP_NUMERIC_SEQAIJ=99,
1662                MATOP_PTAP_SYMBOLIC_MPIAIJ=100,
1663                MATOP_PTAP_NUMERIC_MPIAIJ=101,
1664                MATOP_CONJUGATE=102,
1665                MATOP_SET_SIZES=103,
1666                MATOP_SET_VALUES_ROW=104,
1667                MATOP_REAL_PART=105,
1668                MATOP_IMAG_PART=106,
1669                MATOP_GET_ROW_UTRIANGULAR=107,
1670                MATOP_RESTORE_ROW_UTRIANGULAR=108,
1671                MATOP_MATSOLVE=109,
1672                MATOP_GET_REDUNDANTMATRIX=110,
1673                MATOP_GET_ROW_MIN=111,
1674                MATOP_GET_COLUMN_VEC=112,
1675                MATOP_MISSING_DIAGONAL=113,
1676                MATOP_MATGETSEQNONZEROSTRUCTURE=114,
1677                MATOP_CREATE=115,
1678                MATOP_GET_GHOSTS=116,
1679                MATOP_GET_LOCALSUBMATRIX=117,
1680                MATOP_RESTORE_LOCALSUBMATRIX=118,
1681                MATOP_MULT_DIAGONAL_BLOCK=119,
1682                MATOP_HERMITIANTRANSPOSE=120,
1683                MATOP_MULTHERMITIANTRANSPOSE=121,
1684                MATOP_MULTHERMITIANTRANSPOSEADD=122,
1685                MATOP_GETMULTIPROCBLOCK=123,
1686                MATOP_GETCOLUMNNORMS=125,
1687 	       MATOP_GET_SUBMATRICES_PARALLEL=128,
1688                MATOP_SET_VALUES_BATCH=129,
1689                MATOP_TRANSPOSEMATMULT=130,
1690                MATOP_TRANSPOSEMATMULT_SYMBOLIC=131,
1691                MATOP_TRANSPOSEMATMULT_NUMERIC=132,
1692                MATOP_TRANSPOSECOLORING_CREATE=133,
1693                MATOP_TRANSCOLORING_APPLY_SPTODEN=134,
1694                MATOP_TRANSCOLORING_APPLY_DENTOSP=135,
1695                MATOP_RARt=136,
1696                MATOP_RARt_SYMBOLIC=137,
1697                MATOP_RARt_NUMERIC=138
1698              } MatOperation;
1699 extern PetscErrorCode  MatHasOperation(Mat,MatOperation,PetscBool *);
1700 extern PetscErrorCode  MatShellSetOperation(Mat,MatOperation,void(*)(void));
1701 extern PetscErrorCode  MatShellGetOperation(Mat,MatOperation,void(**)(void));
1702 extern PetscErrorCode  MatShellSetContext(Mat,void*);
1703 
1704 /*
1705    Codes for matrices stored on disk. By default they are
1706    stored in a universal format. By changing the format with
1707    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
1708    be stored in a way natural for the matrix, for example dense matrices
1709    would be stored as dense. Matrices stored this way may only be
1710    read into matrices of the same type.
1711 */
1712 #define MATRIX_BINARY_FORMAT_DENSE -1
1713 
1714 extern PetscErrorCode  MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1715 extern PetscErrorCode  MatISGetLocalMat(Mat,Mat*);
1716 extern PetscErrorCode  MatISSetLocalMat(Mat,Mat);
1717 
1718 /*S
1719      MatNullSpace - Object that removes a null space from a vector, i.e.
1720          orthogonalizes the vector to a subsapce
1721 
1722    Level: advanced
1723 
1724   Concepts: matrix; linear operator, null space
1725 
1726   Users manual sections:
1727 .   sec_singular
1728 
1729 .seealso:  MatNullSpaceCreate()
1730 S*/
1731 typedef struct _p_MatNullSpace* MatNullSpace;
1732 
1733 extern PetscErrorCode  MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1734 extern PetscErrorCode  MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1735 extern PetscErrorCode  MatNullSpaceDestroy(MatNullSpace*);
1736 extern PetscErrorCode  MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1737 extern PetscErrorCode  MatGetNullSpace(Mat, MatNullSpace *);
1738 extern PetscErrorCode  MatSetNullSpace(Mat,MatNullSpace);
1739 extern PetscErrorCode  MatSetNearNullSpace(Mat,MatNullSpace);
1740 extern PetscErrorCode  MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1741 extern PetscErrorCode  MatNullSpaceView(MatNullSpace,PetscViewer);
1742 
1743 extern PetscErrorCode  MatReorderingSeqSBAIJ(Mat,IS);
1744 extern PetscErrorCode  MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1745 extern PetscErrorCode  MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1746 extern PetscErrorCode  MatSeqBAIJInvertBlockDiagonal(Mat);
1747 
1748 extern PetscErrorCode  MatCreateMAIJ(Mat,PetscInt,Mat*);
1749 extern PetscErrorCode  MatMAIJRedimension(Mat,PetscInt,Mat*);
1750 extern PetscErrorCode  MatMAIJGetAIJ(Mat,Mat*);
1751 
1752 extern PetscErrorCode  MatComputeExplicitOperator(Mat,Mat*);
1753 
1754 extern PetscErrorCode  MatDiagonalScaleLocal(Mat,Vec);
1755 
1756 extern PetscErrorCode  MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1757 extern PetscErrorCode  MatMFFDSetBase(Mat,Vec,Vec);
1758 extern PetscErrorCode  MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1759 extern PetscErrorCode  MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1760 extern PetscErrorCode  MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1761 extern PetscErrorCode  MatMFFDAddNullSpace(Mat,MatNullSpace);
1762 extern PetscErrorCode  MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1763 extern PetscErrorCode  MatMFFDResetHHistory(Mat);
1764 extern PetscErrorCode  MatMFFDSetFunctionError(Mat,PetscReal);
1765 extern PetscErrorCode  MatMFFDSetPeriod(Mat,PetscInt);
1766 extern PetscErrorCode  MatMFFDGetH(Mat,PetscScalar *);
1767 extern PetscErrorCode  MatMFFDSetOptionsPrefix(Mat,const char[]);
1768 extern PetscErrorCode  MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1769 extern PetscErrorCode  MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1770 
1771 /*S
1772     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
1773               Jacobian vector products
1774 
1775     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
1776 
1777            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
1778 
1779     Level: developer
1780 
1781 .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
1782 S*/
1783 typedef struct _p_MatMFFD* MatMFFD;
1784 
1785 /*J
1786     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1787 
1788    Level: beginner
1789 
1790 .seealso: MatMFFDSetType(), MatMFFDRegister()
1791 J*/
1792 #define MatMFFDType char*
1793 #define MATMFFD_DS  "ds"
1794 #define MATMFFD_WP  "wp"
1795 
1796 extern PetscErrorCode  MatMFFDSetType(Mat,const MatMFFDType);
1797 extern PetscErrorCode  MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1798 
1799 /*MC
1800    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1801 
1802    Synopsis:
1803    PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1804 
1805    Not Collective
1806 
1807    Input Parameters:
1808 +  name_solver - name of a new user-defined compute-h module
1809 .  path - path (either absolute or relative) the library containing this solver
1810 .  name_create - name of routine to create method context
1811 -  routine_create - routine to create method context
1812 
1813    Level: developer
1814 
1815    Notes:
1816    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1817 
1818    If dynamic libraries are used, then the fourth input argument (routine_create)
1819    is ignored.
1820 
1821    Sample usage:
1822 .vb
1823    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1824                "MyHCreate",MyHCreate);
1825 .ve
1826 
1827    Then, your solver can be chosen with the procedural interface via
1828 $     MatMFFDSetType(mfctx,"my_h")
1829    or at runtime via the option
1830 $     -snes_mf_type my_h
1831 
1832 .keywords: MatMFFD, register
1833 
1834 .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1835 M*/
1836 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1837 #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1838 #else
1839 #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1840 #endif
1841 
1842 extern PetscErrorCode  MatMFFDRegisterAll(const char[]);
1843 extern PetscErrorCode  MatMFFDRegisterDestroy(void);
1844 extern PetscErrorCode  MatMFFDDSSetUmin(Mat,PetscReal);
1845 extern PetscErrorCode  MatMFFDWPSetComputeNormU(Mat,PetscBool );
1846 
1847 
1848 extern PetscErrorCode  PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1849 extern PetscErrorCode  PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
1850 
1851 /*
1852    PETSc interface to MUMPS
1853 */
1854 #ifdef PETSC_HAVE_MUMPS
1855 extern PetscErrorCode  MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1856 #endif
1857 
1858 /*
1859    PETSc interface to SUPERLU
1860 */
1861 #ifdef PETSC_HAVE_SUPERLU
1862 extern PetscErrorCode  MatSuperluSetILUDropTol(Mat,PetscReal);
1863 #endif
1864 
1865 #if defined(PETSC_HAVE_CUSP)
1866 extern PetscErrorCode  MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1867 extern PetscErrorCode  MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1868 #endif
1869 
1870 /*
1871    PETSc interface to FFTW
1872 */
1873 #if defined(PETSC_HAVE_FFTW)
1874 extern PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1875 extern PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
1876 extern PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*);
1877 #endif
1878 
1879 extern PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1880 extern PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1881 extern PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1882 extern PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1883 extern PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1884 extern PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
1885 extern PetscErrorCode MatNestSetVecType(Mat,const VecType);
1886 extern PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1887 extern PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1888 
1889 /*
1890  MatIJ:
1891  An unweighted directed pseudograph
1892  An interpretation of this matrix as a (pseudo)graph allows us to define additional operations on it:
1893  A MatIJ can act on sparse arrays: arrays of indices, or index arrays of integers, scalars, or integer-scalar pairs
1894  by mapping the indices to the indices connected to them by the (pseudo)graph ed
1895  */
1896 typedef enum {MATIJ_LOCAL, MATIJ_GLOBAL} MatIJIndexType;
1897 extern  PetscErrorCode MatIJSetMultivalued(Mat, PetscBool);
1898 extern  PetscErrorCode MatIJGetMultivalued(Mat, PetscBool*);
1899 extern  PetscErrorCode MatIJSetEdges(Mat, PetscInt, const PetscInt*, const PetscInt*);
1900 extern  PetscErrorCode MatIJGetEdges(Mat, PetscInt *, PetscInt **, PetscInt **);
1901 extern  PetscErrorCode MatIJSetEdgesIS(Mat, IS, IS);
1902 extern  PetscErrorCode MatIJGetEdgesIS(Mat, IS*, IS*);
1903 extern  PetscErrorCode MatIJGetRowSizes(Mat, MatIJIndexType, PetscInt, const PetscInt *, PetscInt **);
1904 extern  PetscErrorCode MatIJGetMinRowSize(Mat, PetscInt *);
1905 extern  PetscErrorCode MatIJGetMaxRowSize(Mat, PetscInt *);
1906 extern  PetscErrorCode MatIJGetSupport(Mat,  PetscInt *, PetscInt **);
1907 extern  PetscErrorCode MatIJGetSupportIS(Mat, IS *);
1908 extern  PetscErrorCode MatIJGetImage(Mat, PetscInt*, PetscInt**);
1909 extern  PetscErrorCode MatIJGetImageIS(Mat, IS *);
1910 extern  PetscErrorCode MatIJGetSupportSize(Mat, PetscInt *);
1911 extern  PetscErrorCode MatIJGetImageSize(Mat, PetscInt *);
1912 
1913 extern  PetscErrorCode MatIJBinRenumber(Mat, Mat*);
1914 
1915 extern  PetscErrorCode MatIJMap(Mat, MatIJIndexType, PetscInt,const PetscInt*,const PetscInt*,const PetscScalar*, MatIJIndexType,PetscInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt**);
1916 extern  PetscErrorCode MatIJBin(Mat, MatIJIndexType, PetscInt,const PetscInt*,const PetscInt*,const PetscScalar*,PetscInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt**);
1917 extern  PetscErrorCode MatIJBinMap(Mat,Mat, MatIJIndexType,PetscInt,const PetscInt*,const PetscInt*,const PetscScalar*,MatIJIndexType,PetscInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt**);
1918 
1919 PETSC_EXTERN_CXX_END
1920 #endif
1921