xref: /petsc/include/petscmat.h (revision 1d6018f00ceeb27ddc94c3818e3e9034da2fce00)
12eac72dbSBarry Smith /*
22eac72dbSBarry Smith      Include file for the matrix component of PETSc
32eac72dbSBarry Smith */
40a835dfdSSatish Balay #ifndef __PETSCMAT_H
50a835dfdSSatish Balay #define __PETSCMAT_H
60a835dfdSSatish Balay #include "petscvec.h"
7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN
82eac72dbSBarry Smith 
9d9274352SBarry Smith /*S
10d9274352SBarry Smith      Mat - Abstract PETSc matrix object
112eac72dbSBarry Smith 
12d91e6319SBarry Smith    Level: beginner
13d91e6319SBarry Smith 
14d9274352SBarry Smith   Concepts: matrix; linear operator
15d9274352SBarry Smith 
16d9274352SBarry Smith .seealso:  MatCreate(), MatType, MatSetType()
17d9274352SBarry Smith S*/
18d9274352SBarry Smith typedef struct _p_Mat*           Mat;
19d9274352SBarry Smith 
20d9274352SBarry Smith /*E
21d9274352SBarry Smith     MatType - String with the name of a PETSc matrix or the creation function
22d9274352SBarry Smith        with an optional dynamic library name, for example
23d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mymatcreate()
24d9274352SBarry Smith 
25d9274352SBarry Smith    Level: beginner
26d9274352SBarry Smith 
27c7393fdbSBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage
28d91e6319SBarry Smith E*/
29a313700dSBarry Smith #define MatType char*
30273d9f13SBarry Smith #define MATSAME            "same"
31273d9f13SBarry Smith #define MATSEQMAIJ         "seqmaij"
32273d9f13SBarry Smith #define MATMPIMAIJ         "mpimaij"
33209238afSKris Buschelman #define MATMAIJ            "maij"
34273d9f13SBarry Smith #define MATIS              "is"
35273d9f13SBarry Smith #define MATMPIROWBS        "mpirowbs"
36273d9f13SBarry Smith #define MATSEQAIJ          "seqaij"
37273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
38209238afSKris Buschelman #define MATAIJ             "aij"
39273d9f13SBarry Smith #define MATSHELL           "shell"
40209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
41273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
42209238afSKris Buschelman #define MATDENSE           "dense"
43273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
44273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
45209238afSKris Buschelman #define MATBAIJ            "baij"
46273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
47273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
48273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
49209238afSKris Buschelman #define MATSBAIJ           "sbaij"
50cebc7f6cSBarry Smith #define MATDAAD            "daad"
51cebc7f6cSBarry Smith #define MATMFFD            "mffd"
52c8a8475eSBarry Smith #define MATNORMAL          "normal"
53ab92ecdeSBarry Smith #define MATLRC             "lrc"
544b8d99adSRichard Tran Mills #define MATSEQCSRPERM      "seqcsrperm"
5518def556SRichard Tran Mills #define MATMPICSRPERM      "mpicsrperm"
56814655b8SBarry Smith #define MATCSRPERM         "csrperm"
5781824310SBarry Smith #define MATSEQCRL          "seqcrl"
58c02b24c5SBarry Smith #define MATMPICRL          "mpicrl"
59c02b24c5SBarry Smith #define MATCRL             "crl"
602a6744ebSBarry Smith #define MATSCATTER         "scatter"
61421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
62793850ffSBarry Smith #define MATCOMPOSITE       "composite"
635a7f1df3SHong Zhang #define MATSEQFFTW         "seqfftw"
64557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
6572ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
66*1d6018f0SLisandro Dalcin #define MATPYTHON          "python"
67d91e6319SBarry Smith 
689989ab13SBarry Smith /*E
69c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
709989ab13SBarry Smith 
719989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
729989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
739989ab13SBarry Smith 
749989ab13SBarry Smith 
759989ab13SBarry Smith    Level: beginner
769989ab13SBarry Smith 
775c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
789989ab13SBarry Smith E*/
79c7393fdbSBarry Smith #define MatSolverPackage char*
80b24902e0SBarry Smith #define MAT_SOLVER_SPOOLES      "spooles"
81b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU      "superlu"
82b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU_DIST "superlu_dist"
83b24902e0SBarry Smith #define MAT_SOLVER_UMFPACK      "umfpack"
84b24902e0SBarry Smith #define MAT_SOLVER_ESSL         "essl"
85b24902e0SBarry Smith #define MAT_SOLVER_LUSOL        "lusol"
86b24902e0SBarry Smith #define MAT_SOLVER_MUMPS        "mumps"
873bf14a46SMatthew Knepley #define MAT_SOLVER_PASTIX        "pastix"
88b24902e0SBarry Smith #define MAT_SOLVER_DSCPACK      "dscpack"
89b24902e0SBarry Smith #define MAT_SOLVER_MATLAB       "matlab"
9043244d56SBarry Smith #define MAT_SOLVER_PETSC        "petsc"
91b6806ab0SHong Zhang #define MAT_SOLVER_PLAPACK      "plapack"
92b24902e0SBarry Smith 
93b24902e0SBarry Smith /*E
94b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
95b24902e0SBarry Smith 
96b24902e0SBarry Smith     Level: beginner
97b24902e0SBarry Smith 
98b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
99b24902e0SBarry Smith 
100c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
101b24902e0SBarry Smith E*/
1025c9eb25fSBarry Smith typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC} MatFactorType;
103c7393fdbSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
104db4efbfdSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscTruth*);
10535bd34faSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
106b24902e0SBarry Smith 
1079989ab13SBarry Smith 
108c06d978dSMatthew Knepley /* Logging support */
109552e946dSBarry Smith #define    MAT_FILE_COOKIE 1211216    /* used to indicate matrices in binary files */
110be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE;
111be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
112be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE;
113be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE;
114e884886fSBarry Smith extern PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE;
115c06d978dSMatthew Knepley 
116ceb03754SKris Buschelman /*E
117ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
118ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
119ceb03754SKris Buschelman 
120ceb03754SKris Buschelman     Level: beginner
121ceb03754SKris Buschelman 
122ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
123ceb03754SKris Buschelman 
1240c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
125ceb03754SKris Buschelman E*/
126ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse;
127ceb03754SKris Buschelman 
1285494a064SHong Zhang /*E
1295494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
130829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1315494a064SHong Zhang 
1325494a064SHong Zhang     Level: beginner
1335494a064SHong Zhang 
134829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1355494a064SHong Zhang E*/
1365494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1375494a064SHong Zhang 
138e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]);
139c06d978dSMatthew Knepley 
140f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
141a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
142a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
143f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
144a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,const MatType);
145be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
146be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
147be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
148be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
14930de9b25SBarry Smith 
150f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
151f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
152f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
153f69a0ea3SMatthew Knepley 
15430de9b25SBarry Smith /*MC
15530de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
15630de9b25SBarry Smith 
15730de9b25SBarry Smith    Synopsis:
158c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
15930de9b25SBarry Smith 
16030de9b25SBarry Smith    Not Collective
16130de9b25SBarry Smith 
16230de9b25SBarry Smith    Input Parameters:
16330de9b25SBarry Smith +  name - name of a new user-defined matrix type
16430de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
16530de9b25SBarry Smith .  name_create - name of routine to create method context
16630de9b25SBarry Smith -  routine_create - routine to create method context
16730de9b25SBarry Smith 
16830de9b25SBarry Smith    Notes:
16930de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
17030de9b25SBarry Smith 
17130de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
17230de9b25SBarry Smith    is ignored.
17330de9b25SBarry Smith 
17430de9b25SBarry Smith    Sample usage:
17530de9b25SBarry Smith .vb
17630de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
17730de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
17830de9b25SBarry Smith .ve
17930de9b25SBarry Smith 
18030de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
18130de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
18230de9b25SBarry Smith    or at runtime via the option
18330de9b25SBarry Smith $     -mat_type my_mat
18430de9b25SBarry Smith 
18530de9b25SBarry Smith    Level: advanced
18630de9b25SBarry Smith 
187ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
18830de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
18930de9b25SBarry Smith 
19030de9b25SBarry Smith .keywords: Mat, register
19130de9b25SBarry Smith 
19230de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
19330de9b25SBarry Smith 
19430de9b25SBarry Smith M*/
195273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
196273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
197273d9f13SBarry Smith #else
198273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
19930de9b25SBarry Smith #endif
20030de9b25SBarry Smith 
201273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
202b0a32e0cSBarry Smith extern PetscFList MatList;
20328988994SBarry Smith 
2043b224e63SBarry Smith /*E
2053b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
2063b224e63SBarry Smith 
2073b224e63SBarry Smith     Level: beginner
2083b224e63SBarry Smith 
2093b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
2103b224e63SBarry Smith 
2113b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
2123b224e63SBarry Smith E*/
2133b224e63SBarry Smith typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
2143b224e63SBarry Smith 
215be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
216be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
217be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
218ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
219ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
220ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
221ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
222ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
223ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
224ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
225be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
226ba966639SSatish Balay 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)
227ba966639SSatish Balay 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)
228ba966639SSatish Balay 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)
229ba966639SSatish Balay 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)
230ba966639SSatish Balay 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))
231ba966639SSatish Balay 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))
232ba966639SSatish Balay 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))
233ba966639SSatish Balay 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)
234ba966639SSatish Balay 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)
235ba966639SSatish Balay 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)
236ba966639SSatish Balay 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)
237ba966639SSatish Balay 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))
238ba966639SSatish Balay 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))
239ba966639SSatish Balay 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))
24063ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2418d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2428d7a6e47SBarry Smith 
243be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
244be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
245ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,nz,nnz,&A),Mat,A)
246ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
247ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
248ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
249ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
250ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
251ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
252be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
253ba966639SSatish Balay 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)
254ba966639SSatish Balay 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)
255ba966639SSatish Balay 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)
256ba966639SSatish Balay 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)
257ba966639SSatish Balay 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))
258ba966639SSatish Balay 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))
259ba966639SSatish Balay 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))
260ba966639SSatish Balay 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)
261ba966639SSatish Balay 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)
262ba966639SSatish Balay 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)
263ba966639SSatish Balay 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)
264ba966639SSatish Balay 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))
265ba966639SSatish Balay 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))
266ba966639SSatish Balay 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))
267be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
269ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,nz,nnz,&A),Mat,A)
270ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
271ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
272ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
273ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
274ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
275ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
277ba966639SSatish Balay 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)
278ba966639SSatish Balay 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)
279ba966639SSatish Balay 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)
280ba966639SSatish Balay 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)
281ba966639SSatish Balay 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))
282ba966639SSatish Balay 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))
283ba966639SSatish Balay 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))
284ba966639SSatish Balay 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)
285ba966639SSatish Balay 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)
286ba966639SSatish Balay 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)
287ba966639SSatish Balay 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)
288ba966639SSatish Balay 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))
289ba966639SSatish Balay 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))
290ba966639SSatish Balay 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))
291be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
292ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(comm,m,n,M,N,ctx,&A),Mat,A)
293ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
294be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
295be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
296ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
297ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*);
298284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
2991472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
3001472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
3012a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*);
3022a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter);
3032a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*);
3048cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
305793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat);
306b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat);
307793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
3086d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
3096d7c1e57SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeSetType(Mat,MatCompositeType);
3106d7c1e57SBarry Smith 
3115a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*);
312f20085c4SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateTranspose(Mat,Mat*);
3131472f72bSBarry Smith 
314*1d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreatePython(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const char[],Mat*);
315*1d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPythonSetType(Mat,const char[]);
316*1d6018f0SLisandro Dalcin 
317*1d6018f0SLisandro Dalcin 
318f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
319be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
32021c89e3eSBarry Smith 
3210c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
32299cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat);
32399cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat);
32499cafbc1SBarry Smith 
3258ed539a5SBarry Smith /* ------------------------------------------------------------*/
326be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
327be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
32887d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
329f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
33084cb2905SBarry Smith 
3312ef4de8bSBarry Smith /*S
3322ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
3332ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3342ef4de8bSBarry Smith 
3352ef4de8bSBarry Smith    Level: beginner
3362ef4de8bSBarry Smith 
3372ef4de8bSBarry Smith   Concepts: matrix; linear operator
3382ef4de8bSBarry Smith 
339d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3402ef4de8bSBarry Smith S*/
341435da068SBarry Smith typedef struct {
342c1ac3661SBarry Smith   PetscInt k,j,i,c;
343435da068SBarry Smith } MatStencil;
3442ef4de8bSBarry Smith 
345be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
346be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
347be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
348435da068SBarry Smith 
349be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
350be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
351be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
3523a7fca6bSBarry Smith 
353d91e6319SBarry Smith /*E
354d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
355d91e6319SBarry Smith      to continue to add values to it
356d91e6319SBarry Smith 
357d91e6319SBarry Smith     Level: beginner
358d91e6319SBarry Smith 
359d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
360d91e6319SBarry Smith E*/
3616d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
362be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
363be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
364be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
3654f9c727eSBarry Smith 
3661d73ed98SBarry Smith 
36730de9b25SBarry Smith 
368d91e6319SBarry Smith /*E
369d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
370d91e6319SBarry Smith 
371d91e6319SBarry Smith     Level: beginner
372d91e6319SBarry Smith 
3730a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
374d91e6319SBarry Smith 
375d91e6319SBarry Smith .seealso: MatSetOption()
376d91e6319SBarry Smith E*/
377512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
3784e0d8c25SBarry Smith               MAT_SYMMETRIC,
3794e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
3808047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
3814e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
3824e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
3834e0d8c25SBarry Smith               MAT_KEEP_ZEROED_ROWS,MAT_IGNORE_ZERO_ENTRIES,
3844e0d8c25SBarry Smith               MAT_USE_INODES,
3854e0d8c25SBarry Smith               MAT_HERMITIAN,
3864e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
3874e0d8c25SBarry Smith               MAT_USE_COMPRESSEDROW,
3884e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
38928b2fa4aSMatthew Knepley               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
39028b2fa4aSMatthew Knepley               NUM_MAT_OPTIONS} MatOption;
391290bbb0aSBarry Smith extern const char *MatOptions[];
3924e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth);
393a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,const MatType*);
394a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
39584cb2905SBarry Smith 
396be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
398be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
399f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat);
400f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat);
401be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
402be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
403be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
404be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
405ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
406be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
407be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
408ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
409be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
4107b80b807SBarry Smith 
4111620fd73SBarry Smith 
412be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
413be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
414ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
415be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
416be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
417ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
418e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t)
4191cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*);
420be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
421ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
422be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
423be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
4242bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat);
425f5edf698SHong Zhang 
426d91e6319SBarry Smith /*E
427d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
428d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
429d91e6319SBarry Smith 
430d91e6319SBarry Smith     Level: beginner
431d91e6319SBarry Smith 
432d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
433d91e6319SBarry Smith 
434d91e6319SBarry Smith .seealso: MatDuplicate()
435d91e6319SBarry Smith E*/
4362e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
4372e8a6d31SBarry Smith 
438a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,const MatType,MatReuse,Mat*);
439a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
440be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
441ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
442ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
44394a9d846SBarry Smith 
444cb5b572fSBarry Smith 
445be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
446be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
447be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
448ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
449e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t)
450be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
451ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
4521cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*);
4531cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t)
454be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
455be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
45664c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *);
457a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,const MatType,Mat*);
458a313700dSBarry Smith PetscPolymorphicFunction(MatLoad,(PetscViewer v,const MatType t),(v,t,&a),Mat,a)
4597b80b807SBarry Smith 
4608f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4618f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
4628f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4638f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
464d4fbbf0eSBarry Smith 
465d91e6319SBarry Smith /*S
466d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
467d91e6319SBarry Smith 
468d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
469d91e6319SBarry Smith 
470d91e6319SBarry Smith    Level: intermediate
471d91e6319SBarry Smith 
472d91e6319SBarry Smith   Concepts: matrix^nonzero information
473d91e6319SBarry Smith 
474d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
475d91e6319SBarry Smith S*/
4764e220ebcSLois Curfman McInnes typedef struct {
477b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
478b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
479b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
480b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
481b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
482b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
483b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4844e220ebcSLois Curfman McInnes } MatInfo;
4854e220ebcSLois Curfman McInnes 
486d9274352SBarry Smith /*E
487d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
488d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
489d9274352SBarry Smith 
490d9274352SBarry Smith     Level: beginner
491d9274352SBarry Smith 
492d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
493d9274352SBarry Smith 
494d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
495d9274352SBarry Smith E*/
4967b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
497be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
499be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
500985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]);
501985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]);
502985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
503c87e5d42SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMinAbs(Mat,Vec,PetscInt[]);
50486697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec);
505fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*);
506fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
508ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
509be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
513ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
515be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
517be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
5187b80b807SBarry Smith 
519be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
520ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
521be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
522f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
523f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
524f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
525f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5267b80b807SBarry Smith 
527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
529be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5305ef9f2a5SBarry Smith 
531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5348ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**);
535f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
53672ca8751SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5377b80b807SBarry Smith 
538be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *);
541abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *);
542829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat *[]);
543829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat *[]);
5445494a064SHong Zhang 
545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
547be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
553dd6ea824SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,MatScalar**,Mat*);
55443eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
555cd116e26SSatish Balay #include "petscctable.h"
55643eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
55743eb5e2fSMatthew Knepley #else
55843eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
55943eb5e2fSMatthew Knepley #endif
5608efafbd8SBarry Smith 
561be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5627b80b807SBarry Smith 
563be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
565be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
56622440eb1SKris Buschelman 
567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
568be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
569be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
57022440eb1SKris Buschelman 
571be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
572be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
573be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
574bc011b1eSHong Zhang 
575f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
57604aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure);
577be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat);
5781c741599SHong Zhang 
579f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
580f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5817b80b807SBarry Smith 
582be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
583be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
584f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
585f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
586be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
587be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
588052efed2SBarry Smith 
589be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
590be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
59190f02eecSBarry Smith 
592be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
5930c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
594ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
595be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
596be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
59769db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
598bd481603SBarry Smith 
599bd481603SBarry Smith /*MC
6001620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
6011620fd73SBarry Smith 
6021620fd73SBarry Smith    Not collective
6031620fd73SBarry Smith 
6041620fd73SBarry Smith    Input Parameters:
6051620fd73SBarry Smith +  m - the matrix
6061620fd73SBarry Smith .  row - the row location of the entry
6071620fd73SBarry Smith .  col - the column location of the entry
6081620fd73SBarry Smith .  value - the value to insert
6091620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
6101620fd73SBarry Smith 
6111620fd73SBarry Smith    Notes:
6121620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6131620fd73SBarry Smith    values simultaneously if possible.
6141620fd73SBarry Smith 
6151620fd73SBarry Smith    Level: beginner
6161620fd73SBarry Smith 
6171620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6181620fd73SBarry Smith M*/
6191620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSetValue(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValues(v,1,&i,1,&j,&va,mode);}
6201620fd73SBarry Smith 
6211620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6221620fd73SBarry Smith 
6231620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSetValueLocal(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValuesLocal(v,1,&i,1,&j,&va,mode);}
6241620fd73SBarry Smith 
6251620fd73SBarry Smith /*MC
6260d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
627bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
628bd481603SBarry Smith 
629bd481603SBarry Smith    Synopsis:
630c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
631bd481603SBarry Smith 
632bd481603SBarry Smith    Collective on MPI_Comm
633bd481603SBarry Smith 
634bd481603SBarry Smith    Input Parameters:
635bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
636859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
637859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
638bd481603SBarry Smith 
639bd481603SBarry Smith    Output Parameters:
640bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
641bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
642bd481603SBarry Smith 
643bd481603SBarry Smith 
644bd481603SBarry Smith    Level: intermediate
645bd481603SBarry Smith 
646bd481603SBarry Smith    Notes:
647bd481603SBarry Smith    See the chapter in the users manual on performance for more details
648bd481603SBarry Smith 
6491d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
650bd481603SBarry Smith 
651bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
652bd481603SBarry Smith 
6531620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6541620fd73SBarry Smith 
655bd481603SBarry Smith   Concepts: preallocation^Matrix
656bd481603SBarry Smith 
657bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
658bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
659bd481603SBarry Smith M*/
660c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6617c922b88SBarry Smith { \
662a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
663a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
664a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
665a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
666c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
667a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6687c922b88SBarry Smith 
669bd481603SBarry Smith /*MC
6700d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
671bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
672bd481603SBarry Smith 
673bd481603SBarry Smith    Synopsis:
674c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
675bd481603SBarry Smith 
676bd481603SBarry Smith    Collective on MPI_Comm
677bd481603SBarry Smith 
678bd481603SBarry Smith    Input Parameters:
679bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
680859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
681859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
682bd481603SBarry Smith 
683bd481603SBarry Smith    Output Parameters:
684bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
685bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
686bd481603SBarry Smith 
687bd481603SBarry Smith 
688bd481603SBarry Smith    Level: intermediate
689bd481603SBarry Smith 
690bd481603SBarry Smith    Notes:
691bd481603SBarry Smith    See the chapter in the users manual on performance for more details
692bd481603SBarry Smith 
6931d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
694bd481603SBarry Smith 
6951620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6961620fd73SBarry Smith 
697bd481603SBarry Smith   Concepts: preallocation^Matrix
698bd481603SBarry Smith 
699bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
700bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
701bd481603SBarry Smith M*/
702222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
703222b16d4SBarry Smith { \
704a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
705a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
706a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
707a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
708c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
709a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
710222b16d4SBarry Smith 
711bd481603SBarry Smith /*MC
712bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
713bd481603SBarry Smith        inserted using a local number of the rows and columns
714bd481603SBarry Smith 
715bd481603SBarry Smith    Synopsis:
716c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
717bd481603SBarry Smith 
718bd481603SBarry Smith    Not Collective
719bd481603SBarry Smith 
720bd481603SBarry Smith    Input Parameters:
721bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
722bd481603SBarry Smith .  nrows - the number of rows indicated
7231d73ed98SBarry Smith .  rows - the indices of the rows
724bd481603SBarry Smith .  ncols - the number of columns in the matrix
725bd481603SBarry Smith .  cols - the columns indicated
726bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
727bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
728bd481603SBarry Smith 
729bd481603SBarry Smith 
730bd481603SBarry Smith    Level: intermediate
731bd481603SBarry Smith 
732bd481603SBarry Smith    Notes:
733bd481603SBarry Smith    See the chapter in the users manual on performance for more details
734bd481603SBarry Smith 
7351d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
736bd481603SBarry Smith 
737bd481603SBarry Smith   Concepts: preallocation^Matrix
738bd481603SBarry Smith 
739bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
740bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
741bd481603SBarry Smith M*/
742c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
743c4f061fbSSatish Balay {\
744c1ac3661SBarry Smith   PetscInt __l;\
745ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
746ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
747c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
748ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
749c4f061fbSSatish Balay   }\
750c4f061fbSSatish Balay }
751c4f061fbSSatish Balay 
752bd481603SBarry Smith /*MC
753bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
754bd481603SBarry Smith        inserted using a local number of the rows and columns
755bd481603SBarry Smith 
756bd481603SBarry Smith    Synopsis:
757c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
758bd481603SBarry Smith 
759bd481603SBarry Smith    Not Collective
760bd481603SBarry Smith 
761bd481603SBarry Smith    Input Parameters:
762bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
763bd481603SBarry Smith .  nrows - the number of rows indicated
7641d73ed98SBarry Smith .  rows - the indices of the rows
765bd481603SBarry Smith .  ncols - the number of columns in the matrix
766bd481603SBarry Smith .  cols - the columns indicated
767bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
768bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
769bd481603SBarry Smith 
770bd481603SBarry Smith 
771bd481603SBarry Smith    Level: intermediate
772bd481603SBarry Smith 
773bd481603SBarry Smith    Notes:
774bd481603SBarry Smith    See the chapter in the users manual on performance for more details
775bd481603SBarry Smith 
776bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
777bd481603SBarry Smith 
778bd481603SBarry Smith   Concepts: preallocation^Matrix
779bd481603SBarry Smith 
780bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
781bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
782bd481603SBarry Smith M*/
783d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
784d3d32019SBarry Smith {\
785c1ac3661SBarry Smith   PetscInt __l;\
786d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
787d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
788d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
789d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
790d3d32019SBarry Smith   }\
791d3d32019SBarry Smith }
792d3d32019SBarry Smith 
793bd481603SBarry Smith /*MC
794bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
795bd481603SBarry Smith        inserted using a local number of the rows and columns
796bd481603SBarry Smith 
797bd481603SBarry Smith    Synopsis:
798c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
799bd481603SBarry Smith 
800bd481603SBarry Smith    Not Collective
801bd481603SBarry Smith 
802bd481603SBarry Smith    Input Parameters:
80364075487SBarry Smith +  row - the row
804bd481603SBarry Smith .  ncols - the number of columns in the matrix
805a50f8bf6SHong Zhang -  cols - the columns indicated
806a50f8bf6SHong Zhang 
807a50f8bf6SHong Zhang    Output Parameters:
808a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
809bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
810bd481603SBarry Smith 
811bd481603SBarry Smith 
812bd481603SBarry Smith    Level: intermediate
813bd481603SBarry Smith 
814bd481603SBarry Smith    Notes:
815bd481603SBarry Smith    See the chapter in the users manual on performance for more details
816bd481603SBarry Smith 
817bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
818bd481603SBarry Smith 
8191620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8201620fd73SBarry Smith 
821bd481603SBarry Smith   Concepts: preallocation^Matrix
822bd481603SBarry Smith 
823bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
824bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
825bd481603SBarry Smith M*/
826c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
827c1ac3661SBarry Smith { PetscInt __i; \
8288ee2e534SBarry Smith   if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
829a89bc211SBarry Smith   if (row >= __rstart+__nrows) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D greater than last local row %D",row,__rstart+__nrows-1);\
8307c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
83164075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8327cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8337c922b88SBarry Smith   }\
8347c922b88SBarry Smith }
8357c922b88SBarry Smith 
836bd481603SBarry Smith /*MC
837bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
838bd481603SBarry Smith        inserted using a local number of the rows and columns
839bd481603SBarry Smith 
840bd481603SBarry Smith    Synopsis:
841c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
842bd481603SBarry Smith 
843bd481603SBarry Smith    Not Collective
844bd481603SBarry Smith 
845bd481603SBarry Smith    Input Parameters:
846bd481603SBarry Smith +  nrows - the number of rows indicated
8471d73ed98SBarry Smith .  rows - the indices of the rows
848bd481603SBarry Smith .  ncols - the number of columns in the matrix
849bd481603SBarry Smith .  cols - the columns indicated
850bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
851bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
852bd481603SBarry Smith 
853bd481603SBarry Smith 
854bd481603SBarry Smith    Level: intermediate
855bd481603SBarry Smith 
856bd481603SBarry Smith    Notes:
857bd481603SBarry Smith    See the chapter in the users manual on performance for more details
858bd481603SBarry Smith 
859bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
860bd481603SBarry Smith 
8611620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8621620fd73SBarry Smith 
863bd481603SBarry Smith   Concepts: preallocation^Matrix
864bd481603SBarry Smith 
865bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
866bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
867bd481603SBarry Smith M*/
868d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
869c1ac3661SBarry Smith { PetscInt __i; \
870d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
871d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
872d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
873d3d32019SBarry Smith   }\
874d3d32019SBarry Smith }
875d3d32019SBarry Smith 
876bd481603SBarry Smith /*MC
87716371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
87816371a99SBarry Smith 
87916371a99SBarry Smith    Synopsis:
88016371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
88116371a99SBarry Smith 
88216371a99SBarry Smith    Not Collective
88316371a99SBarry Smith 
88416371a99SBarry Smith    Input Parameters:
88516371a99SBarry Smith .  A - matrix
88616371a99SBarry Smith .  row - row where values exist (must be local to this process)
88716371a99SBarry Smith .  ncols - number of columns
88816371a99SBarry Smith .  cols - columns with nonzeros
88916371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
89016371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
89116371a99SBarry Smith 
89216371a99SBarry Smith 
89316371a99SBarry Smith    Level: intermediate
89416371a99SBarry Smith 
89516371a99SBarry Smith    Notes:
89616371a99SBarry Smith    See the chapter in the users manual on performance for more details
89716371a99SBarry Smith 
89816371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
89916371a99SBarry Smith 
90016371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
90116371a99SBarry Smith 
90216371a99SBarry Smith   Concepts: preallocation^Matrix
90316371a99SBarry Smith 
90416371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
90516371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
90616371a99SBarry Smith M*/
90716371a99SBarry Smith #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);}
90816371a99SBarry Smith 
90916371a99SBarry Smith 
91016371a99SBarry Smith /*MC
9110d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
912bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
913bd481603SBarry Smith 
914bd481603SBarry Smith    Synopsis:
915c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
916bd481603SBarry Smith 
917bd481603SBarry Smith    Collective on MPI_Comm
918bd481603SBarry Smith 
919bd481603SBarry Smith    Input Parameters:
92016371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
921bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
922bd481603SBarry Smith 
923bd481603SBarry Smith 
924bd481603SBarry Smith    Level: intermediate
925bd481603SBarry Smith 
926bd481603SBarry Smith    Notes:
927bd481603SBarry Smith    See the chapter in the users manual on performance for more details
928bd481603SBarry Smith 
929bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
930bd481603SBarry Smith 
9311620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9321620fd73SBarry Smith 
933bd481603SBarry Smith   Concepts: preallocation^Matrix
934bd481603SBarry Smith 
935bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
936bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
937bd481603SBarry Smith M*/
938a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9397c922b88SBarry Smith 
940bd481603SBarry Smith 
941bd481603SBarry Smith 
9427b80b807SBarry Smith /* Routines unique to particular data structures */
943be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
944ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
945698d4c6aSKris Buschelman 
946be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
947be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
948698d4c6aSKris Buschelman 
949be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
950be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
951be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
952c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
953c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9547b80b807SBarry Smith 
955a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
956a96a251dSBarry Smith 
957be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
958ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
959be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
960ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
961be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
962ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
963273d9f13SBarry Smith 
964be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
965ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
966be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
967be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
96853803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
969725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
970be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
971be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
972be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
973be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
974284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
975be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
976be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
977be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
978be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
979273d9f13SBarry Smith 
980be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
981ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*);
9821b807ce4Svictorle 
983be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
984be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
9852e8a6d31SBarry Smith 
986be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
9873a7fca6bSBarry Smith 
9887b80b807SBarry Smith /*
9897b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
99094b7f48cSBarry Smith   done through the KSP and PC interfaces.
9917b80b807SBarry Smith */
9927b80b807SBarry Smith 
993d9274352SBarry Smith /*E
994d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
995d9274352SBarry Smith        with an optional dynamic library name, for example
996d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
997d9274352SBarry Smith 
998d9274352SBarry Smith    Level: beginner
999d9274352SBarry Smith 
1000e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
1001e5a9bf91SBarry Smith 
1002d9274352SBarry Smith .seealso: MatGetOrdering()
1003d9274352SBarry Smith E*/
10043eaad2caSSatish Balay #define MatOrderingType char*
1005b12f92e5SBarry Smith #define MATORDERING_NATURAL   "natural"
1006b12f92e5SBarry Smith #define MATORDERING_ND        "nd"
1007b12f92e5SBarry Smith #define MATORDERING_1WD       "1wd"
1008b12f92e5SBarry Smith #define MATORDERING_RCM       "rcm"
1009b12f92e5SBarry Smith #define MATORDERING_QMD       "qmd"
1010b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength"
101162152c8bSBarry Smith #define MATORDERING_DSC_ND    "dsc_nd"
101262152c8bSBarry Smith #define MATORDERING_DSC_MMD   "dsc_mmd"
101362152c8bSBarry Smith #define MATORDERING_DSC_MDF   "dsc_mdf"
1014c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
1015c06d978dSMatthew Knepley #define MATORDERING_IDENTITY  "identity"
1016c06d978dSMatthew Knepley #define MATORDERING_REVERSE   "reverse"
1017b12f92e5SBarry Smith 
1018f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
1019be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
1020f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
102130de9b25SBarry Smith 
102230de9b25SBarry Smith /*MC
102330de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
102430de9b25SBarry Smith                                matrix package.
102530de9b25SBarry Smith 
102630de9b25SBarry Smith    Synopsis:
1027c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
102830de9b25SBarry Smith 
102930de9b25SBarry Smith    Not Collective
103030de9b25SBarry Smith 
103130de9b25SBarry Smith    Input Parameters:
103230de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
103330de9b25SBarry Smith .  path - location of library where creation routine is
103430de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
103530de9b25SBarry Smith -  function - function pointer that creates the ordering
103630de9b25SBarry Smith 
103730de9b25SBarry Smith    Level: developer
103830de9b25SBarry Smith 
103930de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
104030de9b25SBarry Smith    is ignored.
104130de9b25SBarry Smith 
104230de9b25SBarry Smith    Sample usage:
104330de9b25SBarry Smith .vb
104430de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
104530de9b25SBarry Smith                "MyOrder",MyOrder);
104630de9b25SBarry Smith .ve
104730de9b25SBarry Smith 
104830de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
104930de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
105030de9b25SBarry Smith    or at runtime via the option
10512401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
105230de9b25SBarry Smith 
1053ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
105430de9b25SBarry Smith 
105530de9b25SBarry Smith .keywords: matrix, ordering, register
105630de9b25SBarry Smith 
105730de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
105830de9b25SBarry Smith M*/
1059aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1060f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1061b12f92e5SBarry Smith #else
1062f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1063b12f92e5SBarry Smith #endif
106430de9b25SBarry Smith 
1065be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
1066be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
10672bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
1068b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1069d4fbbf0eSBarry Smith 
1070be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1071a2ce50c7SBarry Smith 
1072d91e6319SBarry Smith /*S
10732401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1074b00f7748SHong Zhang 
107561cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
107661cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1077b00f7748SHong Zhang 
107815e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1079b00f7748SHong Zhang 
108061cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
108161cad860SBarry Smith 
1082b00f7748SHong Zhang    Level: developer
1083b00f7748SHong Zhang 
1084d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1085d7d82daaSBarry Smith           MatFactorInfoInitialize()
1086b00f7748SHong Zhang 
1087b00f7748SHong Zhang S*/
1088b00f7748SHong Zhang typedef struct {
10890a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
1090fbf22428SSatish Balay   PetscReal     shiftpd;        /* if true, shift until positive pivots */
109115e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
109285317021SBarry Smith   PetscReal     usedt;
109315e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1094b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
109515e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
109667e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1097348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1098bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1099bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
110062bba022SBarry Smith   PetscReal     shiftinblocks;  /* if block in block factorization has zero pivot then shift diagonal until non-singular */
110115e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
110215e8a5b3SHong Zhang } MatFactorInfo;
1103ffa6d0a5SLois Curfman McInnes 
1104be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
11050481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
11060481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11070481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
11080481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
11090481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
11100481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11110481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11120481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11130481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,const MatFactorInfo*);
11140481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
11150481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,const MatFactorInfo*,Mat *);
1116be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1117be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1118be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1119be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1120be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1121be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1122be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1123be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
11248ed539a5SBarry Smith 
1125be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1126bb5a7306SBarry Smith 
1127d91e6319SBarry Smith /*E
1128d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1129bb1eb677SSatish Balay 
1130d91e6319SBarry Smith     Level: beginner
1131d91e6319SBarry Smith 
1132d9274352SBarry Smith    May be bitwise ORd together
1133d9274352SBarry Smith 
1134d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1135d91e6319SBarry Smith 
11364e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11374e7234bfSBarry Smith 
1138a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax()
1139d91e6319SBarry Smith E*/
1140ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1141ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1142ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
114384cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1144be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1145be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11468ed539a5SBarry Smith 
1147d4fbbf0eSBarry Smith /*
1148639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1149639f9d9dSBarry Smith */
1150b12f92e5SBarry Smith 
1151d9274352SBarry Smith /*E
1152d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1153d9274352SBarry Smith        with an optional dynamic library name, for example
1154d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1155d9274352SBarry Smith 
1156d9274352SBarry Smith    Level: beginner
1157d9274352SBarry Smith 
1158d9274352SBarry Smith .seealso: MatGetColoring()
1159d9274352SBarry Smith E*/
1160a313700dSBarry Smith #define MatColoringType char*
1161b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1162b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1163b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1164b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1165b12f92e5SBarry Smith 
1166ab970d9bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*);
11672e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
116830de9b25SBarry Smith 
116930de9b25SBarry Smith /*MC
117030de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
117130de9b25SBarry Smith                                matrix package.
117230de9b25SBarry Smith 
117330de9b25SBarry Smith    Synopsis:
1174c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
117530de9b25SBarry Smith 
117630de9b25SBarry Smith    Not Collective
117730de9b25SBarry Smith 
117830de9b25SBarry Smith    Input Parameters:
117930de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
118030de9b25SBarry Smith .  path - location of library where creation routine is
118130de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
118230de9b25SBarry Smith -  function - function pointer that creates the coloring
118330de9b25SBarry Smith 
118430de9b25SBarry Smith    Level: developer
118530de9b25SBarry Smith 
118630de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
118730de9b25SBarry Smith    is ignored.
118830de9b25SBarry Smith 
118930de9b25SBarry Smith    Sample usage:
119030de9b25SBarry Smith .vb
119130de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
119230de9b25SBarry Smith                "MyColor",MyColor);
119330de9b25SBarry Smith .ve
119430de9b25SBarry Smith 
119530de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
119630de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
119730de9b25SBarry Smith    or at runtime via the option
119830de9b25SBarry Smith $     -mat_coloring_type my_color
119930de9b25SBarry Smith 
1200ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
120130de9b25SBarry Smith 
120230de9b25SBarry Smith .keywords: matrix, Coloring, register
120330de9b25SBarry Smith 
120430de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
120530de9b25SBarry Smith M*/
1206aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1207f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1208b12f92e5SBarry Smith #else
1209f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1210b12f92e5SBarry Smith #endif
121130de9b25SBarry Smith 
12122bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1213f1144a30SSatish Balay 
1214be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1215be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1216be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1217639f9d9dSBarry Smith 
1218d9274352SBarry Smith /*S
1219d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1220d9274352SBarry Smith         and coloring
1221639f9d9dSBarry Smith 
1222d9274352SBarry Smith    Level: beginner
1223d9274352SBarry Smith 
1224d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1225d9274352SBarry Smith 
1226d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1227d9274352SBarry Smith S*/
1228e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1229639f9d9dSBarry Smith 
1230be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1231be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1232be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1233be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
12344a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1235be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1236be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1237be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1238be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1239be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1240be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1241639f9d9dSBarry Smith /*
12420752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12433eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12440752156aSBarry Smith */
1245ca161407SBarry Smith 
1246d9274352SBarry Smith /*S
1247d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1248d9274352SBarry Smith 
1249d9274352SBarry Smith    Level: beginner
1250d9274352SBarry Smith 
1251d9274352SBarry Smith   Concepts: partitioning
1252d9274352SBarry Smith 
1253743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1254d9274352SBarry Smith S*/
125591e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1256d9274352SBarry Smith 
1257d9274352SBarry Smith /*E
12585bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1259d9274352SBarry Smith        with an optional dynamic library name, for example
1260d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1261d9274352SBarry Smith 
1262d9274352SBarry Smith    Level: beginner
1263d9274352SBarry Smith 
1264b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1265d9274352SBarry Smith E*/
1266a313700dSBarry Smith #define MatPartitioningType char*
12678ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1268c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
12698ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
127071306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
127171306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
127271306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
127371306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
127471306c60Sjroman 
1275ca161407SBarry Smith 
1276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1277a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1278be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1279be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1280be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1281be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1282be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1283be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
12842aabb6bbSBarry Smith 
1285be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
128630de9b25SBarry Smith 
128730de9b25SBarry Smith /*MC
128830de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
128930de9b25SBarry Smith    matrix package.
129030de9b25SBarry Smith 
129130de9b25SBarry Smith    Synopsis:
1292c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
129330de9b25SBarry Smith 
129430de9b25SBarry Smith    Not Collective
129530de9b25SBarry Smith 
129630de9b25SBarry Smith    Input Parameters:
129730de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
129830de9b25SBarry Smith .  path - location of library where creation routine is
129930de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
130030de9b25SBarry Smith -  function - function pointer that creates the partitioning type
130130de9b25SBarry Smith 
130230de9b25SBarry Smith    Level: developer
130330de9b25SBarry Smith 
130430de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
130530de9b25SBarry Smith    is ignored.
130630de9b25SBarry Smith 
130730de9b25SBarry Smith    Sample usage:
130830de9b25SBarry Smith .vb
130930de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
131030de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
131130de9b25SBarry Smith .ve
131230de9b25SBarry Smith 
131330de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
131430de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
131530de9b25SBarry Smith    or at runtime via the option
131630de9b25SBarry Smith $     -mat_partitioning_type my_part
131730de9b25SBarry Smith 
1318ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
131930de9b25SBarry Smith 
132030de9b25SBarry Smith .keywords: matrix, partitioning, register
132130de9b25SBarry Smith 
132230de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
132330de9b25SBarry Smith M*/
1324aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1325f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13262aabb6bbSBarry Smith #else
1327f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13282aabb6bbSBarry Smith #endif
13292aabb6bbSBarry Smith 
13302bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1331f1144a30SSatish Balay 
1332be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1333be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
13342bad1931SBarry Smith 
1335be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1336be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1337a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1338ca161407SBarry Smith 
1339be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
13400e8af3caSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
13410752156aSBarry Smith 
1342be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1343be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
134471306c60Sjroman 
1345954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
1346be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
134771306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1348be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1349be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
135071306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1351be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1352be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1353be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
135471306c60Sjroman 
135571306c60Sjroman #define MP_PARTY_OPT "opt"
135671306c60Sjroman #define MP_PARTY_LIN "lin"
135771306c60Sjroman #define MP_PARTY_SCA "sca"
135871306c60Sjroman #define MP_PARTY_RAN "ran"
135971306c60Sjroman #define MP_PARTY_GBF "gbf"
136071306c60Sjroman #define MP_PARTY_GCF "gcf"
136171306c60Sjroman #define MP_PARTY_BUB "bub"
136271306c60Sjroman #define MP_PARTY_DEF "def"
1363be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
136471306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
136571306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
136671306c60Sjroman #define MP_PARTY_NONE "no"
1367be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1368be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1370be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
137171306c60Sjroman 
137271306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1373be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1374be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1375be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1376be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
137871306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1379be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
138271306c60Sjroman 
13830752156aSBarry Smith /*
13840a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1385d4fbbf0eSBarry Smith */
13861c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13871c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13881c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13891c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13901c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13917c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13927c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13931c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13941c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13957c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13967c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13971c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13981c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
13991c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
14001c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
14011c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
14021c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
14031c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
14041c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
14051c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
14061c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
14071c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
14081c1c02c0SLois Curfman McInnes                MATOP_COMPRESS=22,
14091c1c02c0SLois Curfman McInnes                MATOP_SET_OPTION=23,
14101c1c02c0SLois Curfman McInnes                MATOP_ZERO_ENTRIES=24,
14111c1c02c0SLois Curfman McInnes                MATOP_ZERO_ROWS=25,
14121c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_SYMBOLIC=26,
14131c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_NUMERIC=27,
14141c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_SYMBOLIC=28,
14151c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_NUMERIC=29,
1416d643ce63SMatthew Knepley                MATOP_SETUP_PREALLOCATION=30,
1417d643ce63SMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=31,
1418d643ce63SMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=32,
1419d643ce63SMatthew Knepley                MATOP_GET_ARRAY=33,
1420d643ce63SMatthew Knepley                MATOP_RESTORE_ARRAY=34,
1421e26f1c3cSBarry Smith                MATOP_DUPLICATE=35,
1422d643ce63SMatthew Knepley                MATOP_FORWARD_SOLVE=36,
1423d643ce63SMatthew Knepley                MATOP_BACKWARD_SOLVE=37,
1424d643ce63SMatthew Knepley                MATOP_ILUFACTOR=38,
1425d643ce63SMatthew Knepley                MATOP_ICCFACTOR=39,
1426d643ce63SMatthew Knepley                MATOP_AXPY=40,
1427d643ce63SMatthew Knepley                MATOP_GET_SUBMATRICES=41,
1428d643ce63SMatthew Knepley                MATOP_INCREASE_OVERLAP=42,
1429d643ce63SMatthew Knepley                MATOP_GET_VALUES=43,
1430d643ce63SMatthew Knepley                MATOP_COPY=44,
1431ff1cced0SMatthew Knepley                MATOP_GET_ROW_MAX=45,
1432d643ce63SMatthew Knepley                MATOP_SCALE=46,
1433d643ce63SMatthew Knepley                MATOP_SHIFT=47,
1434d643ce63SMatthew Knepley                MATOP_DIAGONAL_SHIFT=48,
1435d643ce63SMatthew Knepley                MATOP_ILUDT_FACTOR=49,
1436ff1cced0SMatthew Knepley                MATOP_SET_BLOCK_SIZE=50,
1437d643ce63SMatthew Knepley                MATOP_GET_ROW_IJ=51,
1438d643ce63SMatthew Knepley                MATOP_RESTORE_ROW_IJ=52,
1439d643ce63SMatthew Knepley                MATOP_GET_COLUMN_IJ=53,
1440d643ce63SMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=54,
1441d643ce63SMatthew Knepley                MATOP_FDCOLORING_CREATE=55,
1442d643ce63SMatthew Knepley                MATOP_COLORING_PATCH=56,
1443d643ce63SMatthew Knepley                MATOP_SET_UNFACTORED=57,
1444d643ce63SMatthew Knepley                MATOP_PERMUTE=58,
1445d643ce63SMatthew Knepley                MATOP_SET_VALUES_BLOCKED=59,
1446d643ce63SMatthew Knepley                MATOP_GET_SUBMATRIX=60,
1447d643ce63SMatthew Knepley                MATOP_DESTROY=61,
1448d643ce63SMatthew Knepley                MATOP_VIEW=62,
1449ff1cced0SMatthew Knepley                MATOP_CONVERT_FROM=63,
1450d643ce63SMatthew Knepley                MATOP_USE_SCALED_FORM=64,
1451d643ce63SMatthew Knepley                MATOP_SCALE_SYSTEM=65,
1452d643ce63SMatthew Knepley                MATOP_UNSCALE_SYSTEM=66,
1453ba16a8cbSKris Buschelman                MATOP_SET_LOCAL_TO_GLOBAL_MAP=67,
1454d643ce63SMatthew Knepley                MATOP_SET_VALUES_LOCAL=68,
1455d643ce63SMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=69,
1456ff1cced0SMatthew Knepley                MATOP_GET_ROW_MAX_ABS=70,
1457c87e5d42SMatthew Knepley                MATOP_GET_ROW_MIN_ABS=71,
1458c87e5d42SMatthew Knepley                MATOP_CONVERT=72,
1459c87e5d42SMatthew Knepley                MATOP_SET_COLORING=73,
1460c87e5d42SMatthew Knepley                MATOP_SET_VALUES_ADIC=74,
1461c87e5d42SMatthew Knepley                MATOP_SET_VALUES_ADIFOR=75,
1462c87e5d42SMatthew Knepley                MATOP_FD_COLORING_APPLY=76,
1463c87e5d42SMatthew Knepley                MATOP_SET_FROM_OPTIONS=77,
1464c87e5d42SMatthew Knepley                MATOP_MULT_CON=78,
1465c87e5d42SMatthew Knepley                MATOP_MULT_TRANSPOSE_CON=79,
1466d643ce63SMatthew Knepley                MATOP_PERMUTE_SPARSIFY=80,
1467d643ce63SMatthew Knepley                MATOP_MULT_MULTIPLE=81,
146841acf15aSKris Buschelman                MATOP_SOLVE_MULTIPLE=82,
146941acf15aSKris Buschelman                MATOP_GET_INERTIA=83,
147041acf15aSKris Buschelman                MATOP_LOAD=84,
147141acf15aSKris Buschelman                MATOP_IS_SYMMETRIC=85,
147241acf15aSKris Buschelman                MATOP_IS_HERMITIAN=86,
147341acf15aSKris Buschelman                MATOP_IS_STRUCTURALLY_SYMMETRIC=87,
147441acf15aSKris Buschelman                MATOP_PB_RELAX=88,
147541acf15aSKris Buschelman                MATOP_GET_VECS=89,
147641acf15aSKris Buschelman                MATOP_MAT_MULT=90,
147741acf15aSKris Buschelman                MATOP_MAT_MULT_SYMBOLIC=91,
147841acf15aSKris Buschelman                MATOP_MAT_MULT_NUMERIC=92,
147941acf15aSKris Buschelman                MATOP_PTAP=93,
148041acf15aSKris Buschelman                MATOP_PTAP_SYMBOLIC=94,
1481bc011b1eSHong Zhang                MATOP_PTAP_NUMERIC=95,
1482bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE=96,
1483bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97,
14847ba1a0bfSKris Buschelman                MATOP_MAT_MULTTRANSPOSE_NUMERIC=98,
14857ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_SEQAIJ=99,
14867ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_SEQAIJ=100,
14877ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_MPIAIJ=101,
148887d4246cSBarry Smith                MATOP_PTAP_NUMERIC_MPIAIJ=102,
1489ff1cced0SMatthew Knepley                MATOP_CONJUGATE=103,
1490ff1cced0SMatthew Knepley                MATOP_SET_SIZES=104,
1491f5edf698SHong Zhang                MATOP_SET_VALUES_ROW = 105,
1492ff1cced0SMatthew Knepley                MATOP_REAL_PART=106,
1493ff1cced0SMatthew Knepley                MATOP_IMAG_PART=107,
1494d5c63f83SSatish Balay                MATOP_GET_ROW_UTRIANGULAR=108,
14952bebee5dSHong Zhang                MATOP_RESTORE_ROW_UTRIANGULAR=109,
149669db28dcSHong Zhang                MATOP_MATSOLVE=110,
1497829201f2SHong Zhang                MATOP_GET_REDUNDANTMATRIX=111,
1498ff1cced0SMatthew Knepley                MATOP_GET_ROW_MIN=112,
1499ff1cced0SMatthew Knepley                MATOP_GET_COLUMN_VEC=113,
1500ff1cced0SMatthew Knepley                MATOP_MISSING_DIAGONAL=114,
1501ff1cced0SMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=115,
1502ff1cced0SMatthew Knepley                MATOP_DESTROY_SOLVER=116
1503fae171e0SBarry Smith              } MatOperation;
1504be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1505be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1506be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1508112a2221SBarry Smith 
150990ace30eSBarry Smith /*
151090ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
151190ace30eSBarry Smith  stored in a universal format. By changing the format with
15127973ad04SBarry Smith  PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
151390ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
151490ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
151590ace30eSBarry Smith  read into matrices of the same time.
151690ace30eSBarry Smith */
151790ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
151890ace30eSBarry Smith 
1519be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1520be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
15211f4e1ec7SBarry Smith 
1522d9274352SBarry Smith /*S
1523d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1524d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1525d9274352SBarry Smith 
1526f7a9e4ceSBarry Smith    Level: advanced
1527d9274352SBarry Smith 
1528d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1529d9274352SBarry Smith 
15306e1639daSBarry Smith   Users manual sections:
15316e1639daSBarry Smith .   sec_singular
15326e1639daSBarry Smith 
1533d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1534d9274352SBarry Smith S*/
153574637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1536d9274352SBarry Smith 
1537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1538281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*);
1539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
154295902228SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat,PetscTruth *);
154374637425SBarry Smith 
1544be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
154746129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
15483f1d51d7SBarry Smith 
1549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1552c4f061fbSSatish Balay 
1553be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1554b0a32e0cSBarry Smith 
1555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
155604f1ad80SBarry Smith 
1557fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
15583ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec);
15593ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1560e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1561e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1562e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace);
1563e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1564e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat);
1565e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal);
1566e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt);
1567e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *);
1568e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat);
1569e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1570e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1571e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal);
1572e884886fSBarry Smith 
1573e884886fSBarry Smith 
1574e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1575e884886fSBarry Smith 
1576e884886fSBarry Smith /*E
1577e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1578e884886fSBarry Smith 
1579e884886fSBarry Smith    Level: beginner
1580e884886fSBarry Smith 
1581e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1582e884886fSBarry Smith E*/
1583a313700dSBarry Smith #define MatMFFDType char*
1584e884886fSBarry Smith #define MATMFFD_DS  "ds"
1585e884886fSBarry Smith #define MATMFFD_WP  "wp"
1586e884886fSBarry Smith 
1587a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,const MatMFFDType);
1588e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1589e884886fSBarry Smith 
1590e884886fSBarry Smith /*MC
1591e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1592e884886fSBarry Smith 
1593e884886fSBarry Smith    Synopsis:
1594e884886fSBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1595e884886fSBarry Smith 
1596e884886fSBarry Smith    Not Collective
1597e884886fSBarry Smith 
1598e884886fSBarry Smith    Input Parameters:
1599e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1600e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1601e884886fSBarry Smith .  name_create - name of routine to create method context
1602e884886fSBarry Smith -  routine_create - routine to create method context
1603e884886fSBarry Smith 
1604e884886fSBarry Smith    Level: developer
1605e884886fSBarry Smith 
1606e884886fSBarry Smith    Notes:
1607e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1608e884886fSBarry Smith 
1609e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1610e884886fSBarry Smith    is ignored.
1611e884886fSBarry Smith 
1612e884886fSBarry Smith    Sample usage:
1613e884886fSBarry Smith .vb
1614e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1615e884886fSBarry Smith                "MyHCreate",MyHCreate);
1616e884886fSBarry Smith .ve
1617e884886fSBarry Smith 
1618e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1619e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1620e884886fSBarry Smith    or at runtime via the option
1621e884886fSBarry Smith $     -snes_mf_type my_h
1622e884886fSBarry Smith 
1623e884886fSBarry Smith .keywords: MatMFFD, register
1624e884886fSBarry Smith 
1625e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1626e884886fSBarry Smith M*/
1627e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1628e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1629e884886fSBarry Smith #else
1630e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1631e884886fSBarry Smith #endif
1632e884886fSBarry Smith 
1633e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]);
1634e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void);
1635e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal);
1636e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth);
1637e884886fSBarry Smith 
1638e884886fSBarry Smith 
1639be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1640be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
16417dbadf16SMatthew Knepley 
1642e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
16432eac72dbSBarry Smith #endif
1644