xref: /petsc/include/petscmat.h (revision 85e3dda7fddc385e91bc763ba6bbff3d650f88ca)
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"
64*85e3dda7SBarry Smith #define MATTRANSPOSE       "transpose"
65d91e6319SBarry Smith 
669989ab13SBarry Smith /*E
67c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
689989ab13SBarry Smith 
699989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
709989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
719989ab13SBarry Smith 
729989ab13SBarry Smith 
739989ab13SBarry Smith    Level: beginner
749989ab13SBarry Smith 
755c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
769989ab13SBarry Smith E*/
77c7393fdbSBarry Smith #define MatSolverPackage char*
78b24902e0SBarry Smith #define MAT_SOLVER_SPOOLES      "spooles"
79b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU      "superlu"
80b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU_DIST "superlu_dist"
81b24902e0SBarry Smith #define MAT_SOLVER_UMFPACK      "umfpack"
82b24902e0SBarry Smith #define MAT_SOLVER_ESSL         "essl"
83b24902e0SBarry Smith #define MAT_SOLVER_LUSOL        "lusol"
84b24902e0SBarry Smith #define MAT_SOLVER_MUMPS        "mumps"
85b24902e0SBarry Smith #define MAT_SOLVER_DSCPACK      "dscpack"
86b24902e0SBarry Smith #define MAT_SOLVER_MATLAB       "matlab"
8743244d56SBarry Smith #define MAT_SOLVER_PETSC        "petsc"
88b24902e0SBarry Smith 
89b24902e0SBarry Smith /*E
90b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
91b24902e0SBarry Smith 
92b24902e0SBarry Smith     Level: beginner
93b24902e0SBarry Smith 
94b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
95b24902e0SBarry Smith 
96c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
97b24902e0SBarry Smith E*/
985c9eb25fSBarry Smith typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC} MatFactorType;
99c7393fdbSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
100b24902e0SBarry Smith 
1019989ab13SBarry Smith 
102c06d978dSMatthew Knepley /* Logging support */
103552e946dSBarry Smith #define    MAT_FILE_COOKIE 1211216    /* used to indicate matrices in binary files */
104be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE;
105be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
106be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE;
107be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE;
108e884886fSBarry Smith extern PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE;
109c06d978dSMatthew Knepley 
110ceb03754SKris Buschelman /*E
111ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
112ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
113ceb03754SKris Buschelman 
114ceb03754SKris Buschelman     Level: beginner
115ceb03754SKris Buschelman 
116ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
117ceb03754SKris Buschelman 
1180c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
119ceb03754SKris Buschelman E*/
120ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse;
121ceb03754SKris Buschelman 
1225494a064SHong Zhang /*E
1235494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
124829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1255494a064SHong Zhang 
1265494a064SHong Zhang     Level: beginner
1275494a064SHong Zhang 
128829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1295494a064SHong Zhang E*/
1305494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1315494a064SHong Zhang 
132e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]);
133c06d978dSMatthew Knepley 
134f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
135a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
136a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
137f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
138a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,const MatType);
139be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
140be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
141be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
142be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
14330de9b25SBarry Smith 
144f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
145f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
146f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
147f69a0ea3SMatthew Knepley 
14830de9b25SBarry Smith /*MC
14930de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
15030de9b25SBarry Smith 
15130de9b25SBarry Smith    Synopsis:
152c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
15330de9b25SBarry Smith 
15430de9b25SBarry Smith    Not Collective
15530de9b25SBarry Smith 
15630de9b25SBarry Smith    Input Parameters:
15730de9b25SBarry Smith +  name - name of a new user-defined matrix type
15830de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
15930de9b25SBarry Smith .  name_create - name of routine to create method context
16030de9b25SBarry Smith -  routine_create - routine to create method context
16130de9b25SBarry Smith 
16230de9b25SBarry Smith    Notes:
16330de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
16430de9b25SBarry Smith 
16530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
16630de9b25SBarry Smith    is ignored.
16730de9b25SBarry Smith 
16830de9b25SBarry Smith    Sample usage:
16930de9b25SBarry Smith .vb
17030de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
17130de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
17230de9b25SBarry Smith .ve
17330de9b25SBarry Smith 
17430de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
17530de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
17630de9b25SBarry Smith    or at runtime via the option
17730de9b25SBarry Smith $     -mat_type my_mat
17830de9b25SBarry Smith 
17930de9b25SBarry Smith    Level: advanced
18030de9b25SBarry Smith 
181ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
18230de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
18330de9b25SBarry Smith 
18430de9b25SBarry Smith .keywords: Mat, register
18530de9b25SBarry Smith 
18630de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
18730de9b25SBarry Smith 
18830de9b25SBarry Smith M*/
189273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
190273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
191273d9f13SBarry Smith #else
192273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
19330de9b25SBarry Smith #endif
19430de9b25SBarry Smith 
1959c666560SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
196273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
197b0a32e0cSBarry Smith extern PetscFList MatList;
19828988994SBarry Smith 
199be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
200be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
201be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
202ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
203ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
204ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
205ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
206ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
207ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
208ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
209be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
210ba966639SSatish 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)
211ba966639SSatish 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)
212ba966639SSatish 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)
213ba966639SSatish 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)
214ba966639SSatish 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))
215ba966639SSatish 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))
216ba966639SSatish 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))
217ba966639SSatish 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)
218ba966639SSatish 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)
219ba966639SSatish 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)
220ba966639SSatish 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)
221ba966639SSatish 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))
222ba966639SSatish 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))
223ba966639SSatish 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))
22463ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2258d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2268d7a6e47SBarry Smith 
227be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
228be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
229ba966639SSatish 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)
230ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
231ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
232ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
233ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
234ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
235ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
236be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
237ba966639SSatish 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)
238ba966639SSatish 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)
239ba966639SSatish 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)
240ba966639SSatish 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)
241ba966639SSatish 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))
242ba966639SSatish 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))
243ba966639SSatish 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))
244ba966639SSatish 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)
245ba966639SSatish 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)
246ba966639SSatish 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)
247ba966639SSatish 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)
248ba966639SSatish 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))
249ba966639SSatish 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))
250ba966639SSatish 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))
251be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
252be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
253ba966639SSatish 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)
254ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
255ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
256ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
257ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
258ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
259ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
260be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
261ba966639SSatish 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)
262ba966639SSatish 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)
263ba966639SSatish 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)
264ba966639SSatish 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)
265ba966639SSatish 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))
266ba966639SSatish 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))
267ba966639SSatish 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))
268ba966639SSatish 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)
269ba966639SSatish 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)
270ba966639SSatish 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)
271ba966639SSatish 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)
272ba966639SSatish 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))
273ba966639SSatish 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))
274ba966639SSatish 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))
275be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
276ba966639SSatish 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)
277ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
278be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
279be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
280ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
281ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*);
282284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
2831472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
2841472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
2852a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*);
2862a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter);
2872a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*);
2888cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
289793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat);
290b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat);
291793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2925a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*);
2931472f72bSBarry Smith 
294f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
295be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
29621c89e3eSBarry Smith 
2970c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
29899cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat);
29999cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat);
30099cafbc1SBarry Smith 
3018ed539a5SBarry Smith /* ------------------------------------------------------------*/
302be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
303be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
30487d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
305f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
30684cb2905SBarry Smith 
3072ef4de8bSBarry Smith /*S
3082ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
3092ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3102ef4de8bSBarry Smith 
3112ef4de8bSBarry Smith    Level: beginner
3122ef4de8bSBarry Smith 
3132ef4de8bSBarry Smith   Concepts: matrix; linear operator
3142ef4de8bSBarry Smith 
315d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3162ef4de8bSBarry Smith S*/
317435da068SBarry Smith typedef struct {
318c1ac3661SBarry Smith   PetscInt k,j,i,c;
319435da068SBarry Smith } MatStencil;
3202ef4de8bSBarry Smith 
321be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
322be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
323be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
324435da068SBarry Smith 
325be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
326be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
327be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
3283a7fca6bSBarry Smith 
329d91e6319SBarry Smith /*E
330d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
331d91e6319SBarry Smith      to continue to add values to it
332d91e6319SBarry Smith 
333d91e6319SBarry Smith     Level: beginner
334d91e6319SBarry Smith 
335d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
336d91e6319SBarry Smith E*/
3376d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
338be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
339be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
340be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
3414f9c727eSBarry Smith 
3421d73ed98SBarry Smith 
34330de9b25SBarry Smith 
344d91e6319SBarry Smith /*E
345d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
346d91e6319SBarry Smith 
347d91e6319SBarry Smith     Level: beginner
348d91e6319SBarry Smith 
3490a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
350d91e6319SBarry Smith 
351d91e6319SBarry Smith .seealso: MatSetOption()
352d91e6319SBarry Smith E*/
353512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
3544e0d8c25SBarry Smith               MAT_SYMMETRIC,
3554e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
3568047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
3574e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
3584e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
3594e0d8c25SBarry Smith               MAT_KEEP_ZEROED_ROWS,MAT_IGNORE_ZERO_ENTRIES,
3604e0d8c25SBarry Smith               MAT_USE_INODES,
3614e0d8c25SBarry Smith               MAT_HERMITIAN,
3624e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
3634e0d8c25SBarry Smith               MAT_USE_COMPRESSEDROW,
3644e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
3654e0d8c25SBarry Smith               MAT_GETROW_UPPERTRIANGULAR} MatOption;
366290bbb0aSBarry Smith extern const char *MatOptions[];
3674e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth);
368a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,const MatType*);
369a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
37084cb2905SBarry Smith 
371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
373be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
374f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat);
375f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat);
376be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
378be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
379be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
380ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
382be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
383ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
384be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
3857b80b807SBarry Smith 
3861620fd73SBarry Smith 
387be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
388be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
389ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
390be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
391be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
392ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
393e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t)
3941cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*);
395be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
396ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
398be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
3992bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat);
400f5edf698SHong Zhang 
401d91e6319SBarry Smith /*E
402d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
403d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
404d91e6319SBarry Smith 
405d91e6319SBarry Smith     Level: beginner
406d91e6319SBarry Smith 
407d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
408d91e6319SBarry Smith 
409d91e6319SBarry Smith .seealso: MatDuplicate()
410d91e6319SBarry Smith E*/
4112e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
4122e8a6d31SBarry Smith 
413a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,const MatType,MatReuse,Mat*);
414a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
415be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
416ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
417ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
41894a9d846SBarry Smith 
419d91e6319SBarry Smith /*E
420d91e6319SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
421d91e6319SBarry Smith 
422d91e6319SBarry Smith     Level: beginner
423d91e6319SBarry Smith 
424d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
425d91e6319SBarry Smith 
42694b7f48cSBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
427d91e6319SBarry Smith E*/
428c537a176SHong Zhang typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
429cb5b572fSBarry Smith 
430be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
431be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
432be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
433ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
434e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t)
435be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
436ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
4371cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*);
4381cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t)
439be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
440be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
44164c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *);
442a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,const MatType,Mat*);
443a313700dSBarry Smith PetscPolymorphicFunction(MatLoad,(PetscViewer v,const MatType t),(v,t,&a),Mat,a)
4447b80b807SBarry Smith 
4458f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4468f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
4478f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4488f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
449d4fbbf0eSBarry Smith 
450d91e6319SBarry Smith /*S
451d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
452d91e6319SBarry Smith 
453d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
454d91e6319SBarry Smith 
455d91e6319SBarry Smith    Level: intermediate
456d91e6319SBarry Smith 
457d91e6319SBarry Smith   Concepts: matrix^nonzero information
458d91e6319SBarry Smith 
459d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
460d91e6319SBarry Smith S*/
4614e220ebcSLois Curfman McInnes typedef struct {
462b0a32e0cSBarry Smith   PetscLogDouble rows_global,columns_global;         /* number of global rows and columns */
463b0a32e0cSBarry Smith   PetscLogDouble rows_local,columns_local;           /* number of local rows and columns */
464b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
465b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
466b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
467b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
468b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
469b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
470b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4714e220ebcSLois Curfman McInnes } MatInfo;
4724e220ebcSLois Curfman McInnes 
473d9274352SBarry Smith /*E
474d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
475d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
476d9274352SBarry Smith 
477d9274352SBarry Smith     Level: beginner
478d9274352SBarry Smith 
479d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
480d9274352SBarry Smith 
481d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
482d9274352SBarry Smith E*/
4837b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
484be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
485be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
486be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
487985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]);
488985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]);
489985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
49086697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec);
491fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*);
492fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
493be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
494ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
495be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
496be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
497be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
499ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
501be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
502be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
503be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
5047b80b807SBarry Smith 
505be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
506ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
508f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
509f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
510f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
511f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5127b80b807SBarry Smith 
513be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
515be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5165ef9f2a5SBarry Smith 
517be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
519be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5208ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**);
521f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
5227b80b807SBarry Smith 
523be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
524be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
525be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *);
526abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *);
527829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat *[]);
528829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat *[]);
5295494a064SHong Zhang 
530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
535be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
536be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
538dd6ea824SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,MatScalar**,Mat*);
53943eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
540cd116e26SSatish Balay #include "petscctable.h"
54143eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
54243eb5e2fSMatthew Knepley #else
54343eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
54443eb5e2fSMatthew Knepley #endif
5458efafbd8SBarry Smith 
546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5477b80b807SBarry Smith 
548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
55122440eb1SKris Buschelman 
552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
553be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
55522440eb1SKris Buschelman 
556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
557be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
558be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
559bc011b1eSHong Zhang 
560f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
56104aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure);
562be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat);
5631c741599SHong Zhang 
564f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
565f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5667b80b807SBarry Smith 
567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
568be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
569f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
570f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
571be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
572be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
573052efed2SBarry Smith 
574be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
575be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
57690f02eecSBarry Smith 
577be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
5780c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
579ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
580be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
581be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
58269db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
583bd481603SBarry Smith 
584bd481603SBarry Smith /*MC
5851620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5861620fd73SBarry Smith 
5871620fd73SBarry Smith    Not collective
5881620fd73SBarry Smith 
5891620fd73SBarry Smith    Input Parameters:
5901620fd73SBarry Smith +  m - the matrix
5911620fd73SBarry Smith .  row - the row location of the entry
5921620fd73SBarry Smith .  col - the column location of the entry
5931620fd73SBarry Smith .  value - the value to insert
5941620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5951620fd73SBarry Smith 
5961620fd73SBarry Smith    Notes:
5971620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
5981620fd73SBarry Smith    values simultaneously if possible.
5991620fd73SBarry Smith 
6001620fd73SBarry Smith    Level: beginner
6011620fd73SBarry Smith 
6021620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6031620fd73SBarry Smith M*/
6041620fd73SBarry 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);}
6051620fd73SBarry Smith 
6061620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6071620fd73SBarry Smith 
6081620fd73SBarry 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);}
6091620fd73SBarry Smith 
6101620fd73SBarry Smith /*MC
6110d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
612bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
613bd481603SBarry Smith 
614bd481603SBarry Smith    Synopsis:
615c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
616bd481603SBarry Smith 
617bd481603SBarry Smith    Collective on MPI_Comm
618bd481603SBarry Smith 
619bd481603SBarry Smith    Input Parameters:
620bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
621859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
622859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
623bd481603SBarry Smith 
624bd481603SBarry Smith    Output Parameters:
625bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
626bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
627bd481603SBarry Smith 
628bd481603SBarry Smith 
629bd481603SBarry Smith    Level: intermediate
630bd481603SBarry Smith 
631bd481603SBarry Smith    Notes:
632bd481603SBarry Smith    See the chapter in the users manual on performance for more details
633bd481603SBarry Smith 
6341d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
635bd481603SBarry Smith 
636bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
637bd481603SBarry Smith 
6381620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6391620fd73SBarry Smith 
640bd481603SBarry Smith   Concepts: preallocation^Matrix
641bd481603SBarry Smith 
642bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
643bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
644bd481603SBarry Smith M*/
645c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6467c922b88SBarry Smith { \
647a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
648a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
649a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
650a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
651c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
652a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6537c922b88SBarry Smith 
654bd481603SBarry Smith /*MC
6550d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
656bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
657bd481603SBarry Smith 
658bd481603SBarry Smith    Synopsis:
659c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
660bd481603SBarry Smith 
661bd481603SBarry Smith    Collective on MPI_Comm
662bd481603SBarry Smith 
663bd481603SBarry Smith    Input Parameters:
664bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
665859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
666859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
667bd481603SBarry Smith 
668bd481603SBarry Smith    Output Parameters:
669bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
670bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
671bd481603SBarry Smith 
672bd481603SBarry Smith 
673bd481603SBarry Smith    Level: intermediate
674bd481603SBarry Smith 
675bd481603SBarry Smith    Notes:
676bd481603SBarry Smith    See the chapter in the users manual on performance for more details
677bd481603SBarry Smith 
6781d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
679bd481603SBarry Smith 
6801620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6811620fd73SBarry Smith 
682bd481603SBarry Smith   Concepts: preallocation^Matrix
683bd481603SBarry Smith 
684bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
685bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
686bd481603SBarry Smith M*/
687222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
688222b16d4SBarry Smith { \
689a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
690a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
691a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
692a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
693c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
694a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
695222b16d4SBarry Smith 
696bd481603SBarry Smith /*MC
697bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
698bd481603SBarry Smith        inserted using a local number of the rows and columns
699bd481603SBarry Smith 
700bd481603SBarry Smith    Synopsis:
701c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
702bd481603SBarry Smith 
703bd481603SBarry Smith    Not Collective
704bd481603SBarry Smith 
705bd481603SBarry Smith    Input Parameters:
706bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
707bd481603SBarry Smith .  nrows - the number of rows indicated
7081d73ed98SBarry Smith .  rows - the indices of the rows
709bd481603SBarry Smith .  ncols - the number of columns in the matrix
710bd481603SBarry Smith .  cols - the columns indicated
711bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
712bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
713bd481603SBarry Smith 
714bd481603SBarry Smith 
715bd481603SBarry Smith    Level: intermediate
716bd481603SBarry Smith 
717bd481603SBarry Smith    Notes:
718bd481603SBarry Smith    See the chapter in the users manual on performance for more details
719bd481603SBarry Smith 
7201d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
721bd481603SBarry Smith 
722bd481603SBarry Smith   Concepts: preallocation^Matrix
723bd481603SBarry Smith 
724bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
725bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
726bd481603SBarry Smith M*/
727c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
728c4f061fbSSatish Balay {\
729c1ac3661SBarry Smith   PetscInt __l;\
730ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
731ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
732c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
733ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
734c4f061fbSSatish Balay   }\
735c4f061fbSSatish Balay }
736c4f061fbSSatish Balay 
737bd481603SBarry Smith /*MC
738bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
739bd481603SBarry Smith        inserted using a local number of the rows and columns
740bd481603SBarry Smith 
741bd481603SBarry Smith    Synopsis:
742c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
743bd481603SBarry Smith 
744bd481603SBarry Smith    Not Collective
745bd481603SBarry Smith 
746bd481603SBarry Smith    Input Parameters:
747bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
748bd481603SBarry Smith .  nrows - the number of rows indicated
7491d73ed98SBarry Smith .  rows - the indices of the rows
750bd481603SBarry Smith .  ncols - the number of columns in the matrix
751bd481603SBarry Smith .  cols - the columns indicated
752bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
753bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
754bd481603SBarry Smith 
755bd481603SBarry Smith 
756bd481603SBarry Smith    Level: intermediate
757bd481603SBarry Smith 
758bd481603SBarry Smith    Notes:
759bd481603SBarry Smith    See the chapter in the users manual on performance for more details
760bd481603SBarry Smith 
761bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
762bd481603SBarry Smith 
763bd481603SBarry Smith   Concepts: preallocation^Matrix
764bd481603SBarry Smith 
765bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
766bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
767bd481603SBarry Smith M*/
768d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
769d3d32019SBarry Smith {\
770c1ac3661SBarry Smith   PetscInt __l;\
771d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
772d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
773d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
774d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
775d3d32019SBarry Smith   }\
776d3d32019SBarry Smith }
777d3d32019SBarry Smith 
778bd481603SBarry Smith /*MC
779bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
780bd481603SBarry Smith        inserted using a local number of the rows and columns
781bd481603SBarry Smith 
782bd481603SBarry Smith    Synopsis:
783c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
784bd481603SBarry Smith 
785bd481603SBarry Smith    Not Collective
786bd481603SBarry Smith 
787bd481603SBarry Smith    Input Parameters:
78864075487SBarry Smith +  row - the row
789bd481603SBarry Smith .  ncols - the number of columns in the matrix
790a50f8bf6SHong Zhang -  cols - the columns indicated
791a50f8bf6SHong Zhang 
792a50f8bf6SHong Zhang    Output Parameters:
793a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
794bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
795bd481603SBarry Smith 
796bd481603SBarry Smith 
797bd481603SBarry Smith    Level: intermediate
798bd481603SBarry Smith 
799bd481603SBarry Smith    Notes:
800bd481603SBarry Smith    See the chapter in the users manual on performance for more details
801bd481603SBarry Smith 
802bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
803bd481603SBarry Smith 
8041620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8051620fd73SBarry Smith 
806bd481603SBarry Smith   Concepts: preallocation^Matrix
807bd481603SBarry Smith 
808bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
809bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
810bd481603SBarry Smith M*/
811c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
812c1ac3661SBarry Smith { PetscInt __i; \
8138ee2e534SBarry Smith   if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
814a89bc211SBarry 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);\
8157c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
81664075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8177cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8187c922b88SBarry Smith   }\
8197c922b88SBarry Smith }
8207c922b88SBarry Smith 
821bd481603SBarry Smith /*MC
822bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
823bd481603SBarry Smith        inserted using a local number of the rows and columns
824bd481603SBarry Smith 
825bd481603SBarry Smith    Synopsis:
826c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
827bd481603SBarry Smith 
828bd481603SBarry Smith    Not Collective
829bd481603SBarry Smith 
830bd481603SBarry Smith    Input Parameters:
831bd481603SBarry Smith +  nrows - the number of rows indicated
8321d73ed98SBarry Smith .  rows - the indices of the rows
833bd481603SBarry Smith .  ncols - the number of columns in the matrix
834bd481603SBarry Smith .  cols - the columns indicated
835bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
836bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
837bd481603SBarry Smith 
838bd481603SBarry Smith 
839bd481603SBarry Smith    Level: intermediate
840bd481603SBarry Smith 
841bd481603SBarry Smith    Notes:
842bd481603SBarry Smith    See the chapter in the users manual on performance for more details
843bd481603SBarry Smith 
844bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
845bd481603SBarry Smith 
8461620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8471620fd73SBarry Smith 
848bd481603SBarry Smith   Concepts: preallocation^Matrix
849bd481603SBarry Smith 
850bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
851bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
852bd481603SBarry Smith M*/
853d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
854c1ac3661SBarry Smith { PetscInt __i; \
855d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
856d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
857d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
858d3d32019SBarry Smith   }\
859d3d32019SBarry Smith }
860d3d32019SBarry Smith 
861bd481603SBarry Smith /*MC
86216371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
86316371a99SBarry Smith 
86416371a99SBarry Smith    Synopsis:
86516371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
86616371a99SBarry Smith 
86716371a99SBarry Smith    Not Collective
86816371a99SBarry Smith 
86916371a99SBarry Smith    Input Parameters:
87016371a99SBarry Smith .  A - matrix
87116371a99SBarry Smith .  row - row where values exist (must be local to this process)
87216371a99SBarry Smith .  ncols - number of columns
87316371a99SBarry Smith .  cols - columns with nonzeros
87416371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
87516371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
87616371a99SBarry Smith 
87716371a99SBarry Smith 
87816371a99SBarry Smith    Level: intermediate
87916371a99SBarry Smith 
88016371a99SBarry Smith    Notes:
88116371a99SBarry Smith    See the chapter in the users manual on performance for more details
88216371a99SBarry Smith 
88316371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
88416371a99SBarry Smith 
88516371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
88616371a99SBarry Smith 
88716371a99SBarry Smith   Concepts: preallocation^Matrix
88816371a99SBarry Smith 
88916371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
89016371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
89116371a99SBarry Smith M*/
89216371a99SBarry 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);}
89316371a99SBarry Smith 
89416371a99SBarry Smith 
89516371a99SBarry Smith /*MC
8960d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
897bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
898bd481603SBarry Smith 
899bd481603SBarry Smith    Synopsis:
900c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
901bd481603SBarry Smith 
902bd481603SBarry Smith    Collective on MPI_Comm
903bd481603SBarry Smith 
904bd481603SBarry Smith    Input Parameters:
90516371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
906bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
907bd481603SBarry Smith 
908bd481603SBarry Smith 
909bd481603SBarry Smith    Level: intermediate
910bd481603SBarry Smith 
911bd481603SBarry Smith    Notes:
912bd481603SBarry Smith    See the chapter in the users manual on performance for more details
913bd481603SBarry Smith 
914bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
915bd481603SBarry Smith 
9161620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9171620fd73SBarry Smith 
918bd481603SBarry Smith   Concepts: preallocation^Matrix
919bd481603SBarry Smith 
920bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
921bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
922bd481603SBarry Smith M*/
923a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9247c922b88SBarry Smith 
925bd481603SBarry Smith 
926bd481603SBarry Smith 
9277b80b807SBarry Smith /* Routines unique to particular data structures */
928be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
929ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
930698d4c6aSKris Buschelman 
931be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
932be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
933698d4c6aSKris Buschelman 
934be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
935be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
936be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
937c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
938c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9397b80b807SBarry Smith 
940a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
941a96a251dSBarry Smith 
942be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
943ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
944be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
945ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
946be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
947ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
948be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDensePreallocation(Mat,PetscScalar[]);
949273d9f13SBarry Smith 
950be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
951ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
952be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
953be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
95453803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
955725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
956be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
957be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
958be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
959be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
960284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
961be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
962be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
963be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
964be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
965273d9f13SBarry Smith 
966be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
967ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*);
9681b807ce4Svictorle 
969be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
970be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
9712e8a6d31SBarry Smith 
972be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
9733a7fca6bSBarry Smith 
9747b80b807SBarry Smith /*
9757b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
97694b7f48cSBarry Smith   done through the KSP and PC interfaces.
9777b80b807SBarry Smith */
9787b80b807SBarry Smith 
979d9274352SBarry Smith /*E
980d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
981d9274352SBarry Smith        with an optional dynamic library name, for example
982d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
983d9274352SBarry Smith 
984d9274352SBarry Smith    Level: beginner
985d9274352SBarry Smith 
986e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
987e5a9bf91SBarry Smith 
988d9274352SBarry Smith .seealso: MatGetOrdering()
989d9274352SBarry Smith E*/
9903eaad2caSSatish Balay #define MatOrderingType char*
991b12f92e5SBarry Smith #define MATORDERING_NATURAL   "natural"
992b12f92e5SBarry Smith #define MATORDERING_ND        "nd"
993b12f92e5SBarry Smith #define MATORDERING_1WD       "1wd"
994b12f92e5SBarry Smith #define MATORDERING_RCM       "rcm"
995b12f92e5SBarry Smith #define MATORDERING_QMD       "qmd"
996b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength"
99762152c8bSBarry Smith #define MATORDERING_DSC_ND    "dsc_nd"
99862152c8bSBarry Smith #define MATORDERING_DSC_MMD   "dsc_mmd"
99962152c8bSBarry Smith #define MATORDERING_DSC_MDF   "dsc_mdf"
1000c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
1001c06d978dSMatthew Knepley #define MATORDERING_IDENTITY  "identity"
1002c06d978dSMatthew Knepley #define MATORDERING_REVERSE   "reverse"
1003b12f92e5SBarry Smith 
1004f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
1005be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
1006f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
100730de9b25SBarry Smith 
100830de9b25SBarry Smith /*MC
100930de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
101030de9b25SBarry Smith                                matrix package.
101130de9b25SBarry Smith 
101230de9b25SBarry Smith    Synopsis:
1013c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
101430de9b25SBarry Smith 
101530de9b25SBarry Smith    Not Collective
101630de9b25SBarry Smith 
101730de9b25SBarry Smith    Input Parameters:
101830de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
101930de9b25SBarry Smith .  path - location of library where creation routine is
102030de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
102130de9b25SBarry Smith -  function - function pointer that creates the ordering
102230de9b25SBarry Smith 
102330de9b25SBarry Smith    Level: developer
102430de9b25SBarry Smith 
102530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
102630de9b25SBarry Smith    is ignored.
102730de9b25SBarry Smith 
102830de9b25SBarry Smith    Sample usage:
102930de9b25SBarry Smith .vb
103030de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
103130de9b25SBarry Smith                "MyOrder",MyOrder);
103230de9b25SBarry Smith .ve
103330de9b25SBarry Smith 
103430de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
103530de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
103630de9b25SBarry Smith    or at runtime via the option
10372401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
103830de9b25SBarry Smith 
1039ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
104030de9b25SBarry Smith 
104130de9b25SBarry Smith .keywords: matrix, ordering, register
104230de9b25SBarry Smith 
104330de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
104430de9b25SBarry Smith M*/
1045aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1046f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1047b12f92e5SBarry Smith #else
1048f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1049b12f92e5SBarry Smith #endif
105030de9b25SBarry Smith 
1051be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
1052be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
10532bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
1054b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1055d4fbbf0eSBarry Smith 
1056be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1057a2ce50c7SBarry Smith 
1058d91e6319SBarry Smith /*S
10592401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1060b00f7748SHong Zhang 
106161cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
106261cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1063b00f7748SHong Zhang 
106415e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1065b00f7748SHong Zhang 
106661cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
106761cad860SBarry Smith 
1068b00f7748SHong Zhang    Level: developer
1069b00f7748SHong Zhang 
1070d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1071d7d82daaSBarry Smith           MatFactorInfoInitialize()
1072b00f7748SHong Zhang 
1073b00f7748SHong Zhang S*/
1074b00f7748SHong Zhang typedef struct {
10750a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
1076fbf22428SSatish Balay   PetscReal     shiftpd;        /* if true, shift until positive pivots */
10772cea7109SBarry Smith   PetscReal     shift_fraction; /* record shift fraction taken */
107815e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
107915e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1080b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
108115e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
108267e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1083348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1084bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1085bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
108662bba022SBarry Smith   PetscReal     shiftinblocks;  /* if block in block factorization has zero pivot then shift diagonal until non-singular */
108715e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
108815e8a5b3SHong Zhang } MatFactorInfo;
1089ffa6d0a5SLois Curfman McInnes 
1090be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
1091be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,MatFactorInfo*);
1092be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1093be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,MatFactorInfo*,Mat*);
1094be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,MatFactorInfo*);
1095be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,MatFactorInfo*);
1096be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1097be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1098be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1099be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,MatFactorInfo*);
1100be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,MatFactorInfo*,Mat*);
1101be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,MatFactorInfo*,Mat *);
1102be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1103be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1104be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1105be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1106be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1107be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1108be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1109be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
11108ed539a5SBarry Smith 
1111be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1112bb5a7306SBarry Smith 
1113d91e6319SBarry Smith /*E
1114d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1115bb1eb677SSatish Balay 
1116d91e6319SBarry Smith     Level: beginner
1117d91e6319SBarry Smith 
1118d9274352SBarry Smith    May be bitwise ORd together
1119d9274352SBarry Smith 
1120d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1121d91e6319SBarry Smith 
11224e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11234e7234bfSBarry Smith 
1124a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax()
1125d91e6319SBarry Smith E*/
1126ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1127ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1128ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
112984cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1130be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1131be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11328ed539a5SBarry Smith 
1133d4fbbf0eSBarry Smith /*
1134639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1135639f9d9dSBarry Smith */
1136b12f92e5SBarry Smith 
1137d9274352SBarry Smith /*E
1138d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1139d9274352SBarry Smith        with an optional dynamic library name, for example
1140d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1141d9274352SBarry Smith 
1142d9274352SBarry Smith    Level: beginner
1143d9274352SBarry Smith 
1144d9274352SBarry Smith .seealso: MatGetColoring()
1145d9274352SBarry Smith E*/
1146a313700dSBarry Smith #define MatColoringType char*
1147b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1148b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1149b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1150b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1151b12f92e5SBarry Smith 
1152ab970d9bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*);
11532e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
115430de9b25SBarry Smith 
115530de9b25SBarry Smith /*MC
115630de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
115730de9b25SBarry Smith                                matrix package.
115830de9b25SBarry Smith 
115930de9b25SBarry Smith    Synopsis:
1160c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
116130de9b25SBarry Smith 
116230de9b25SBarry Smith    Not Collective
116330de9b25SBarry Smith 
116430de9b25SBarry Smith    Input Parameters:
116530de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
116630de9b25SBarry Smith .  path - location of library where creation routine is
116730de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
116830de9b25SBarry Smith -  function - function pointer that creates the coloring
116930de9b25SBarry Smith 
117030de9b25SBarry Smith    Level: developer
117130de9b25SBarry Smith 
117230de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
117330de9b25SBarry Smith    is ignored.
117430de9b25SBarry Smith 
117530de9b25SBarry Smith    Sample usage:
117630de9b25SBarry Smith .vb
117730de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
117830de9b25SBarry Smith                "MyColor",MyColor);
117930de9b25SBarry Smith .ve
118030de9b25SBarry Smith 
118130de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
118230de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
118330de9b25SBarry Smith    or at runtime via the option
118430de9b25SBarry Smith $     -mat_coloring_type my_color
118530de9b25SBarry Smith 
1186ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
118730de9b25SBarry Smith 
118830de9b25SBarry Smith .keywords: matrix, Coloring, register
118930de9b25SBarry Smith 
119030de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
119130de9b25SBarry Smith M*/
1192aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1193f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1194b12f92e5SBarry Smith #else
1195f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1196b12f92e5SBarry Smith #endif
119730de9b25SBarry Smith 
11982bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1199f1144a30SSatish Balay 
1200be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1201be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1202be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1203639f9d9dSBarry Smith 
1204d9274352SBarry Smith /*S
1205d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1206d9274352SBarry Smith         and coloring
1207639f9d9dSBarry Smith 
1208d9274352SBarry Smith    Level: beginner
1209d9274352SBarry Smith 
1210d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1211d9274352SBarry Smith 
1212d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1213d9274352SBarry Smith S*/
1214e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1215639f9d9dSBarry Smith 
1216be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1217be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1218be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1219be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
12204a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1221be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1222be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring,PetscInt);
1223be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring,PetscInt*);
1224be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1225be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1226be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1227be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring);
1228be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1229be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1230639f9d9dSBarry Smith /*
12310752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12323eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12330752156aSBarry Smith */
1234ca161407SBarry Smith 
1235d9274352SBarry Smith /*S
1236d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1237d9274352SBarry Smith 
1238d9274352SBarry Smith    Level: beginner
1239d9274352SBarry Smith 
1240d9274352SBarry Smith   Concepts: partitioning
1241d9274352SBarry Smith 
1242743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1243d9274352SBarry Smith S*/
124491e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1245d9274352SBarry Smith 
1246d9274352SBarry Smith /*E
12475bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1248d9274352SBarry Smith        with an optional dynamic library name, for example
1249d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1250d9274352SBarry Smith 
1251d9274352SBarry Smith    Level: beginner
1252d9274352SBarry Smith 
1253b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1254d9274352SBarry Smith E*/
1255a313700dSBarry Smith #define MatPartitioningType char*
12568ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1257c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
12588ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
125971306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
126071306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
126171306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
126271306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
126371306c60Sjroman 
1264ca161407SBarry Smith 
1265be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1266a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1267be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1269be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1270be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1271be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1272be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
12732aabb6bbSBarry Smith 
1274be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
127530de9b25SBarry Smith 
127630de9b25SBarry Smith /*MC
127730de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
127830de9b25SBarry Smith    matrix package.
127930de9b25SBarry Smith 
128030de9b25SBarry Smith    Synopsis:
1281c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
128230de9b25SBarry Smith 
128330de9b25SBarry Smith    Not Collective
128430de9b25SBarry Smith 
128530de9b25SBarry Smith    Input Parameters:
128630de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
128730de9b25SBarry Smith .  path - location of library where creation routine is
128830de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
128930de9b25SBarry Smith -  function - function pointer that creates the partitioning type
129030de9b25SBarry Smith 
129130de9b25SBarry Smith    Level: developer
129230de9b25SBarry Smith 
129330de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
129430de9b25SBarry Smith    is ignored.
129530de9b25SBarry Smith 
129630de9b25SBarry Smith    Sample usage:
129730de9b25SBarry Smith .vb
129830de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
129930de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
130030de9b25SBarry Smith .ve
130130de9b25SBarry Smith 
130230de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
130330de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
130430de9b25SBarry Smith    or at runtime via the option
130530de9b25SBarry Smith $     -mat_partitioning_type my_part
130630de9b25SBarry Smith 
1307ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
130830de9b25SBarry Smith 
130930de9b25SBarry Smith .keywords: matrix, partitioning, register
131030de9b25SBarry Smith 
131130de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
131230de9b25SBarry Smith M*/
1313aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1314f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13152aabb6bbSBarry Smith #else
1316f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13172aabb6bbSBarry Smith #endif
13182aabb6bbSBarry Smith 
13192bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1320f1144a30SSatish Balay 
1321be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1322be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
13232bad1931SBarry Smith 
1324be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1325be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1326a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1327ca161407SBarry Smith 
1328be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
13290752156aSBarry Smith 
1330be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1331be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
133271306c60Sjroman 
1333954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
1334be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
133571306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1336be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1337be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
133871306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1339be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1340be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1341be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
134271306c60Sjroman 
134371306c60Sjroman #define MP_PARTY_OPT "opt"
134471306c60Sjroman #define MP_PARTY_LIN "lin"
134571306c60Sjroman #define MP_PARTY_SCA "sca"
134671306c60Sjroman #define MP_PARTY_RAN "ran"
134771306c60Sjroman #define MP_PARTY_GBF "gbf"
134871306c60Sjroman #define MP_PARTY_GCF "gcf"
134971306c60Sjroman #define MP_PARTY_BUB "bub"
135071306c60Sjroman #define MP_PARTY_DEF "def"
1351be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
135271306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
135371306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
135471306c60Sjroman #define MP_PARTY_NONE "no"
1355be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1356be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1357be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1358be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
135971306c60Sjroman 
136071306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1361be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1362be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1363be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1364be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1365be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
136671306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1367be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1368be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
137071306c60Sjroman 
13710752156aSBarry Smith /*
13720a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1373d4fbbf0eSBarry Smith */
13741c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13751c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13761c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13771c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13781c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13797c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13807c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13811c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13821c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13837c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13847c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13851c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13861c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
13871c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
13881c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13891c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13901c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13911c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13921c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13931c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13941c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13951c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
13961c1c02c0SLois Curfman McInnes                MATOP_COMPRESS=22,
13971c1c02c0SLois Curfman McInnes                MATOP_SET_OPTION=23,
13981c1c02c0SLois Curfman McInnes                MATOP_ZERO_ENTRIES=24,
13991c1c02c0SLois Curfman McInnes                MATOP_ZERO_ROWS=25,
14001c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_SYMBOLIC=26,
14011c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_NUMERIC=27,
14021c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_SYMBOLIC=28,
14031c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_NUMERIC=29,
1404d643ce63SMatthew Knepley                MATOP_SETUP_PREALLOCATION=30,
1405d643ce63SMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=31,
1406d643ce63SMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=32,
1407d643ce63SMatthew Knepley                MATOP_GET_ARRAY=33,
1408d643ce63SMatthew Knepley                MATOP_RESTORE_ARRAY=34,
1409e26f1c3cSBarry Smith                MATOP_DUPLICATE=35,
1410d643ce63SMatthew Knepley                MATOP_FORWARD_SOLVE=36,
1411d643ce63SMatthew Knepley                MATOP_BACKWARD_SOLVE=37,
1412d643ce63SMatthew Knepley                MATOP_ILUFACTOR=38,
1413d643ce63SMatthew Knepley                MATOP_ICCFACTOR=39,
1414d643ce63SMatthew Knepley                MATOP_AXPY=40,
1415d643ce63SMatthew Knepley                MATOP_GET_SUBMATRICES=41,
1416d643ce63SMatthew Knepley                MATOP_INCREASE_OVERLAP=42,
1417d643ce63SMatthew Knepley                MATOP_GET_VALUES=43,
1418d643ce63SMatthew Knepley                MATOP_COPY=44,
1419ff1cced0SMatthew Knepley                MATOP_GET_ROW_MAX=45,
1420d643ce63SMatthew Knepley                MATOP_SCALE=46,
1421d643ce63SMatthew Knepley                MATOP_SHIFT=47,
1422d643ce63SMatthew Knepley                MATOP_DIAGONAL_SHIFT=48,
1423d643ce63SMatthew Knepley                MATOP_ILUDT_FACTOR=49,
1424ff1cced0SMatthew Knepley                MATOP_SET_BLOCK_SIZE=50,
1425d643ce63SMatthew Knepley                MATOP_GET_ROW_IJ=51,
1426d643ce63SMatthew Knepley                MATOP_RESTORE_ROW_IJ=52,
1427d643ce63SMatthew Knepley                MATOP_GET_COLUMN_IJ=53,
1428d643ce63SMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=54,
1429d643ce63SMatthew Knepley                MATOP_FDCOLORING_CREATE=55,
1430d643ce63SMatthew Knepley                MATOP_COLORING_PATCH=56,
1431d643ce63SMatthew Knepley                MATOP_SET_UNFACTORED=57,
1432d643ce63SMatthew Knepley                MATOP_PERMUTE=58,
1433d643ce63SMatthew Knepley                MATOP_SET_VALUES_BLOCKED=59,
1434d643ce63SMatthew Knepley                MATOP_GET_SUBMATRIX=60,
1435d643ce63SMatthew Knepley                MATOP_DESTROY=61,
1436d643ce63SMatthew Knepley                MATOP_VIEW=62,
1437ff1cced0SMatthew Knepley                MATOP_CONVERT_FROM=63,
1438d643ce63SMatthew Knepley                MATOP_USE_SCALED_FORM=64,
1439d643ce63SMatthew Knepley                MATOP_SCALE_SYSTEM=65,
1440d643ce63SMatthew Knepley                MATOP_UNSCALE_SYSTEM=66,
1441ba16a8cbSKris Buschelman                MATOP_SET_LOCAL_TO_GLOBAL_MAP=67,
1442d643ce63SMatthew Knepley                MATOP_SET_VALUES_LOCAL=68,
1443d643ce63SMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=69,
1444ff1cced0SMatthew Knepley                MATOP_GET_ROW_MAX_ABS=70,
1445d643ce63SMatthew Knepley                MATOP_CONVERT=71,
1446d643ce63SMatthew Knepley                MATOP_SET_COLORING=72,
1447d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1448d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1449d643ce63SMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1450d643ce63SMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1451ba16a8cbSKris Buschelman                MATOP_MULT_CON=77,
1452ba16a8cbSKris Buschelman                MATOP_MULT_TRANSPOSE_CON=78,
1453ba16a8cbSKris Buschelman                MATOP_ILU_FACTOR_SYMBOLIC_CON=79,
1454d643ce63SMatthew Knepley                MATOP_PERMUTE_SPARSIFY=80,
1455d643ce63SMatthew Knepley                MATOP_MULT_MULTIPLE=81,
145641acf15aSKris Buschelman                MATOP_SOLVE_MULTIPLE=82,
145741acf15aSKris Buschelman                MATOP_GET_INERTIA=83,
145841acf15aSKris Buschelman                MATOP_LOAD=84,
145941acf15aSKris Buschelman                MATOP_IS_SYMMETRIC=85,
146041acf15aSKris Buschelman                MATOP_IS_HERMITIAN=86,
146141acf15aSKris Buschelman                MATOP_IS_STRUCTURALLY_SYMMETRIC=87,
146241acf15aSKris Buschelman                MATOP_PB_RELAX=88,
146341acf15aSKris Buschelman                MATOP_GET_VECS=89,
146441acf15aSKris Buschelman                MATOP_MAT_MULT=90,
146541acf15aSKris Buschelman                MATOP_MAT_MULT_SYMBOLIC=91,
146641acf15aSKris Buschelman                MATOP_MAT_MULT_NUMERIC=92,
146741acf15aSKris Buschelman                MATOP_PTAP=93,
146841acf15aSKris Buschelman                MATOP_PTAP_SYMBOLIC=94,
1469bc011b1eSHong Zhang                MATOP_PTAP_NUMERIC=95,
1470bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE=96,
1471bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97,
14727ba1a0bfSKris Buschelman                MATOP_MAT_MULTTRANSPOSE_NUMERIC=98,
14737ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_SEQAIJ=99,
14747ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_SEQAIJ=100,
14757ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_MPIAIJ=101,
147687d4246cSBarry Smith                MATOP_PTAP_NUMERIC_MPIAIJ=102,
1477ff1cced0SMatthew Knepley                MATOP_CONJUGATE=103,
1478ff1cced0SMatthew Knepley                MATOP_SET_SIZES=104,
1479f5edf698SHong Zhang                MATOP_SET_VALUES_ROW = 105,
1480ff1cced0SMatthew Knepley                MATOP_REAL_PART=106,
1481ff1cced0SMatthew Knepley                MATOP_IMAG_PART=107,
1482d5c63f83SSatish Balay                MATOP_GET_ROW_UTRIANGULAR=108,
14832bebee5dSHong Zhang                MATOP_RESTORE_ROW_UTRIANGULAR=109,
148469db28dcSHong Zhang                MATOP_MATSOLVE=110,
1485829201f2SHong Zhang                MATOP_GET_REDUNDANTMATRIX=111,
1486ff1cced0SMatthew Knepley                MATOP_GET_ROW_MIN=112,
1487ff1cced0SMatthew Knepley                MATOP_GET_COLUMN_VEC=113,
1488ff1cced0SMatthew Knepley                MATOP_MISSING_DIAGONAL=114,
1489ff1cced0SMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=115,
1490ff1cced0SMatthew Knepley                MATOP_DESTROY_SOLVER=116
1491fae171e0SBarry Smith              } MatOperation;
1492be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1493be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1494be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1495be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1496112a2221SBarry Smith 
149790ace30eSBarry Smith /*
149890ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
149990ace30eSBarry Smith  stored in a universal format. By changing the format with
1500fb9695e5SSatish Balay  PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will
150190ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
150290ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
150390ace30eSBarry Smith  read into matrices of the same time.
150490ace30eSBarry Smith */
150590ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
150690ace30eSBarry Smith 
1507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1508be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
15091f4e1ec7SBarry Smith 
1510d9274352SBarry Smith /*S
1511d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1512d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1513d9274352SBarry Smith 
1514f7a9e4ceSBarry Smith    Level: advanced
1515d9274352SBarry Smith 
1516d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1517d9274352SBarry Smith 
15186e1639daSBarry Smith   Users manual sections:
15196e1639daSBarry Smith .   sec_singular
15206e1639daSBarry Smith 
1521d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1522d9274352SBarry Smith S*/
152374637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1524d9274352SBarry Smith 
1525be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1526281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*);
1527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1529be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
1530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat);
153174637425SBarry Smith 
1532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
153546129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
15363f1d51d7SBarry Smith 
1537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1538be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1540c4f061fbSSatish Balay 
1541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1542b0a32e0cSBarry Smith 
1543be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
154404f1ad80SBarry Smith 
1545fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
15463ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec);
15473ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1548e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1549e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1550e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace);
1551e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1552e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat);
1553e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal);
1554e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt);
1555e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *);
1556e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat);
1557e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1558e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1559e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal);
1560e884886fSBarry Smith 
1561e884886fSBarry Smith 
1562e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1563e884886fSBarry Smith 
1564e884886fSBarry Smith /*E
1565e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1566e884886fSBarry Smith 
1567e884886fSBarry Smith    Level: beginner
1568e884886fSBarry Smith 
1569e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1570e884886fSBarry Smith E*/
1571a313700dSBarry Smith #define MatMFFDType char*
1572e884886fSBarry Smith #define MATMFFD_DS  "ds"
1573e884886fSBarry Smith #define MATMFFD_WP  "wp"
1574e884886fSBarry Smith 
1575a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,const MatMFFDType);
1576e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1577e884886fSBarry Smith 
1578e884886fSBarry Smith /*MC
1579e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1580e884886fSBarry Smith 
1581e884886fSBarry Smith    Synopsis:
1582e884886fSBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1583e884886fSBarry Smith 
1584e884886fSBarry Smith    Not Collective
1585e884886fSBarry Smith 
1586e884886fSBarry Smith    Input Parameters:
1587e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1588e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1589e884886fSBarry Smith .  name_create - name of routine to create method context
1590e884886fSBarry Smith -  routine_create - routine to create method context
1591e884886fSBarry Smith 
1592e884886fSBarry Smith    Level: developer
1593e884886fSBarry Smith 
1594e884886fSBarry Smith    Notes:
1595e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1596e884886fSBarry Smith 
1597e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1598e884886fSBarry Smith    is ignored.
1599e884886fSBarry Smith 
1600e884886fSBarry Smith    Sample usage:
1601e884886fSBarry Smith .vb
1602e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1603e884886fSBarry Smith                "MyHCreate",MyHCreate);
1604e884886fSBarry Smith .ve
1605e884886fSBarry Smith 
1606e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1607e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1608e884886fSBarry Smith    or at runtime via the option
1609e884886fSBarry Smith $     -snes_mf_type my_h
1610e884886fSBarry Smith 
1611e884886fSBarry Smith .keywords: MatMFFD, register
1612e884886fSBarry Smith 
1613e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1614e884886fSBarry Smith M*/
1615e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1616e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1617e884886fSBarry Smith #else
1618e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1619e884886fSBarry Smith #endif
1620e884886fSBarry Smith 
1621e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]);
1622e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void);
1623e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal);
1624e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth);
1625e884886fSBarry Smith 
1626e884886fSBarry Smith 
1627be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1628be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
16297dbadf16SMatthew Knepley 
1630e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
16312eac72dbSBarry Smith #endif
1632