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