xref: /petsc/include/petscmat.h (revision 28b2fa4ae375c43e60412405d68670f082d116bf)
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"
66d91e6319SBarry Smith 
679989ab13SBarry Smith /*E
68c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
699989ab13SBarry Smith 
709989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
719989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
729989ab13SBarry Smith 
739989ab13SBarry Smith 
749989ab13SBarry Smith    Level: beginner
759989ab13SBarry Smith 
765c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
779989ab13SBarry Smith E*/
78c7393fdbSBarry Smith #define MatSolverPackage char*
79b24902e0SBarry Smith #define MAT_SOLVER_SPOOLES      "spooles"
80b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU      "superlu"
81b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU_DIST "superlu_dist"
82b24902e0SBarry Smith #define MAT_SOLVER_UMFPACK      "umfpack"
83b24902e0SBarry Smith #define MAT_SOLVER_ESSL         "essl"
84b24902e0SBarry Smith #define MAT_SOLVER_LUSOL        "lusol"
85b24902e0SBarry Smith #define MAT_SOLVER_MUMPS        "mumps"
86b24902e0SBarry Smith #define MAT_SOLVER_DSCPACK      "dscpack"
87b24902e0SBarry Smith #define MAT_SOLVER_MATLAB       "matlab"
8843244d56SBarry Smith #define MAT_SOLVER_PETSC        "petsc"
89b24902e0SBarry Smith 
90b24902e0SBarry Smith /*E
91b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
92b24902e0SBarry Smith 
93b24902e0SBarry Smith     Level: beginner
94b24902e0SBarry Smith 
95b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
96b24902e0SBarry Smith 
97c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
98b24902e0SBarry Smith E*/
995c9eb25fSBarry Smith typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC} MatFactorType;
100c7393fdbSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
101db4efbfdSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscTruth*);
10235bd34faSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
103b24902e0SBarry Smith 
1049989ab13SBarry Smith 
105c06d978dSMatthew Knepley /* Logging support */
106552e946dSBarry Smith #define    MAT_FILE_COOKIE 1211216    /* used to indicate matrices in binary files */
107be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE;
108be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
109be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE;
110be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE;
111e884886fSBarry Smith extern PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE;
112c06d978dSMatthew Knepley 
113ceb03754SKris Buschelman /*E
114ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
115ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
116ceb03754SKris Buschelman 
117ceb03754SKris Buschelman     Level: beginner
118ceb03754SKris Buschelman 
119ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
120ceb03754SKris Buschelman 
1210c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
122ceb03754SKris Buschelman E*/
123ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse;
124ceb03754SKris Buschelman 
1255494a064SHong Zhang /*E
1265494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
127829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1285494a064SHong Zhang 
1295494a064SHong Zhang     Level: beginner
1305494a064SHong Zhang 
131829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1325494a064SHong Zhang E*/
1335494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1345494a064SHong Zhang 
135e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]);
136c06d978dSMatthew Knepley 
137f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
138a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
139a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
140f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
141a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,const MatType);
142be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
143be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
144be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
145be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
14630de9b25SBarry Smith 
147f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
148f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
149f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
150f69a0ea3SMatthew Knepley 
15130de9b25SBarry Smith /*MC
15230de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
15330de9b25SBarry Smith 
15430de9b25SBarry Smith    Synopsis:
155c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
15630de9b25SBarry Smith 
15730de9b25SBarry Smith    Not Collective
15830de9b25SBarry Smith 
15930de9b25SBarry Smith    Input Parameters:
16030de9b25SBarry Smith +  name - name of a new user-defined matrix type
16130de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
16230de9b25SBarry Smith .  name_create - name of routine to create method context
16330de9b25SBarry Smith -  routine_create - routine to create method context
16430de9b25SBarry Smith 
16530de9b25SBarry Smith    Notes:
16630de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
16730de9b25SBarry Smith 
16830de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
16930de9b25SBarry Smith    is ignored.
17030de9b25SBarry Smith 
17130de9b25SBarry Smith    Sample usage:
17230de9b25SBarry Smith .vb
17330de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
17430de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
17530de9b25SBarry Smith .ve
17630de9b25SBarry Smith 
17730de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
17830de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
17930de9b25SBarry Smith    or at runtime via the option
18030de9b25SBarry Smith $     -mat_type my_mat
18130de9b25SBarry Smith 
18230de9b25SBarry Smith    Level: advanced
18330de9b25SBarry Smith 
184ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
18530de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
18630de9b25SBarry Smith 
18730de9b25SBarry Smith .keywords: Mat, register
18830de9b25SBarry Smith 
18930de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
19030de9b25SBarry Smith 
19130de9b25SBarry Smith M*/
192273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
193273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
194273d9f13SBarry Smith #else
195273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
19630de9b25SBarry Smith #endif
19730de9b25SBarry Smith 
198273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
199b0a32e0cSBarry Smith extern PetscFList MatList;
20028988994SBarry Smith 
2013b224e63SBarry Smith /*E
2023b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
2033b224e63SBarry Smith 
2043b224e63SBarry Smith     Level: beginner
2053b224e63SBarry Smith 
2063b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
2073b224e63SBarry Smith 
2083b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
2093b224e63SBarry Smith E*/
2103b224e63SBarry Smith typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
2113b224e63SBarry Smith 
212be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
213be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
214be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
215ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
216ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
217ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
218ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
219ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
220ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
221ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
222be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
223ba966639SSatish 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)
224ba966639SSatish 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)
225ba966639SSatish 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)
226ba966639SSatish 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)
227ba966639SSatish 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))
228ba966639SSatish 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))
229ba966639SSatish 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))
230ba966639SSatish 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)
231ba966639SSatish 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)
232ba966639SSatish 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)
233ba966639SSatish 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)
234ba966639SSatish 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))
235ba966639SSatish 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))
236ba966639SSatish 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))
23763ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2388d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2398d7a6e47SBarry Smith 
240be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
241be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
242ba966639SSatish 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)
243ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
244ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
245ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
246ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
247ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
248ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
249be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
250ba966639SSatish 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)
251ba966639SSatish 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)
252ba966639SSatish 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)
253ba966639SSatish 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)
254ba966639SSatish 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))
255ba966639SSatish 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))
256ba966639SSatish 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))
257ba966639SSatish 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)
258ba966639SSatish 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)
259ba966639SSatish 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)
260ba966639SSatish 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)
261ba966639SSatish 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))
262ba966639SSatish 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))
263ba966639SSatish 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))
264be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
265be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
266ba966639SSatish 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)
267ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
268ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
269ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
270ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
271ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
272ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
273be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
274ba966639SSatish 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)
275ba966639SSatish 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)
276ba966639SSatish 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)
277ba966639SSatish 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)
278ba966639SSatish 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))
279ba966639SSatish 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))
280ba966639SSatish 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))
281ba966639SSatish 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)
282ba966639SSatish 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)
283ba966639SSatish 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)
284ba966639SSatish 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)
285ba966639SSatish 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))
286ba966639SSatish 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))
287ba966639SSatish 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))
288be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
289ba966639SSatish 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)
290ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
291be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
292be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
293ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
294ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*);
295284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
2961472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
2971472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
2982a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*);
2992a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter);
3002a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*);
3018cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
302793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat);
303b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat);
304793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
3056d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
3066d7c1e57SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeSetType(Mat,MatCompositeType);
3076d7c1e57SBarry Smith 
3085a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*);
309f20085c4SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateTranspose(Mat,Mat*);
3101472f72bSBarry Smith 
311f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
312be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
31321c89e3eSBarry Smith 
3140c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
31599cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat);
31699cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat);
31799cafbc1SBarry Smith 
3188ed539a5SBarry Smith /* ------------------------------------------------------------*/
319be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
320be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
32187d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
322f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
32384cb2905SBarry Smith 
3242ef4de8bSBarry Smith /*S
3252ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
3262ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3272ef4de8bSBarry Smith 
3282ef4de8bSBarry Smith    Level: beginner
3292ef4de8bSBarry Smith 
3302ef4de8bSBarry Smith   Concepts: matrix; linear operator
3312ef4de8bSBarry Smith 
332d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3332ef4de8bSBarry Smith S*/
334435da068SBarry Smith typedef struct {
335c1ac3661SBarry Smith   PetscInt k,j,i,c;
336435da068SBarry Smith } MatStencil;
3372ef4de8bSBarry Smith 
338be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
339be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
340be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
341435da068SBarry Smith 
342be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
343be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
344be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
3453a7fca6bSBarry Smith 
346d91e6319SBarry Smith /*E
347d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
348d91e6319SBarry Smith      to continue to add values to it
349d91e6319SBarry Smith 
350d91e6319SBarry Smith     Level: beginner
351d91e6319SBarry Smith 
352d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
353d91e6319SBarry Smith E*/
3546d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
355be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
356be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
357be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
3584f9c727eSBarry Smith 
3591d73ed98SBarry Smith 
36030de9b25SBarry Smith 
361d91e6319SBarry Smith /*E
362d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
363d91e6319SBarry Smith 
364d91e6319SBarry Smith     Level: beginner
365d91e6319SBarry Smith 
3660a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
367d91e6319SBarry Smith 
368d91e6319SBarry Smith .seealso: MatSetOption()
369d91e6319SBarry Smith E*/
370512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
3714e0d8c25SBarry Smith               MAT_SYMMETRIC,
3724e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
3738047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
3744e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
3754e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
3764e0d8c25SBarry Smith               MAT_KEEP_ZEROED_ROWS,MAT_IGNORE_ZERO_ENTRIES,
3774e0d8c25SBarry Smith               MAT_USE_INODES,
3784e0d8c25SBarry Smith               MAT_HERMITIAN,
3794e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
3804e0d8c25SBarry Smith               MAT_USE_COMPRESSEDROW,
3814e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
382*28b2fa4aSMatthew Knepley               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
383*28b2fa4aSMatthew Knepley               NUM_MAT_OPTIONS} MatOption;
384290bbb0aSBarry Smith extern const char *MatOptions[];
3854e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth);
386a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,const MatType*);
387a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
38884cb2905SBarry Smith 
389be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
390be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
391be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
392f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat);
393f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat);
394be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
395be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
396be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
398ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
400be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
401ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
402be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
4037b80b807SBarry Smith 
4041620fd73SBarry Smith 
405be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
406be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
407ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
408be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
409be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
410ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
411e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t)
4121cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*);
413be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
414ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
415be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
416be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
4172bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat);
418f5edf698SHong Zhang 
419d91e6319SBarry Smith /*E
420d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
421d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
422d91e6319SBarry Smith 
423d91e6319SBarry Smith     Level: beginner
424d91e6319SBarry Smith 
425d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
426d91e6319SBarry Smith 
427d91e6319SBarry Smith .seealso: MatDuplicate()
428d91e6319SBarry Smith E*/
4292e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
4302e8a6d31SBarry Smith 
431a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,const MatType,MatReuse,Mat*);
432a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
433be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
434ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
435ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
43694a9d846SBarry Smith 
437cb5b572fSBarry Smith 
438be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
439be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
440be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
441ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
442e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t)
443be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
444ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
4451cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*);
4461cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t)
447be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
44964c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *);
450a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,const MatType,Mat*);
451a313700dSBarry Smith PetscPolymorphicFunction(MatLoad,(PetscViewer v,const MatType t),(v,t,&a),Mat,a)
4527b80b807SBarry Smith 
4538f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4548f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
4558f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4568f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
457d4fbbf0eSBarry Smith 
458d91e6319SBarry Smith /*S
459d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
460d91e6319SBarry Smith 
461d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
462d91e6319SBarry Smith 
463d91e6319SBarry Smith    Level: intermediate
464d91e6319SBarry Smith 
465d91e6319SBarry Smith   Concepts: matrix^nonzero information
466d91e6319SBarry Smith 
467d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
468d91e6319SBarry Smith S*/
4694e220ebcSLois Curfman McInnes typedef struct {
470b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
471b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
472b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
473b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
474b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
475b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
476b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4774e220ebcSLois Curfman McInnes } MatInfo;
4784e220ebcSLois Curfman McInnes 
479d9274352SBarry Smith /*E
480d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
481d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
482d9274352SBarry Smith 
483d9274352SBarry Smith     Level: beginner
484d9274352SBarry Smith 
485d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
486d9274352SBarry Smith 
487d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
488d9274352SBarry Smith E*/
4897b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
490be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
491be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
492be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
493985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]);
494985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]);
495985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
496c87e5d42SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMinAbs(Mat,Vec,PetscInt[]);
49786697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec);
498fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*);
499fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
501ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
502be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
503be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
504be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
505be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
506ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
508be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
509be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
5117b80b807SBarry Smith 
512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
513ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
515f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
516f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
517f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
518f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5197b80b807SBarry Smith 
520be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
521be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
522be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5235ef9f2a5SBarry Smith 
524be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
525be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
526be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5278ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**);
528f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
52972ca8751SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5307b80b807SBarry Smith 
531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *);
534abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *);
535829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat *[]);
536829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat *[]);
5375494a064SHong Zhang 
538be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
542be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
543be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
544be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
546dd6ea824SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,MatScalar**,Mat*);
54743eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
548cd116e26SSatish Balay #include "petscctable.h"
54943eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
55043eb5e2fSMatthew Knepley #else
55143eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
55243eb5e2fSMatthew Knepley #endif
5538efafbd8SBarry Smith 
554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5557b80b807SBarry Smith 
556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
557be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
558be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
55922440eb1SKris Buschelman 
560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
561be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
562be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
56322440eb1SKris Buschelman 
564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
565be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
566be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
567bc011b1eSHong Zhang 
568f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
56904aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure);
570be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat);
5711c741599SHong Zhang 
572f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
573f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5747b80b807SBarry Smith 
575be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
576be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
577f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
578f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
579be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
580be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
581052efed2SBarry Smith 
582be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
583be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
58490f02eecSBarry Smith 
585be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
5860c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
587ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
588be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
589be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
59069db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
591bd481603SBarry Smith 
592bd481603SBarry Smith /*MC
5931620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5941620fd73SBarry Smith 
5951620fd73SBarry Smith    Not collective
5961620fd73SBarry Smith 
5971620fd73SBarry Smith    Input Parameters:
5981620fd73SBarry Smith +  m - the matrix
5991620fd73SBarry Smith .  row - the row location of the entry
6001620fd73SBarry Smith .  col - the column location of the entry
6011620fd73SBarry Smith .  value - the value to insert
6021620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
6031620fd73SBarry Smith 
6041620fd73SBarry Smith    Notes:
6051620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6061620fd73SBarry Smith    values simultaneously if possible.
6071620fd73SBarry Smith 
6081620fd73SBarry Smith    Level: beginner
6091620fd73SBarry Smith 
6101620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6111620fd73SBarry Smith M*/
6121620fd73SBarry 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);}
6131620fd73SBarry Smith 
6141620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6151620fd73SBarry Smith 
6161620fd73SBarry 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);}
6171620fd73SBarry Smith 
6181620fd73SBarry Smith /*MC
6190d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
620bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
621bd481603SBarry Smith 
622bd481603SBarry Smith    Synopsis:
623c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
624bd481603SBarry Smith 
625bd481603SBarry Smith    Collective on MPI_Comm
626bd481603SBarry Smith 
627bd481603SBarry Smith    Input Parameters:
628bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
629859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
630859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
631bd481603SBarry Smith 
632bd481603SBarry Smith    Output Parameters:
633bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
634bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
635bd481603SBarry Smith 
636bd481603SBarry Smith 
637bd481603SBarry Smith    Level: intermediate
638bd481603SBarry Smith 
639bd481603SBarry Smith    Notes:
640bd481603SBarry Smith    See the chapter in the users manual on performance for more details
641bd481603SBarry Smith 
6421d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
643bd481603SBarry Smith 
644bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
645bd481603SBarry Smith 
6461620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6471620fd73SBarry Smith 
648bd481603SBarry Smith   Concepts: preallocation^Matrix
649bd481603SBarry Smith 
650bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
651bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
652bd481603SBarry Smith M*/
653c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6547c922b88SBarry Smith { \
655a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
656a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
657a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
658a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
659c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
660a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6617c922b88SBarry Smith 
662bd481603SBarry Smith /*MC
6630d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
664bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
665bd481603SBarry Smith 
666bd481603SBarry Smith    Synopsis:
667c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
668bd481603SBarry Smith 
669bd481603SBarry Smith    Collective on MPI_Comm
670bd481603SBarry Smith 
671bd481603SBarry Smith    Input Parameters:
672bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
673859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
674859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
675bd481603SBarry Smith 
676bd481603SBarry Smith    Output Parameters:
677bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
678bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
679bd481603SBarry Smith 
680bd481603SBarry Smith 
681bd481603SBarry Smith    Level: intermediate
682bd481603SBarry Smith 
683bd481603SBarry Smith    Notes:
684bd481603SBarry Smith    See the chapter in the users manual on performance for more details
685bd481603SBarry Smith 
6861d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
687bd481603SBarry Smith 
6881620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6891620fd73SBarry Smith 
690bd481603SBarry Smith   Concepts: preallocation^Matrix
691bd481603SBarry Smith 
692bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
693bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
694bd481603SBarry Smith M*/
695222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
696222b16d4SBarry Smith { \
697a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
698a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
699a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
700a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
701c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
702a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
703222b16d4SBarry Smith 
704bd481603SBarry Smith /*MC
705bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
706bd481603SBarry Smith        inserted using a local number of the rows and columns
707bd481603SBarry Smith 
708bd481603SBarry Smith    Synopsis:
709c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
710bd481603SBarry Smith 
711bd481603SBarry Smith    Not Collective
712bd481603SBarry Smith 
713bd481603SBarry Smith    Input Parameters:
714bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
715bd481603SBarry Smith .  nrows - the number of rows indicated
7161d73ed98SBarry Smith .  rows - the indices of the rows
717bd481603SBarry Smith .  ncols - the number of columns in the matrix
718bd481603SBarry Smith .  cols - the columns indicated
719bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
720bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
721bd481603SBarry Smith 
722bd481603SBarry Smith 
723bd481603SBarry Smith    Level: intermediate
724bd481603SBarry Smith 
725bd481603SBarry Smith    Notes:
726bd481603SBarry Smith    See the chapter in the users manual on performance for more details
727bd481603SBarry Smith 
7281d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
729bd481603SBarry Smith 
730bd481603SBarry Smith   Concepts: preallocation^Matrix
731bd481603SBarry Smith 
732bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
733bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
734bd481603SBarry Smith M*/
735c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
736c4f061fbSSatish Balay {\
737c1ac3661SBarry Smith   PetscInt __l;\
738ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
739ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
740c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
741ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
742c4f061fbSSatish Balay   }\
743c4f061fbSSatish Balay }
744c4f061fbSSatish Balay 
745bd481603SBarry Smith /*MC
746bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
747bd481603SBarry Smith        inserted using a local number of the rows and columns
748bd481603SBarry Smith 
749bd481603SBarry Smith    Synopsis:
750c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
751bd481603SBarry Smith 
752bd481603SBarry Smith    Not Collective
753bd481603SBarry Smith 
754bd481603SBarry Smith    Input Parameters:
755bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
756bd481603SBarry Smith .  nrows - the number of rows indicated
7571d73ed98SBarry Smith .  rows - the indices of the rows
758bd481603SBarry Smith .  ncols - the number of columns in the matrix
759bd481603SBarry Smith .  cols - the columns indicated
760bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
761bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
762bd481603SBarry Smith 
763bd481603SBarry Smith 
764bd481603SBarry Smith    Level: intermediate
765bd481603SBarry Smith 
766bd481603SBarry Smith    Notes:
767bd481603SBarry Smith    See the chapter in the users manual on performance for more details
768bd481603SBarry Smith 
769bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
770bd481603SBarry Smith 
771bd481603SBarry Smith   Concepts: preallocation^Matrix
772bd481603SBarry Smith 
773bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
774bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
775bd481603SBarry Smith M*/
776d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
777d3d32019SBarry Smith {\
778c1ac3661SBarry Smith   PetscInt __l;\
779d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
780d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
781d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
782d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
783d3d32019SBarry Smith   }\
784d3d32019SBarry Smith }
785d3d32019SBarry Smith 
786bd481603SBarry Smith /*MC
787bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
788bd481603SBarry Smith        inserted using a local number of the rows and columns
789bd481603SBarry Smith 
790bd481603SBarry Smith    Synopsis:
791c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
792bd481603SBarry Smith 
793bd481603SBarry Smith    Not Collective
794bd481603SBarry Smith 
795bd481603SBarry Smith    Input Parameters:
79664075487SBarry Smith +  row - the row
797bd481603SBarry Smith .  ncols - the number of columns in the matrix
798a50f8bf6SHong Zhang -  cols - the columns indicated
799a50f8bf6SHong Zhang 
800a50f8bf6SHong Zhang    Output Parameters:
801a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
802bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
803bd481603SBarry Smith 
804bd481603SBarry Smith 
805bd481603SBarry Smith    Level: intermediate
806bd481603SBarry Smith 
807bd481603SBarry Smith    Notes:
808bd481603SBarry Smith    See the chapter in the users manual on performance for more details
809bd481603SBarry Smith 
810bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
811bd481603SBarry Smith 
8121620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8131620fd73SBarry Smith 
814bd481603SBarry Smith   Concepts: preallocation^Matrix
815bd481603SBarry Smith 
816bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
817bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
818bd481603SBarry Smith M*/
819c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
820c1ac3661SBarry Smith { PetscInt __i; \
8218ee2e534SBarry Smith   if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
822a89bc211SBarry 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);\
8237c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
82464075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8257cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8267c922b88SBarry Smith   }\
8277c922b88SBarry Smith }
8287c922b88SBarry Smith 
829bd481603SBarry Smith /*MC
830bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
831bd481603SBarry Smith        inserted using a local number of the rows and columns
832bd481603SBarry Smith 
833bd481603SBarry Smith    Synopsis:
834c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
835bd481603SBarry Smith 
836bd481603SBarry Smith    Not Collective
837bd481603SBarry Smith 
838bd481603SBarry Smith    Input Parameters:
839bd481603SBarry Smith +  nrows - the number of rows indicated
8401d73ed98SBarry Smith .  rows - the indices of the rows
841bd481603SBarry Smith .  ncols - the number of columns in the matrix
842bd481603SBarry Smith .  cols - the columns indicated
843bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
844bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
845bd481603SBarry Smith 
846bd481603SBarry Smith 
847bd481603SBarry Smith    Level: intermediate
848bd481603SBarry Smith 
849bd481603SBarry Smith    Notes:
850bd481603SBarry Smith    See the chapter in the users manual on performance for more details
851bd481603SBarry Smith 
852bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
853bd481603SBarry Smith 
8541620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8551620fd73SBarry Smith 
856bd481603SBarry Smith   Concepts: preallocation^Matrix
857bd481603SBarry Smith 
858bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
859bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
860bd481603SBarry Smith M*/
861d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
862c1ac3661SBarry Smith { PetscInt __i; \
863d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
864d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
865d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
866d3d32019SBarry Smith   }\
867d3d32019SBarry Smith }
868d3d32019SBarry Smith 
869bd481603SBarry Smith /*MC
87016371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
87116371a99SBarry Smith 
87216371a99SBarry Smith    Synopsis:
87316371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
87416371a99SBarry Smith 
87516371a99SBarry Smith    Not Collective
87616371a99SBarry Smith 
87716371a99SBarry Smith    Input Parameters:
87816371a99SBarry Smith .  A - matrix
87916371a99SBarry Smith .  row - row where values exist (must be local to this process)
88016371a99SBarry Smith .  ncols - number of columns
88116371a99SBarry Smith .  cols - columns with nonzeros
88216371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
88316371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
88416371a99SBarry Smith 
88516371a99SBarry Smith 
88616371a99SBarry Smith    Level: intermediate
88716371a99SBarry Smith 
88816371a99SBarry Smith    Notes:
88916371a99SBarry Smith    See the chapter in the users manual on performance for more details
89016371a99SBarry Smith 
89116371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
89216371a99SBarry Smith 
89316371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
89416371a99SBarry Smith 
89516371a99SBarry Smith   Concepts: preallocation^Matrix
89616371a99SBarry Smith 
89716371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
89816371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
89916371a99SBarry Smith M*/
90016371a99SBarry 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);}
90116371a99SBarry Smith 
90216371a99SBarry Smith 
90316371a99SBarry Smith /*MC
9040d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
905bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
906bd481603SBarry Smith 
907bd481603SBarry Smith    Synopsis:
908c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
909bd481603SBarry Smith 
910bd481603SBarry Smith    Collective on MPI_Comm
911bd481603SBarry Smith 
912bd481603SBarry Smith    Input Parameters:
91316371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
914bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
915bd481603SBarry Smith 
916bd481603SBarry Smith 
917bd481603SBarry Smith    Level: intermediate
918bd481603SBarry Smith 
919bd481603SBarry Smith    Notes:
920bd481603SBarry Smith    See the chapter in the users manual on performance for more details
921bd481603SBarry Smith 
922bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
923bd481603SBarry Smith 
9241620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9251620fd73SBarry Smith 
926bd481603SBarry Smith   Concepts: preallocation^Matrix
927bd481603SBarry Smith 
928bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
929bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
930bd481603SBarry Smith M*/
931a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9327c922b88SBarry Smith 
933bd481603SBarry Smith 
934bd481603SBarry Smith 
9357b80b807SBarry Smith /* Routines unique to particular data structures */
936be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
937ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
938698d4c6aSKris Buschelman 
939be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
940be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
941698d4c6aSKris Buschelman 
942be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
943be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
944be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
945c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
946c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9477b80b807SBarry Smith 
948a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
949a96a251dSBarry Smith 
950be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
951ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
952be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
953ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
954be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
955ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
956273d9f13SBarry Smith 
957be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
958ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
959be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
960be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
96153803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
962725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
963be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
964be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
965be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
966be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
967284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
968be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
969be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
970be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
971be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
972273d9f13SBarry Smith 
973be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
974ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*);
9751b807ce4Svictorle 
976be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
977be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
9782e8a6d31SBarry Smith 
979be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
9803a7fca6bSBarry Smith 
9817b80b807SBarry Smith /*
9827b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
98394b7f48cSBarry Smith   done through the KSP and PC interfaces.
9847b80b807SBarry Smith */
9857b80b807SBarry Smith 
986d9274352SBarry Smith /*E
987d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
988d9274352SBarry Smith        with an optional dynamic library name, for example
989d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
990d9274352SBarry Smith 
991d9274352SBarry Smith    Level: beginner
992d9274352SBarry Smith 
993e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
994e5a9bf91SBarry Smith 
995d9274352SBarry Smith .seealso: MatGetOrdering()
996d9274352SBarry Smith E*/
9973eaad2caSSatish Balay #define MatOrderingType char*
998b12f92e5SBarry Smith #define MATORDERING_NATURAL   "natural"
999b12f92e5SBarry Smith #define MATORDERING_ND        "nd"
1000b12f92e5SBarry Smith #define MATORDERING_1WD       "1wd"
1001b12f92e5SBarry Smith #define MATORDERING_RCM       "rcm"
1002b12f92e5SBarry Smith #define MATORDERING_QMD       "qmd"
1003b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength"
100462152c8bSBarry Smith #define MATORDERING_DSC_ND    "dsc_nd"
100562152c8bSBarry Smith #define MATORDERING_DSC_MMD   "dsc_mmd"
100662152c8bSBarry Smith #define MATORDERING_DSC_MDF   "dsc_mdf"
1007c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
1008c06d978dSMatthew Knepley #define MATORDERING_IDENTITY  "identity"
1009c06d978dSMatthew Knepley #define MATORDERING_REVERSE   "reverse"
1010b12f92e5SBarry Smith 
1011f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
1012be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
1013f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
101430de9b25SBarry Smith 
101530de9b25SBarry Smith /*MC
101630de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
101730de9b25SBarry Smith                                matrix package.
101830de9b25SBarry Smith 
101930de9b25SBarry Smith    Synopsis:
1020c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
102130de9b25SBarry Smith 
102230de9b25SBarry Smith    Not Collective
102330de9b25SBarry Smith 
102430de9b25SBarry Smith    Input Parameters:
102530de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
102630de9b25SBarry Smith .  path - location of library where creation routine is
102730de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
102830de9b25SBarry Smith -  function - function pointer that creates the ordering
102930de9b25SBarry Smith 
103030de9b25SBarry Smith    Level: developer
103130de9b25SBarry Smith 
103230de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
103330de9b25SBarry Smith    is ignored.
103430de9b25SBarry Smith 
103530de9b25SBarry Smith    Sample usage:
103630de9b25SBarry Smith .vb
103730de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
103830de9b25SBarry Smith                "MyOrder",MyOrder);
103930de9b25SBarry Smith .ve
104030de9b25SBarry Smith 
104130de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
104230de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
104330de9b25SBarry Smith    or at runtime via the option
10442401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
104530de9b25SBarry Smith 
1046ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
104730de9b25SBarry Smith 
104830de9b25SBarry Smith .keywords: matrix, ordering, register
104930de9b25SBarry Smith 
105030de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
105130de9b25SBarry Smith M*/
1052aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1053f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1054b12f92e5SBarry Smith #else
1055f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1056b12f92e5SBarry Smith #endif
105730de9b25SBarry Smith 
1058be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
1059be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
10602bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
1061b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1062d4fbbf0eSBarry Smith 
1063be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1064a2ce50c7SBarry Smith 
1065d91e6319SBarry Smith /*S
10662401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1067b00f7748SHong Zhang 
106861cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
106961cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1070b00f7748SHong Zhang 
107115e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1072b00f7748SHong Zhang 
107361cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
107461cad860SBarry Smith 
1075b00f7748SHong Zhang    Level: developer
1076b00f7748SHong Zhang 
1077d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1078d7d82daaSBarry Smith           MatFactorInfoInitialize()
1079b00f7748SHong Zhang 
1080b00f7748SHong Zhang S*/
1081b00f7748SHong Zhang typedef struct {
10820a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
1083fbf22428SSatish Balay   PetscReal     shiftpd;        /* if true, shift until positive pivots */
108415e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
108585317021SBarry Smith   PetscReal     usedt;
108615e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1087b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
108815e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
108967e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1090348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1091bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1092bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
109362bba022SBarry Smith   PetscReal     shiftinblocks;  /* if block in block factorization has zero pivot then shift diagonal until non-singular */
109415e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
109515e8a5b3SHong Zhang } MatFactorInfo;
1096ffa6d0a5SLois Curfman McInnes 
1097be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
10980481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
10990481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11000481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
11010481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
11020481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
11030481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11040481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11050481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11060481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,const MatFactorInfo*);
11070481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
11080481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,const MatFactorInfo*,Mat *);
1109be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1110be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1111be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1112be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1113be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1114be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1115be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1116be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
11178ed539a5SBarry Smith 
1118be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1119bb5a7306SBarry Smith 
1120d91e6319SBarry Smith /*E
1121d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1122bb1eb677SSatish Balay 
1123d91e6319SBarry Smith     Level: beginner
1124d91e6319SBarry Smith 
1125d9274352SBarry Smith    May be bitwise ORd together
1126d9274352SBarry Smith 
1127d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1128d91e6319SBarry Smith 
11294e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11304e7234bfSBarry Smith 
1131a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax()
1132d91e6319SBarry Smith E*/
1133ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1134ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1135ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
113684cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1137be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1138be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11398ed539a5SBarry Smith 
1140d4fbbf0eSBarry Smith /*
1141639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1142639f9d9dSBarry Smith */
1143b12f92e5SBarry Smith 
1144d9274352SBarry Smith /*E
1145d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1146d9274352SBarry Smith        with an optional dynamic library name, for example
1147d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1148d9274352SBarry Smith 
1149d9274352SBarry Smith    Level: beginner
1150d9274352SBarry Smith 
1151d9274352SBarry Smith .seealso: MatGetColoring()
1152d9274352SBarry Smith E*/
1153a313700dSBarry Smith #define MatColoringType char*
1154b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1155b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1156b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1157b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1158b12f92e5SBarry Smith 
1159ab970d9bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*);
11602e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
116130de9b25SBarry Smith 
116230de9b25SBarry Smith /*MC
116330de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
116430de9b25SBarry Smith                                matrix package.
116530de9b25SBarry Smith 
116630de9b25SBarry Smith    Synopsis:
1167c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
116830de9b25SBarry Smith 
116930de9b25SBarry Smith    Not Collective
117030de9b25SBarry Smith 
117130de9b25SBarry Smith    Input Parameters:
117230de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
117330de9b25SBarry Smith .  path - location of library where creation routine is
117430de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
117530de9b25SBarry Smith -  function - function pointer that creates the coloring
117630de9b25SBarry Smith 
117730de9b25SBarry Smith    Level: developer
117830de9b25SBarry Smith 
117930de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
118030de9b25SBarry Smith    is ignored.
118130de9b25SBarry Smith 
118230de9b25SBarry Smith    Sample usage:
118330de9b25SBarry Smith .vb
118430de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
118530de9b25SBarry Smith                "MyColor",MyColor);
118630de9b25SBarry Smith .ve
118730de9b25SBarry Smith 
118830de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
118930de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
119030de9b25SBarry Smith    or at runtime via the option
119130de9b25SBarry Smith $     -mat_coloring_type my_color
119230de9b25SBarry Smith 
1193ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
119430de9b25SBarry Smith 
119530de9b25SBarry Smith .keywords: matrix, Coloring, register
119630de9b25SBarry Smith 
119730de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
119830de9b25SBarry Smith M*/
1199aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1200f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1201b12f92e5SBarry Smith #else
1202f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1203b12f92e5SBarry Smith #endif
120430de9b25SBarry Smith 
12052bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1206f1144a30SSatish Balay 
1207be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1208be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1209be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1210639f9d9dSBarry Smith 
1211d9274352SBarry Smith /*S
1212d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1213d9274352SBarry Smith         and coloring
1214639f9d9dSBarry Smith 
1215d9274352SBarry Smith    Level: beginner
1216d9274352SBarry Smith 
1217d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1218d9274352SBarry Smith 
1219d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1220d9274352SBarry Smith S*/
1221e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1222639f9d9dSBarry Smith 
1223be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1224be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1225be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1226be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
12274a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1228be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1229be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1230be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1231be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1232be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1233be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1234639f9d9dSBarry Smith /*
12350752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12363eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12370752156aSBarry Smith */
1238ca161407SBarry Smith 
1239d9274352SBarry Smith /*S
1240d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1241d9274352SBarry Smith 
1242d9274352SBarry Smith    Level: beginner
1243d9274352SBarry Smith 
1244d9274352SBarry Smith   Concepts: partitioning
1245d9274352SBarry Smith 
1246743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1247d9274352SBarry Smith S*/
124891e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1249d9274352SBarry Smith 
1250d9274352SBarry Smith /*E
12515bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1252d9274352SBarry Smith        with an optional dynamic library name, for example
1253d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1254d9274352SBarry Smith 
1255d9274352SBarry Smith    Level: beginner
1256d9274352SBarry Smith 
1257b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1258d9274352SBarry Smith E*/
1259a313700dSBarry Smith #define MatPartitioningType char*
12608ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1261c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
12628ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
126371306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
126471306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
126571306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
126671306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
126771306c60Sjroman 
1268ca161407SBarry Smith 
1269be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1270a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1271be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1272be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1273be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1274be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1275be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
12772aabb6bbSBarry Smith 
1278be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
127930de9b25SBarry Smith 
128030de9b25SBarry Smith /*MC
128130de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
128230de9b25SBarry Smith    matrix package.
128330de9b25SBarry Smith 
128430de9b25SBarry Smith    Synopsis:
1285c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
128630de9b25SBarry Smith 
128730de9b25SBarry Smith    Not Collective
128830de9b25SBarry Smith 
128930de9b25SBarry Smith    Input Parameters:
129030de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
129130de9b25SBarry Smith .  path - location of library where creation routine is
129230de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
129330de9b25SBarry Smith -  function - function pointer that creates the partitioning type
129430de9b25SBarry Smith 
129530de9b25SBarry Smith    Level: developer
129630de9b25SBarry Smith 
129730de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
129830de9b25SBarry Smith    is ignored.
129930de9b25SBarry Smith 
130030de9b25SBarry Smith    Sample usage:
130130de9b25SBarry Smith .vb
130230de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
130330de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
130430de9b25SBarry Smith .ve
130530de9b25SBarry Smith 
130630de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
130730de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
130830de9b25SBarry Smith    or at runtime via the option
130930de9b25SBarry Smith $     -mat_partitioning_type my_part
131030de9b25SBarry Smith 
1311ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
131230de9b25SBarry Smith 
131330de9b25SBarry Smith .keywords: matrix, partitioning, register
131430de9b25SBarry Smith 
131530de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
131630de9b25SBarry Smith M*/
1317aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1318f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13192aabb6bbSBarry Smith #else
1320f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13212aabb6bbSBarry Smith #endif
13222aabb6bbSBarry Smith 
13232bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1324f1144a30SSatish Balay 
1325be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1326be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
13272bad1931SBarry Smith 
1328be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1329be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1330a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1331ca161407SBarry Smith 
1332be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
13330e8af3caSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
13340752156aSBarry Smith 
1335be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1336be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
133771306c60Sjroman 
1338954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
1339be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
134071306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1341be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1342be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
134371306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1344be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1345be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1346be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
134771306c60Sjroman 
134871306c60Sjroman #define MP_PARTY_OPT "opt"
134971306c60Sjroman #define MP_PARTY_LIN "lin"
135071306c60Sjroman #define MP_PARTY_SCA "sca"
135171306c60Sjroman #define MP_PARTY_RAN "ran"
135271306c60Sjroman #define MP_PARTY_GBF "gbf"
135371306c60Sjroman #define MP_PARTY_GCF "gcf"
135471306c60Sjroman #define MP_PARTY_BUB "bub"
135571306c60Sjroman #define MP_PARTY_DEF "def"
1356be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
135771306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
135871306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
135971306c60Sjroman #define MP_PARTY_NONE "no"
1360be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1361be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1362be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1363be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
136471306c60Sjroman 
136571306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1366be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1367be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1368be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1370be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
137171306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1373be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1374be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
137571306c60Sjroman 
13760752156aSBarry Smith /*
13770a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1378d4fbbf0eSBarry Smith */
13791c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13801c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13811c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13821c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13831c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13847c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13857c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13861c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13871c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13887c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13897c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13901c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13911c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
13921c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
13931c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13941c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13951c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13961c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13971c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13981c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13991c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
14001c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
14011c1c02c0SLois Curfman McInnes                MATOP_COMPRESS=22,
14021c1c02c0SLois Curfman McInnes                MATOP_SET_OPTION=23,
14031c1c02c0SLois Curfman McInnes                MATOP_ZERO_ENTRIES=24,
14041c1c02c0SLois Curfman McInnes                MATOP_ZERO_ROWS=25,
14051c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_SYMBOLIC=26,
14061c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_NUMERIC=27,
14071c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_SYMBOLIC=28,
14081c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_NUMERIC=29,
1409d643ce63SMatthew Knepley                MATOP_SETUP_PREALLOCATION=30,
1410d643ce63SMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=31,
1411d643ce63SMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=32,
1412d643ce63SMatthew Knepley                MATOP_GET_ARRAY=33,
1413d643ce63SMatthew Knepley                MATOP_RESTORE_ARRAY=34,
1414e26f1c3cSBarry Smith                MATOP_DUPLICATE=35,
1415d643ce63SMatthew Knepley                MATOP_FORWARD_SOLVE=36,
1416d643ce63SMatthew Knepley                MATOP_BACKWARD_SOLVE=37,
1417d643ce63SMatthew Knepley                MATOP_ILUFACTOR=38,
1418d643ce63SMatthew Knepley                MATOP_ICCFACTOR=39,
1419d643ce63SMatthew Knepley                MATOP_AXPY=40,
1420d643ce63SMatthew Knepley                MATOP_GET_SUBMATRICES=41,
1421d643ce63SMatthew Knepley                MATOP_INCREASE_OVERLAP=42,
1422d643ce63SMatthew Knepley                MATOP_GET_VALUES=43,
1423d643ce63SMatthew Knepley                MATOP_COPY=44,
1424ff1cced0SMatthew Knepley                MATOP_GET_ROW_MAX=45,
1425d643ce63SMatthew Knepley                MATOP_SCALE=46,
1426d643ce63SMatthew Knepley                MATOP_SHIFT=47,
1427d643ce63SMatthew Knepley                MATOP_DIAGONAL_SHIFT=48,
1428d643ce63SMatthew Knepley                MATOP_ILUDT_FACTOR=49,
1429ff1cced0SMatthew Knepley                MATOP_SET_BLOCK_SIZE=50,
1430d643ce63SMatthew Knepley                MATOP_GET_ROW_IJ=51,
1431d643ce63SMatthew Knepley                MATOP_RESTORE_ROW_IJ=52,
1432d643ce63SMatthew Knepley                MATOP_GET_COLUMN_IJ=53,
1433d643ce63SMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=54,
1434d643ce63SMatthew Knepley                MATOP_FDCOLORING_CREATE=55,
1435d643ce63SMatthew Knepley                MATOP_COLORING_PATCH=56,
1436d643ce63SMatthew Knepley                MATOP_SET_UNFACTORED=57,
1437d643ce63SMatthew Knepley                MATOP_PERMUTE=58,
1438d643ce63SMatthew Knepley                MATOP_SET_VALUES_BLOCKED=59,
1439d643ce63SMatthew Knepley                MATOP_GET_SUBMATRIX=60,
1440d643ce63SMatthew Knepley                MATOP_DESTROY=61,
1441d643ce63SMatthew Knepley                MATOP_VIEW=62,
1442ff1cced0SMatthew Knepley                MATOP_CONVERT_FROM=63,
1443d643ce63SMatthew Knepley                MATOP_USE_SCALED_FORM=64,
1444d643ce63SMatthew Knepley                MATOP_SCALE_SYSTEM=65,
1445d643ce63SMatthew Knepley                MATOP_UNSCALE_SYSTEM=66,
1446ba16a8cbSKris Buschelman                MATOP_SET_LOCAL_TO_GLOBAL_MAP=67,
1447d643ce63SMatthew Knepley                MATOP_SET_VALUES_LOCAL=68,
1448d643ce63SMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=69,
1449ff1cced0SMatthew Knepley                MATOP_GET_ROW_MAX_ABS=70,
1450c87e5d42SMatthew Knepley                MATOP_GET_ROW_MIN_ABS=71,
1451c87e5d42SMatthew Knepley                MATOP_CONVERT=72,
1452c87e5d42SMatthew Knepley                MATOP_SET_COLORING=73,
1453c87e5d42SMatthew Knepley                MATOP_SET_VALUES_ADIC=74,
1454c87e5d42SMatthew Knepley                MATOP_SET_VALUES_ADIFOR=75,
1455c87e5d42SMatthew Knepley                MATOP_FD_COLORING_APPLY=76,
1456c87e5d42SMatthew Knepley                MATOP_SET_FROM_OPTIONS=77,
1457c87e5d42SMatthew Knepley                MATOP_MULT_CON=78,
1458c87e5d42SMatthew Knepley                MATOP_MULT_TRANSPOSE_CON=79,
1459d643ce63SMatthew Knepley                MATOP_PERMUTE_SPARSIFY=80,
1460d643ce63SMatthew Knepley                MATOP_MULT_MULTIPLE=81,
146141acf15aSKris Buschelman                MATOP_SOLVE_MULTIPLE=82,
146241acf15aSKris Buschelman                MATOP_GET_INERTIA=83,
146341acf15aSKris Buschelman                MATOP_LOAD=84,
146441acf15aSKris Buschelman                MATOP_IS_SYMMETRIC=85,
146541acf15aSKris Buschelman                MATOP_IS_HERMITIAN=86,
146641acf15aSKris Buschelman                MATOP_IS_STRUCTURALLY_SYMMETRIC=87,
146741acf15aSKris Buschelman                MATOP_PB_RELAX=88,
146841acf15aSKris Buschelman                MATOP_GET_VECS=89,
146941acf15aSKris Buschelman                MATOP_MAT_MULT=90,
147041acf15aSKris Buschelman                MATOP_MAT_MULT_SYMBOLIC=91,
147141acf15aSKris Buschelman                MATOP_MAT_MULT_NUMERIC=92,
147241acf15aSKris Buschelman                MATOP_PTAP=93,
147341acf15aSKris Buschelman                MATOP_PTAP_SYMBOLIC=94,
1474bc011b1eSHong Zhang                MATOP_PTAP_NUMERIC=95,
1475bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE=96,
1476bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97,
14777ba1a0bfSKris Buschelman                MATOP_MAT_MULTTRANSPOSE_NUMERIC=98,
14787ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_SEQAIJ=99,
14797ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_SEQAIJ=100,
14807ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_MPIAIJ=101,
148187d4246cSBarry Smith                MATOP_PTAP_NUMERIC_MPIAIJ=102,
1482ff1cced0SMatthew Knepley                MATOP_CONJUGATE=103,
1483ff1cced0SMatthew Knepley                MATOP_SET_SIZES=104,
1484f5edf698SHong Zhang                MATOP_SET_VALUES_ROW = 105,
1485ff1cced0SMatthew Knepley                MATOP_REAL_PART=106,
1486ff1cced0SMatthew Knepley                MATOP_IMAG_PART=107,
1487d5c63f83SSatish Balay                MATOP_GET_ROW_UTRIANGULAR=108,
14882bebee5dSHong Zhang                MATOP_RESTORE_ROW_UTRIANGULAR=109,
148969db28dcSHong Zhang                MATOP_MATSOLVE=110,
1490829201f2SHong Zhang                MATOP_GET_REDUNDANTMATRIX=111,
1491ff1cced0SMatthew Knepley                MATOP_GET_ROW_MIN=112,
1492ff1cced0SMatthew Knepley                MATOP_GET_COLUMN_VEC=113,
1493ff1cced0SMatthew Knepley                MATOP_MISSING_DIAGONAL=114,
1494ff1cced0SMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=115,
1495ff1cced0SMatthew Knepley                MATOP_DESTROY_SOLVER=116
1496fae171e0SBarry Smith              } MatOperation;
1497be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1499be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1501112a2221SBarry Smith 
150290ace30eSBarry Smith /*
150390ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
150490ace30eSBarry Smith  stored in a universal format. By changing the format with
15057973ad04SBarry Smith  PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
150690ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
150790ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
150890ace30eSBarry Smith  read into matrices of the same time.
150990ace30eSBarry Smith */
151090ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
151190ace30eSBarry Smith 
1512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1513be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
15141f4e1ec7SBarry Smith 
1515d9274352SBarry Smith /*S
1516d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1517d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1518d9274352SBarry Smith 
1519f7a9e4ceSBarry Smith    Level: advanced
1520d9274352SBarry Smith 
1521d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1522d9274352SBarry Smith 
15236e1639daSBarry Smith   Users manual sections:
15246e1639daSBarry Smith .   sec_singular
15256e1639daSBarry Smith 
1526d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1527d9274352SBarry Smith S*/
152874637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1529d9274352SBarry Smith 
1530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1531281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*);
1532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
153595902228SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat,PetscTruth *);
153674637425SBarry Smith 
1537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1538be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
154046129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
15413f1d51d7SBarry Smith 
1542be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1543be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1544be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1545c4f061fbSSatish Balay 
1546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1547b0a32e0cSBarry Smith 
1548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
154904f1ad80SBarry Smith 
1550fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
15513ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec);
15523ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1553e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1554e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1555e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace);
1556e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1557e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat);
1558e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal);
1559e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt);
1560e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *);
1561e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat);
1562e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1563e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1564e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal);
1565e884886fSBarry Smith 
1566e884886fSBarry Smith 
1567e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1568e884886fSBarry Smith 
1569e884886fSBarry Smith /*E
1570e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1571e884886fSBarry Smith 
1572e884886fSBarry Smith    Level: beginner
1573e884886fSBarry Smith 
1574e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1575e884886fSBarry Smith E*/
1576a313700dSBarry Smith #define MatMFFDType char*
1577e884886fSBarry Smith #define MATMFFD_DS  "ds"
1578e884886fSBarry Smith #define MATMFFD_WP  "wp"
1579e884886fSBarry Smith 
1580a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,const MatMFFDType);
1581e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1582e884886fSBarry Smith 
1583e884886fSBarry Smith /*MC
1584e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1585e884886fSBarry Smith 
1586e884886fSBarry Smith    Synopsis:
1587e884886fSBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1588e884886fSBarry Smith 
1589e884886fSBarry Smith    Not Collective
1590e884886fSBarry Smith 
1591e884886fSBarry Smith    Input Parameters:
1592e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1593e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1594e884886fSBarry Smith .  name_create - name of routine to create method context
1595e884886fSBarry Smith -  routine_create - routine to create method context
1596e884886fSBarry Smith 
1597e884886fSBarry Smith    Level: developer
1598e884886fSBarry Smith 
1599e884886fSBarry Smith    Notes:
1600e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1601e884886fSBarry Smith 
1602e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1603e884886fSBarry Smith    is ignored.
1604e884886fSBarry Smith 
1605e884886fSBarry Smith    Sample usage:
1606e884886fSBarry Smith .vb
1607e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1608e884886fSBarry Smith                "MyHCreate",MyHCreate);
1609e884886fSBarry Smith .ve
1610e884886fSBarry Smith 
1611e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1612e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1613e884886fSBarry Smith    or at runtime via the option
1614e884886fSBarry Smith $     -snes_mf_type my_h
1615e884886fSBarry Smith 
1616e884886fSBarry Smith .keywords: MatMFFD, register
1617e884886fSBarry Smith 
1618e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1619e884886fSBarry Smith M*/
1620e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1621e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1622e884886fSBarry Smith #else
1623e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1624e884886fSBarry Smith #endif
1625e884886fSBarry Smith 
1626e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]);
1627e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void);
1628e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal);
1629e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth);
1630e884886fSBarry Smith 
1631e884886fSBarry Smith 
1632be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1633be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
16347dbadf16SMatthew Knepley 
1635e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
16362eac72dbSBarry Smith #endif
1637