xref: /petsc/include/petscmat.h (revision 5494a0643c6f9624505f1aa1d78ebbf006b2a7a9)
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 
27d9274352SBarry Smith .seealso: MatSetType(), Mat
28d91e6319SBarry Smith E*/
29e5a9bf91SBarry Smith #define MatType const 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"
40273d9f13SBarry Smith #define MATSEQBDIAG        "seqbdiag"
41273d9f13SBarry Smith #define MATMPIBDIAG        "mpibdiag"
42209238afSKris Buschelman #define MATBDIAG           "bdiag"
43209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
44273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
45209238afSKris Buschelman #define MATDENSE           "dense"
46273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
47273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
48209238afSKris Buschelman #define MATBAIJ            "baij"
49273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
50273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
51273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
52209238afSKris Buschelman #define MATSBAIJ           "sbaij"
53cebc7f6cSBarry Smith #define MATDAAD            "daad"
54cebc7f6cSBarry Smith #define MATMFFD            "mffd"
55c8a8475eSBarry Smith #define MATNORMAL          "normal"
56ab92ecdeSBarry Smith #define MATLRC             "lrc"
57b3a1e11cSKris Buschelman #define MATSEQAIJSPOOLES   "seqaijspooles"
58d10c748bSKris Buschelman #define MATMPIAIJSPOOLES   "mpiaijspooles"
599abb65ffSKris Buschelman #define MATSEQSBAIJSPOOLES "seqsbaijspooles"
6022191285SKris Buschelman #define MATMPISBAIJSPOOLES "mpisbaijspooles"
61bb4d4166SHong Zhang #define MATAIJSPOOLES      "aijspooles"
62bb4d4166SHong Zhang #define MATSBAIJSPOOLES    "sbaijspooles"
6314b396bbSKris Buschelman #define MATSUPERLU         "superlu"
641d96aa28SKris Buschelman #define MATSUPERLU_DIST    "superlu_dist"
651677d0b8SKris Buschelman #define MATUMFPACK         "umfpack"
66e8aa55a4SKris Buschelman #define MATESSL            "essl"
674eb8e494SKris Buschelman #define MATLUSOL           "lusol"
68397b6df1SKris Buschelman #define MATAIJMUMPS        "aijmumps"
69397b6df1SKris Buschelman #define MATSBAIJMUMPS      "sbaijmumps"
708da957c5SKris Buschelman #define MATDSCPACK         "dscpack"
717065e2aaSKris Buschelman #define MATMATLAB          "matlab"
724b8d99adSRichard Tran Mills #define MATSEQCSRPERM      "seqcsrperm"
7318def556SRichard Tran Mills #define MATMPICSRPERM      "mpicsrperm"
74814655b8SBarry Smith #define MATCSRPERM         "csrperm"
7581824310SBarry Smith #define MATSEQCRL          "seqcrl"
76c02b24c5SBarry Smith #define MATMPICRL          "mpicrl"
77c02b24c5SBarry Smith #define MATCRL             "crl"
78c0aa2d19SHong Zhang #define MATPLAPACK         "plapack"
792a6744ebSBarry Smith #define MATSCATTER         "scatter"
80421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
81793850ffSBarry Smith #define MATCOMPOSITE       "composite"
825a7f1df3SHong Zhang #define MATSEQFFTW         "seqfftw"
83d91e6319SBarry Smith 
84c06d978dSMatthew Knepley /* Logging support */
85552e946dSBarry Smith #define    MAT_FILE_COOKIE 1211216    /* used to indicate matrices in binary files */
86be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE;
87be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
88be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE;
89be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE;
90e884886fSBarry Smith extern PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE;
91c06d978dSMatthew Knepley 
92ceb03754SKris Buschelman /*E
93ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
94ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
95ceb03754SKris Buschelman 
96ceb03754SKris Buschelman     Level: beginner
97ceb03754SKris Buschelman 
98ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
99ceb03754SKris Buschelman 
1000c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
101ceb03754SKris Buschelman E*/
102ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse;
103ceb03754SKris Buschelman 
104*5494a064SHong Zhang /*E
105*5494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
106*5494a064SHong Zhang      include the matrix values. Currently it is only used by MatGetSequentialNonzerostructure().
107*5494a064SHong Zhang 
108*5494a064SHong Zhang     Level: beginner
109*5494a064SHong Zhang 
110*5494a064SHong Zhang .seealso: MatGetSequentialNonzerostructure()
111*5494a064SHong Zhang E*/
112*5494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
113*5494a064SHong Zhang 
114e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]);
115c06d978dSMatthew Knepley 
116f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
117a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
118a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
119f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
120f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,MatType);
121be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
122be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
123be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
124be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
12530de9b25SBarry Smith 
126f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
127f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
128f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
129f69a0ea3SMatthew Knepley 
13030de9b25SBarry Smith /*MC
13130de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
13230de9b25SBarry Smith 
13330de9b25SBarry Smith    Synopsis:
134c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
13530de9b25SBarry Smith 
13630de9b25SBarry Smith    Not Collective
13730de9b25SBarry Smith 
13830de9b25SBarry Smith    Input Parameters:
13930de9b25SBarry Smith +  name - name of a new user-defined matrix type
14030de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
14130de9b25SBarry Smith .  name_create - name of routine to create method context
14230de9b25SBarry Smith -  routine_create - routine to create method context
14330de9b25SBarry Smith 
14430de9b25SBarry Smith    Notes:
14530de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
14630de9b25SBarry Smith 
14730de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
14830de9b25SBarry Smith    is ignored.
14930de9b25SBarry Smith 
15030de9b25SBarry Smith    Sample usage:
15130de9b25SBarry Smith .vb
15230de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
15330de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
15430de9b25SBarry Smith .ve
15530de9b25SBarry Smith 
15630de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
15730de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
15830de9b25SBarry Smith    or at runtime via the option
15930de9b25SBarry Smith $     -mat_type my_mat
16030de9b25SBarry Smith 
16130de9b25SBarry Smith    Level: advanced
16230de9b25SBarry Smith 
163ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
16430de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
16530de9b25SBarry Smith 
16630de9b25SBarry Smith .keywords: Mat, register
16730de9b25SBarry Smith 
16830de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
16930de9b25SBarry Smith 
17030de9b25SBarry Smith M*/
171273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
172273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
173273d9f13SBarry Smith #else
174273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
17530de9b25SBarry Smith #endif
17630de9b25SBarry Smith 
1779c666560SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
178273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
179b0a32e0cSBarry Smith extern PetscFList MatList;
18028988994SBarry Smith 
181be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
182be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
183be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
184ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
185ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
186ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
187ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
188ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
189ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
190ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
191be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
192ba966639SSatish 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)
193ba966639SSatish 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)
194ba966639SSatish 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)
195ba966639SSatish 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)
196ba966639SSatish 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))
197ba966639SSatish 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))
198ba966639SSatish 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))
199ba966639SSatish 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)
200ba966639SSatish 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)
201ba966639SSatish 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)
202ba966639SSatish 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)
203ba966639SSatish 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))
204ba966639SSatish 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))
205ba966639SSatish 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))
20663ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2078d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2088d7a6e47SBarry Smith 
209be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
210be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*);
211be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*);
212be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
213ba966639SSatish 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)
214ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
215ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
216ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
217ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
218ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
219ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
220be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
221ba966639SSatish 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)
222ba966639SSatish 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)
223ba966639SSatish 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)
224ba966639SSatish 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)
225ba966639SSatish 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))
226ba966639SSatish 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))
227ba966639SSatish 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))
228ba966639SSatish 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)
229ba966639SSatish 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)
230ba966639SSatish 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)
231ba966639SSatish 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)
232ba966639SSatish 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))
233ba966639SSatish 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))
234ba966639SSatish 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))
235be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
236be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
237ba966639SSatish 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)
238ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
239ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
240ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
241ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
242ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
243ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
244be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
245ba966639SSatish 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)
246ba966639SSatish 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)
247ba966639SSatish 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)
248ba966639SSatish 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)
249ba966639SSatish 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))
250ba966639SSatish 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))
251ba966639SSatish 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))
252ba966639SSatish 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)
253ba966639SSatish 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)
254ba966639SSatish 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)
255ba966639SSatish 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)
256ba966639SSatish 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))
257ba966639SSatish 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))
258ba966639SSatish 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))
259be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
260ba966639SSatish 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)
261ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
262be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
263be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
264ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
265ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*);
266284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
2671472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
2681472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
2692a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*);
2702a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter);
2712a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*);
2728cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
273793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat);
274b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat);
275793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2765a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*);
2771472f72bSBarry Smith 
278f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
279be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
28021c89e3eSBarry Smith 
2810c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
28299cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat);
28399cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat);
28499cafbc1SBarry Smith 
2858ed539a5SBarry Smith /* ------------------------------------------------------------*/
286be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
287be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
28887d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
289f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
29084cb2905SBarry Smith 
2912ef4de8bSBarry Smith /*S
2922ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
2932ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
2942ef4de8bSBarry Smith 
2952ef4de8bSBarry Smith    Level: beginner
2962ef4de8bSBarry Smith 
2972ef4de8bSBarry Smith   Concepts: matrix; linear operator
2982ef4de8bSBarry Smith 
299d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3002ef4de8bSBarry Smith S*/
301435da068SBarry Smith typedef struct {
302c1ac3661SBarry Smith   PetscInt k,j,i,c;
303435da068SBarry Smith } MatStencil;
3042ef4de8bSBarry Smith 
305be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
306be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
307be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
308435da068SBarry Smith 
309be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
310be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
311be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
3123a7fca6bSBarry Smith 
313d91e6319SBarry Smith /*E
314d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
315d91e6319SBarry Smith      to continue to add values to it
316d91e6319SBarry Smith 
317d91e6319SBarry Smith     Level: beginner
318d91e6319SBarry Smith 
319d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
320d91e6319SBarry Smith E*/
3216d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
322be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
323be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
324be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
3254f9c727eSBarry Smith 
3261d73ed98SBarry Smith 
32730de9b25SBarry Smith 
328d91e6319SBarry Smith /*E
329d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
330d91e6319SBarry Smith 
331d91e6319SBarry Smith     Level: beginner
332d91e6319SBarry Smith 
3330a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
334d91e6319SBarry Smith 
335d91e6319SBarry Smith .seealso: MatSetOption()
336d91e6319SBarry Smith E*/
337512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
3384e0d8c25SBarry Smith               MAT_SYMMETRIC,
3394e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
3408047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
3414e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
3424e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
3434e0d8c25SBarry Smith               MAT_KEEP_ZEROED_ROWS,MAT_IGNORE_ZERO_ENTRIES,
3444e0d8c25SBarry Smith               MAT_USE_INODES,
3454e0d8c25SBarry Smith               MAT_HERMITIAN,
3464e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
3474e0d8c25SBarry Smith               MAT_USE_COMPRESSEDROW,
3484e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
3494e0d8c25SBarry Smith               MAT_GETROW_UPPERTRIANGULAR} MatOption;
350290bbb0aSBarry Smith extern const char *MatOptions[];
3514e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth);
352be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,MatType*);
353ba966639SSatish Balay PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),MatType,t)
35484cb2905SBarry Smith 
355be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
356be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
357be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
358f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat);
359f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat);
360be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
361be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
362be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
363be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
364ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
365be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
366be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
367ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
368be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
3697b80b807SBarry Smith 
3701620fd73SBarry Smith 
371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
373ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
374be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
375be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
376ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
377e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t)
3781cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*);
379be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
380ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
382be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
3832bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat);
384f5edf698SHong Zhang 
385d91e6319SBarry Smith /*E
386d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
387d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
388d91e6319SBarry Smith 
389d91e6319SBarry Smith     Level: beginner
390d91e6319SBarry Smith 
391d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
392d91e6319SBarry Smith 
393d91e6319SBarry Smith .seealso: MatDuplicate()
394d91e6319SBarry Smith E*/
3952e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
3962e8a6d31SBarry Smith 
397f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,MatType,MatReuse,Mat*);
398ba966639SSatish Balay PetscPolymorphicFunction(MatConvert,(Mat A,MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
400ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
401ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
40294a9d846SBarry Smith 
403d91e6319SBarry Smith /*E
404d91e6319SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
405d91e6319SBarry Smith 
406d91e6319SBarry Smith     Level: beginner
407d91e6319SBarry Smith 
408d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
409d91e6319SBarry Smith 
41094b7f48cSBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
411d91e6319SBarry Smith E*/
412c537a176SHong Zhang typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
413cb5b572fSBarry Smith 
414be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
415be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
416be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
417ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
418e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t)
419be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
420ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
4211cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*);
4221cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t)
423be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
424be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
42564c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *);
426f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,MatType,Mat*);
427ba966639SSatish Balay PetscPolymorphicFunction(MatLoad,(PetscViewer v,MatType t),(v,t,&a),Mat,a)
4287b80b807SBarry Smith 
4298f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4308f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
4318f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4328f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
433d4fbbf0eSBarry Smith 
434d91e6319SBarry Smith /*S
435d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
436d91e6319SBarry Smith 
437d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
438d91e6319SBarry Smith 
439d91e6319SBarry Smith    Level: intermediate
440d91e6319SBarry Smith 
441d91e6319SBarry Smith   Concepts: matrix^nonzero information
442d91e6319SBarry Smith 
443d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
444d91e6319SBarry Smith S*/
4454e220ebcSLois Curfman McInnes typedef struct {
446b0a32e0cSBarry Smith   PetscLogDouble rows_global,columns_global;         /* number of global rows and columns */
447b0a32e0cSBarry Smith   PetscLogDouble rows_local,columns_local;           /* number of local rows and columns */
448b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
449b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
450b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
451b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
452b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
453b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
454b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4554e220ebcSLois Curfman McInnes } MatInfo;
4564e220ebcSLois Curfman McInnes 
457d9274352SBarry Smith /*E
458d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
459d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
460d9274352SBarry Smith 
461d9274352SBarry Smith     Level: beginner
462d9274352SBarry Smith 
463d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
464d9274352SBarry Smith 
465d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
466d9274352SBarry Smith E*/
4677b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
468be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
469be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
470be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
471985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]);
472985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]);
473985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
47486697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec);
475be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,Mat*);
476ba966639SSatish Balay PetscPolymorphicFunction(MatTranspose,(Mat A),(A,&t),Mat,t)
477be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
478ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
479be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
480be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
481be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
482be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
483ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
484be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
485be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
486be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
487be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
4887b80b807SBarry Smith 
489be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
490ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
491be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
492f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
493f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
494f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
495f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
4967b80b807SBarry Smith 
497be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
499be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5005ef9f2a5SBarry Smith 
501be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
502be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
503be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5048ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**);
505f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
5067b80b807SBarry Smith 
507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
508be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
509be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *);
510abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *);
511*5494a064SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSequentialNonzeroStructure(Mat,Mat *[]);
512*5494a064SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySequentialNonzeroStructure(Mat *[]);
513*5494a064SHong Zhang 
514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
515be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
517be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
519be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
520be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
521be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
522be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscScalar**,Mat*);
52343eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
524cd116e26SSatish Balay #include "petscctable.h"
52543eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
52643eb5e2fSMatthew Knepley #else
52743eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
52843eb5e2fSMatthew Knepley #endif
5298efafbd8SBarry Smith 
530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5317b80b807SBarry Smith 
532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
53522440eb1SKris Buschelman 
536be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
538be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
53922440eb1SKris Buschelman 
540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
542be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
543bc011b1eSHong Zhang 
544f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
54504aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure);
546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat);
5471c741599SHong Zhang 
548f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
549f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5507b80b807SBarry Smith 
551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
553f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
554f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
557052efed2SBarry Smith 
558be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
559be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
56090f02eecSBarry Smith 
561be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
5620c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
563ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
565be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
56669db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
567bd481603SBarry Smith 
568bd481603SBarry Smith /*MC
5691620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5701620fd73SBarry Smith 
5711620fd73SBarry Smith    Not collective
5721620fd73SBarry Smith 
5731620fd73SBarry Smith    Input Parameters:
5741620fd73SBarry Smith +  m - the matrix
5751620fd73SBarry Smith .  row - the row location of the entry
5761620fd73SBarry Smith .  col - the column location of the entry
5771620fd73SBarry Smith .  value - the value to insert
5781620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5791620fd73SBarry Smith 
5801620fd73SBarry Smith    Notes:
5811620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
5821620fd73SBarry Smith    values simultaneously if possible.
5831620fd73SBarry Smith 
5841620fd73SBarry Smith    Level: beginner
5851620fd73SBarry Smith 
5861620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
5871620fd73SBarry Smith M*/
5881620fd73SBarry 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);}
5891620fd73SBarry Smith 
5901620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
5911620fd73SBarry Smith 
5921620fd73SBarry 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);}
5931620fd73SBarry Smith 
5941620fd73SBarry Smith /*MC
5950d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
596bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
597bd481603SBarry Smith 
598bd481603SBarry Smith    Synopsis:
599c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
600bd481603SBarry Smith 
601bd481603SBarry Smith    Collective on MPI_Comm
602bd481603SBarry Smith 
603bd481603SBarry Smith    Input Parameters:
604bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
605859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
606859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
607bd481603SBarry Smith 
608bd481603SBarry Smith    Output Parameters:
609bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
610bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
611bd481603SBarry Smith 
612bd481603SBarry Smith 
613bd481603SBarry Smith    Level: intermediate
614bd481603SBarry Smith 
615bd481603SBarry Smith    Notes:
616bd481603SBarry Smith    See the chapter in the users manual on performance for more details
617bd481603SBarry Smith 
6181d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
619bd481603SBarry Smith 
620bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
621bd481603SBarry Smith 
6221620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6231620fd73SBarry Smith 
624bd481603SBarry Smith   Concepts: preallocation^Matrix
625bd481603SBarry Smith 
626bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
627bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
628bd481603SBarry Smith M*/
629c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6307c922b88SBarry Smith { \
631c1ac3661SBarry Smith   PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
632c1ac3661SBarry Smith   _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\
633c1ac3661SBarry Smith   _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\
634c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
635c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp;
6367c922b88SBarry Smith 
637bd481603SBarry Smith /*MC
6380d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
639bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
640bd481603SBarry Smith 
641bd481603SBarry Smith    Synopsis:
642c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
643bd481603SBarry Smith 
644bd481603SBarry Smith    Collective on MPI_Comm
645bd481603SBarry Smith 
646bd481603SBarry Smith    Input Parameters:
647bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
648859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
649859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
650bd481603SBarry Smith 
651bd481603SBarry Smith    Output Parameters:
652bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
653bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
654bd481603SBarry Smith 
655bd481603SBarry Smith 
656bd481603SBarry Smith    Level: intermediate
657bd481603SBarry Smith 
658bd481603SBarry Smith    Notes:
659bd481603SBarry Smith    See the chapter in the users manual on performance for more details
660bd481603SBarry Smith 
6611d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
662bd481603SBarry Smith 
6631620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6641620fd73SBarry Smith 
665bd481603SBarry Smith   Concepts: preallocation^Matrix
666bd481603SBarry Smith 
667bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
668bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
669bd481603SBarry Smith M*/
670222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
671222b16d4SBarry Smith { \
672c1ac3661SBarry Smith   PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__end; \
673c1ac3661SBarry Smith   _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\
674c1ac3661SBarry Smith   _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\
675c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
676c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp;
677222b16d4SBarry Smith 
678bd481603SBarry Smith /*MC
679bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
680bd481603SBarry Smith        inserted using a local number of the rows and columns
681bd481603SBarry Smith 
682bd481603SBarry Smith    Synopsis:
683c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
684bd481603SBarry Smith 
685bd481603SBarry Smith    Not Collective
686bd481603SBarry Smith 
687bd481603SBarry Smith    Input Parameters:
688bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
689bd481603SBarry Smith .  nrows - the number of rows indicated
6901d73ed98SBarry Smith .  rows - the indices of the rows
691bd481603SBarry Smith .  ncols - the number of columns in the matrix
692bd481603SBarry Smith .  cols - the columns indicated
693bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
694bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
695bd481603SBarry Smith 
696bd481603SBarry Smith 
697bd481603SBarry Smith    Level: intermediate
698bd481603SBarry Smith 
699bd481603SBarry Smith    Notes:
700bd481603SBarry Smith    See the chapter in the users manual on performance for more details
701bd481603SBarry Smith 
7021d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
703bd481603SBarry Smith 
704bd481603SBarry Smith   Concepts: preallocation^Matrix
705bd481603SBarry Smith 
706bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
707bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
708bd481603SBarry Smith M*/
709c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
710c4f061fbSSatish Balay {\
711c1ac3661SBarry Smith   PetscInt __l;\
712ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
713ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
714c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
715ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
716c4f061fbSSatish Balay   }\
717c4f061fbSSatish Balay }
718c4f061fbSSatish Balay 
719bd481603SBarry Smith /*MC
720bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
721bd481603SBarry Smith        inserted using a local number of the rows and columns
722bd481603SBarry Smith 
723bd481603SBarry Smith    Synopsis:
724c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
725bd481603SBarry Smith 
726bd481603SBarry Smith    Not Collective
727bd481603SBarry Smith 
728bd481603SBarry Smith    Input Parameters:
729bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
730bd481603SBarry Smith .  nrows - the number of rows indicated
7311d73ed98SBarry Smith .  rows - the indices of the rows
732bd481603SBarry Smith .  ncols - the number of columns in the matrix
733bd481603SBarry Smith .  cols - the columns indicated
734bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
735bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
736bd481603SBarry Smith 
737bd481603SBarry Smith 
738bd481603SBarry Smith    Level: intermediate
739bd481603SBarry Smith 
740bd481603SBarry Smith    Notes:
741bd481603SBarry Smith    See the chapter in the users manual on performance for more details
742bd481603SBarry Smith 
743bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
744bd481603SBarry Smith 
745bd481603SBarry Smith   Concepts: preallocation^Matrix
746bd481603SBarry Smith 
747bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
748bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
749bd481603SBarry Smith M*/
750d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
751d3d32019SBarry Smith {\
752c1ac3661SBarry Smith   PetscInt __l;\
753d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
754d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
755d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
756d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
757d3d32019SBarry Smith   }\
758d3d32019SBarry Smith }
759d3d32019SBarry Smith 
760bd481603SBarry Smith /*MC
761bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
762bd481603SBarry Smith        inserted using a local number of the rows and columns
763bd481603SBarry Smith 
764bd481603SBarry Smith    Synopsis:
765c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
766bd481603SBarry Smith 
767bd481603SBarry Smith    Not Collective
768bd481603SBarry Smith 
769bd481603SBarry Smith    Input Parameters:
77064075487SBarry Smith +  row - the row
771bd481603SBarry Smith .  ncols - the number of columns in the matrix
772a50f8bf6SHong Zhang -  cols - the columns indicated
773a50f8bf6SHong Zhang 
774a50f8bf6SHong Zhang    Output Parameters:
775a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
776bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
777bd481603SBarry Smith 
778bd481603SBarry Smith 
779bd481603SBarry Smith    Level: intermediate
780bd481603SBarry Smith 
781bd481603SBarry Smith    Notes:
782bd481603SBarry Smith    See the chapter in the users manual on performance for more details
783bd481603SBarry Smith 
784bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
785bd481603SBarry Smith 
7861620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7871620fd73SBarry Smith 
788bd481603SBarry Smith   Concepts: preallocation^Matrix
789bd481603SBarry Smith 
790bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
791bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
792bd481603SBarry Smith M*/
793c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
794c1ac3661SBarry Smith { PetscInt __i; \
7958ee2e534SBarry Smith   if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
7968ee2e534SBarry Smith   if (row >= __rstart+__tmp) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D greater than last local row %D",row,__rstart+__tmp-1);\
7977c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
79864075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
7997cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8007c922b88SBarry Smith   }\
8017c922b88SBarry Smith }
8027c922b88SBarry Smith 
803bd481603SBarry Smith /*MC
804bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
805bd481603SBarry Smith        inserted using a local number of the rows and columns
806bd481603SBarry Smith 
807bd481603SBarry Smith    Synopsis:
808c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
809bd481603SBarry Smith 
810bd481603SBarry Smith    Not Collective
811bd481603SBarry Smith 
812bd481603SBarry Smith    Input Parameters:
813bd481603SBarry Smith +  nrows - the number of rows indicated
8141d73ed98SBarry Smith .  rows - the indices of the rows
815bd481603SBarry Smith .  ncols - the number of columns in the matrix
816bd481603SBarry Smith .  cols - the columns indicated
817bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
818bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
819bd481603SBarry Smith 
820bd481603SBarry Smith 
821bd481603SBarry Smith    Level: intermediate
822bd481603SBarry Smith 
823bd481603SBarry Smith    Notes:
824bd481603SBarry Smith    See the chapter in the users manual on performance for more details
825bd481603SBarry Smith 
826bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
827bd481603SBarry Smith 
8281620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8291620fd73SBarry Smith 
830bd481603SBarry Smith   Concepts: preallocation^Matrix
831bd481603SBarry Smith 
832bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
833bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
834bd481603SBarry Smith M*/
835d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
836c1ac3661SBarry Smith { PetscInt __i; \
837d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
838d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
839d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
840d3d32019SBarry Smith   }\
841d3d32019SBarry Smith }
842d3d32019SBarry Smith 
843bd481603SBarry Smith /*MC
8440d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
845bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
846bd481603SBarry Smith 
847bd481603SBarry Smith    Synopsis:
848c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
849bd481603SBarry Smith 
850bd481603SBarry Smith    Collective on MPI_Comm
851bd481603SBarry Smith 
852bd481603SBarry Smith    Input Parameters:
853bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
854bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
855bd481603SBarry Smith 
856bd481603SBarry Smith 
857bd481603SBarry Smith    Level: intermediate
858bd481603SBarry Smith 
859bd481603SBarry Smith    Notes:
860bd481603SBarry Smith    See the chapter in the users manual on performance for more details
861bd481603SBarry Smith 
862bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
863bd481603SBarry Smith 
8641620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
8651620fd73SBarry Smith 
866bd481603SBarry Smith   Concepts: preallocation^Matrix
867bd481603SBarry Smith 
868bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
869bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
870bd481603SBarry Smith M*/
871ef66eb69SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree(dnz);CHKERRQ(_4_ierr);}
8727c922b88SBarry Smith 
873bd481603SBarry Smith 
874bd481603SBarry Smith 
8757b80b807SBarry Smith /* Routines unique to particular data structures */
876be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
877ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
878be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBDiagGetData(Mat,PetscInt*,PetscInt*,PetscInt*[],PetscInt*[],PetscScalar***);
879698d4c6aSKris Buschelman 
880be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
881be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
882698d4c6aSKris Buschelman 
883be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
884be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
885be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
886c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
887c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
8887b80b807SBarry Smith 
889a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
890a96a251dSBarry Smith 
891be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
892ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
893be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
894ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
895be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
896ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
897be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDensePreallocation(Mat,PetscScalar[]);
898be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]);
899273d9f13SBarry Smith 
900be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
901ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
902be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
903be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
90453803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
905be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
906be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
907be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]);
908be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
909be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
910284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
911be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
912be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
913be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
914be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
915273d9f13SBarry Smith 
916be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
917ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*);
9181b807ce4Svictorle 
919be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
920be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
9212e8a6d31SBarry Smith 
922be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
9233a7fca6bSBarry Smith 
9247b80b807SBarry Smith /*
9257b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
92694b7f48cSBarry Smith   done through the KSP and PC interfaces.
9277b80b807SBarry Smith */
9287b80b807SBarry Smith 
929d9274352SBarry Smith /*E
930d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
931d9274352SBarry Smith        with an optional dynamic library name, for example
932d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
933d9274352SBarry Smith 
934d9274352SBarry Smith    Level: beginner
935d9274352SBarry Smith 
936e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
937e5a9bf91SBarry Smith 
938d9274352SBarry Smith .seealso: MatGetOrdering()
939d9274352SBarry Smith E*/
94049773a63SBarry Smith #define MatOrderingType char*
941b12f92e5SBarry Smith #define MATORDERING_NATURAL   "natural"
942b12f92e5SBarry Smith #define MATORDERING_ND        "nd"
943b12f92e5SBarry Smith #define MATORDERING_1WD       "1wd"
944b12f92e5SBarry Smith #define MATORDERING_RCM       "rcm"
945b12f92e5SBarry Smith #define MATORDERING_QMD       "qmd"
946b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength"
94762152c8bSBarry Smith #define MATORDERING_DSC_ND    "dsc_nd"
94862152c8bSBarry Smith #define MATORDERING_DSC_MMD   "dsc_mmd"
94962152c8bSBarry Smith #define MATORDERING_DSC_MDF   "dsc_mdf"
950c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
951c06d978dSMatthew Knepley #define MATORDERING_IDENTITY  "identity"
952c06d978dSMatthew Knepley #define MATORDERING_REVERSE   "reverse"
953b12f92e5SBarry Smith 
954be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
955be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
956be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
95730de9b25SBarry Smith 
95830de9b25SBarry Smith /*MC
95930de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
96030de9b25SBarry Smith                                matrix package.
96130de9b25SBarry Smith 
96230de9b25SBarry Smith    Synopsis:
963c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
96430de9b25SBarry Smith 
96530de9b25SBarry Smith    Not Collective
96630de9b25SBarry Smith 
96730de9b25SBarry Smith    Input Parameters:
96830de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
96930de9b25SBarry Smith .  path - location of library where creation routine is
97030de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
97130de9b25SBarry Smith -  function - function pointer that creates the ordering
97230de9b25SBarry Smith 
97330de9b25SBarry Smith    Level: developer
97430de9b25SBarry Smith 
97530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
97630de9b25SBarry Smith    is ignored.
97730de9b25SBarry Smith 
97830de9b25SBarry Smith    Sample usage:
97930de9b25SBarry Smith .vb
98030de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
98130de9b25SBarry Smith                "MyOrder",MyOrder);
98230de9b25SBarry Smith .ve
98330de9b25SBarry Smith 
98430de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
98530de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
98630de9b25SBarry Smith    or at runtime via the option
9872401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
98830de9b25SBarry Smith 
989ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
99030de9b25SBarry Smith 
99130de9b25SBarry Smith .keywords: matrix, ordering, register
99230de9b25SBarry Smith 
99330de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
99430de9b25SBarry Smith M*/
995aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
996f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
997b12f92e5SBarry Smith #else
998f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
999b12f92e5SBarry Smith #endif
100030de9b25SBarry Smith 
1001be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
1002be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
10032bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
1004b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1005d4fbbf0eSBarry Smith 
1006be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1007a2ce50c7SBarry Smith 
1008d91e6319SBarry Smith /*S
10092401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1010b00f7748SHong Zhang 
101161cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
101261cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1013b00f7748SHong Zhang 
101415e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1015b00f7748SHong Zhang 
101661cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
101761cad860SBarry Smith 
1018b00f7748SHong Zhang    Level: developer
1019b00f7748SHong Zhang 
1020d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1021d7d82daaSBarry Smith           MatFactorInfoInitialize()
1022b00f7748SHong Zhang 
1023b00f7748SHong Zhang S*/
1024b00f7748SHong Zhang typedef struct {
10250a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
1026fbf22428SSatish Balay   PetscReal     shiftpd;        /* if true, shift until positive pivots */
10272cea7109SBarry Smith   PetscReal     shift_fraction; /* record shift fraction taken */
102815e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
102915e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1030b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
103115e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
103267e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1033348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1034bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1035bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
103615e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
103715e8a5b3SHong Zhang } MatFactorInfo;
1038ffa6d0a5SLois Curfman McInnes 
1039be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
1040be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,MatFactorInfo*);
1041be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1042be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,MatFactorInfo*,Mat*);
1043be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,MatFactorInfo*);
1044be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,MatFactorInfo*);
1045be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1046be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1047be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1048be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,MatFactorInfo*);
1049be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,MatFactorInfo*,Mat*);
1050be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,MatFactorInfo*,Mat *);
1051be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1052be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1053be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1054be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1055be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1056be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1057be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1058be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
10598ed539a5SBarry Smith 
1060be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1061bb5a7306SBarry Smith 
1062d91e6319SBarry Smith /*E
1063d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1064bb1eb677SSatish Balay 
1065d91e6319SBarry Smith     Level: beginner
1066d91e6319SBarry Smith 
1067d9274352SBarry Smith    May be bitwise ORd together
1068d9274352SBarry Smith 
1069d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1070d91e6319SBarry Smith 
10714e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10724e7234bfSBarry Smith 
1073a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax()
1074d91e6319SBarry Smith E*/
1075ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1076ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1077ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
107884cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1079be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1080be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10818ed539a5SBarry Smith 
1082d4fbbf0eSBarry Smith /*
1083639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1084639f9d9dSBarry Smith */
1085b12f92e5SBarry Smith 
1086d9274352SBarry Smith /*E
1087d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1088d9274352SBarry Smith        with an optional dynamic library name, for example
1089d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1090d9274352SBarry Smith 
1091d9274352SBarry Smith    Level: beginner
1092d9274352SBarry Smith 
1093d9274352SBarry Smith .seealso: MatGetColoring()
1094d9274352SBarry Smith E*/
1095e5a9bf91SBarry Smith #define MatColoringType const char*
1096b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1097b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1098b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1099b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1100b12f92e5SBarry Smith 
11012e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,MatColoringType,ISColoring*);
11022e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
110330de9b25SBarry Smith 
110430de9b25SBarry Smith /*MC
110530de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
110630de9b25SBarry Smith                                matrix package.
110730de9b25SBarry Smith 
110830de9b25SBarry Smith    Synopsis:
1109c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
111030de9b25SBarry Smith 
111130de9b25SBarry Smith    Not Collective
111230de9b25SBarry Smith 
111330de9b25SBarry Smith    Input Parameters:
111430de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
111530de9b25SBarry Smith .  path - location of library where creation routine is
111630de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
111730de9b25SBarry Smith -  function - function pointer that creates the coloring
111830de9b25SBarry Smith 
111930de9b25SBarry Smith    Level: developer
112030de9b25SBarry Smith 
112130de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
112230de9b25SBarry Smith    is ignored.
112330de9b25SBarry Smith 
112430de9b25SBarry Smith    Sample usage:
112530de9b25SBarry Smith .vb
112630de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
112730de9b25SBarry Smith                "MyColor",MyColor);
112830de9b25SBarry Smith .ve
112930de9b25SBarry Smith 
113030de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
113130de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
113230de9b25SBarry Smith    or at runtime via the option
113330de9b25SBarry Smith $     -mat_coloring_type my_color
113430de9b25SBarry Smith 
1135ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
113630de9b25SBarry Smith 
113730de9b25SBarry Smith .keywords: matrix, Coloring, register
113830de9b25SBarry Smith 
113930de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
114030de9b25SBarry Smith M*/
1141aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1142f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1143b12f92e5SBarry Smith #else
1144f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1145b12f92e5SBarry Smith #endif
114630de9b25SBarry Smith 
11472bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1148f1144a30SSatish Balay 
1149be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1150be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1151be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1152639f9d9dSBarry Smith 
1153d9274352SBarry Smith /*S
1154d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1155d9274352SBarry Smith         and coloring
1156639f9d9dSBarry Smith 
1157d9274352SBarry Smith    Level: beginner
1158d9274352SBarry Smith 
1159d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1160d9274352SBarry Smith 
1161d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1162d9274352SBarry Smith S*/
1163e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1164639f9d9dSBarry Smith 
1165be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1166be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1167be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1168be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
11694a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1170be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1171be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring,PetscInt);
1172be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring,PetscInt*);
1173be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1174be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1175be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1176be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring);
1177be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1178be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1179639f9d9dSBarry Smith /*
11800752156aSBarry Smith     These routines are for partitioning matrices: currently used only
11813eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
11820752156aSBarry Smith */
1183ca161407SBarry Smith 
1184d9274352SBarry Smith /*S
1185d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1186d9274352SBarry Smith 
1187d9274352SBarry Smith    Level: beginner
1188d9274352SBarry Smith 
1189d9274352SBarry Smith   Concepts: partitioning
1190d9274352SBarry Smith 
1191743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1192d9274352SBarry Smith S*/
119391e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1194d9274352SBarry Smith 
1195d9274352SBarry Smith /*E
11965bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1197d9274352SBarry Smith        with an optional dynamic library name, for example
1198d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1199d9274352SBarry Smith 
1200d9274352SBarry Smith    Level: beginner
1201d9274352SBarry Smith 
1202b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1203d9274352SBarry Smith E*/
1204e5a9bf91SBarry Smith #define MatPartitioningType const char*
12058ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1206c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
12078ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
120871306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
120971306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
121071306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
121171306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
121271306c60Sjroman 
1213ca161407SBarry Smith 
1214be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
12152e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1216be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1217be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1218be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1219be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1220be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1221be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
12222aabb6bbSBarry Smith 
1223be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
122430de9b25SBarry Smith 
122530de9b25SBarry Smith /*MC
122630de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
122730de9b25SBarry Smith    matrix package.
122830de9b25SBarry Smith 
122930de9b25SBarry Smith    Synopsis:
1230c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
123130de9b25SBarry Smith 
123230de9b25SBarry Smith    Not Collective
123330de9b25SBarry Smith 
123430de9b25SBarry Smith    Input Parameters:
123530de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
123630de9b25SBarry Smith .  path - location of library where creation routine is
123730de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
123830de9b25SBarry Smith -  function - function pointer that creates the partitioning type
123930de9b25SBarry Smith 
124030de9b25SBarry Smith    Level: developer
124130de9b25SBarry Smith 
124230de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
124330de9b25SBarry Smith    is ignored.
124430de9b25SBarry Smith 
124530de9b25SBarry Smith    Sample usage:
124630de9b25SBarry Smith .vb
124730de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
124830de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
124930de9b25SBarry Smith .ve
125030de9b25SBarry Smith 
125130de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
125230de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
125330de9b25SBarry Smith    or at runtime via the option
125430de9b25SBarry Smith $     -mat_partitioning_type my_part
125530de9b25SBarry Smith 
1256ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
125730de9b25SBarry Smith 
125830de9b25SBarry Smith .keywords: matrix, partitioning, register
125930de9b25SBarry Smith 
126030de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
126130de9b25SBarry Smith M*/
1262aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1263f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
12642aabb6bbSBarry Smith #else
1265f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
12662aabb6bbSBarry Smith #endif
12672aabb6bbSBarry Smith 
12682bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1269f1144a30SSatish Balay 
1270be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1271be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
12722bad1931SBarry Smith 
1273be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1274be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1275be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1276ca161407SBarry Smith 
1277be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
12780752156aSBarry Smith 
1279be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1280be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
128171306c60Sjroman 
1282954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
1283be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
128471306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1285be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1286be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
128771306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1288be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1289be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1290be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
129171306c60Sjroman 
129271306c60Sjroman #define MP_PARTY_OPT "opt"
129371306c60Sjroman #define MP_PARTY_LIN "lin"
129471306c60Sjroman #define MP_PARTY_SCA "sca"
129571306c60Sjroman #define MP_PARTY_RAN "ran"
129671306c60Sjroman #define MP_PARTY_GBF "gbf"
129771306c60Sjroman #define MP_PARTY_GCF "gcf"
129871306c60Sjroman #define MP_PARTY_BUB "bub"
129971306c60Sjroman #define MP_PARTY_DEF "def"
1300be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
130171306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
130271306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
130371306c60Sjroman #define MP_PARTY_NONE "no"
1304be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1305be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1306be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1307be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
130871306c60Sjroman 
130971306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1310be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1311be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1312be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1313be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1314be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
131571306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1316be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1317be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1318be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
131971306c60Sjroman 
13200752156aSBarry Smith /*
13210a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1322d4fbbf0eSBarry Smith */
13231c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13241c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13251c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13261c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13271c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13287c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13297c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13301c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13311c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13327c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13337c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13341c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13351c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
13361c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
13371c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13381c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13391c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13401c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13411c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13421c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13431c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13441c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
13451c1c02c0SLois Curfman McInnes                MATOP_COMPRESS=22,
13461c1c02c0SLois Curfman McInnes                MATOP_SET_OPTION=23,
13471c1c02c0SLois Curfman McInnes                MATOP_ZERO_ENTRIES=24,
13481c1c02c0SLois Curfman McInnes                MATOP_ZERO_ROWS=25,
13491c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_SYMBOLIC=26,
13501c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_NUMERIC=27,
13511c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_SYMBOLIC=28,
13521c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_NUMERIC=29,
1353d643ce63SMatthew Knepley                MATOP_SETUP_PREALLOCATION=30,
1354d643ce63SMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=31,
1355d643ce63SMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=32,
1356d643ce63SMatthew Knepley                MATOP_GET_ARRAY=33,
1357d643ce63SMatthew Knepley                MATOP_RESTORE_ARRAY=34,
1358e26f1c3cSBarry Smith                MATOP_DUPLICATE=35,
1359d643ce63SMatthew Knepley                MATOP_FORWARD_SOLVE=36,
1360d643ce63SMatthew Knepley                MATOP_BACKWARD_SOLVE=37,
1361d643ce63SMatthew Knepley                MATOP_ILUFACTOR=38,
1362d643ce63SMatthew Knepley                MATOP_ICCFACTOR=39,
1363d643ce63SMatthew Knepley                MATOP_AXPY=40,
1364d643ce63SMatthew Knepley                MATOP_GET_SUBMATRICES=41,
1365d643ce63SMatthew Knepley                MATOP_INCREASE_OVERLAP=42,
1366d643ce63SMatthew Knepley                MATOP_GET_VALUES=43,
1367d643ce63SMatthew Knepley                MATOP_COPY=44,
1368d643ce63SMatthew Knepley                MATOP_PRINT_HELP=45,
1369d643ce63SMatthew Knepley                MATOP_SCALE=46,
1370d643ce63SMatthew Knepley                MATOP_SHIFT=47,
1371d643ce63SMatthew Knepley                MATOP_DIAGONAL_SHIFT=48,
1372d643ce63SMatthew Knepley                MATOP_ILUDT_FACTOR=49,
1373d643ce63SMatthew Knepley                MATOP_GET_BLOCK_SIZE=50,
1374d643ce63SMatthew Knepley                MATOP_GET_ROW_IJ=51,
1375d643ce63SMatthew Knepley                MATOP_RESTORE_ROW_IJ=52,
1376d643ce63SMatthew Knepley                MATOP_GET_COLUMN_IJ=53,
1377d643ce63SMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=54,
1378d643ce63SMatthew Knepley                MATOP_FDCOLORING_CREATE=55,
1379d643ce63SMatthew Knepley                MATOP_COLORING_PATCH=56,
1380d643ce63SMatthew Knepley                MATOP_SET_UNFACTORED=57,
1381d643ce63SMatthew Knepley                MATOP_PERMUTE=58,
1382d643ce63SMatthew Knepley                MATOP_SET_VALUES_BLOCKED=59,
1383d643ce63SMatthew Knepley                MATOP_GET_SUBMATRIX=60,
1384d643ce63SMatthew Knepley                MATOP_DESTROY=61,
1385d643ce63SMatthew Knepley                MATOP_VIEW=62,
1386d643ce63SMatthew Knepley                MATOP_GET_MAPS=63,
1387d643ce63SMatthew Knepley                MATOP_USE_SCALED_FORM=64,
1388d643ce63SMatthew Knepley                MATOP_SCALE_SYSTEM=65,
1389d643ce63SMatthew Knepley                MATOP_UNSCALE_SYSTEM=66,
1390ba16a8cbSKris Buschelman                MATOP_SET_LOCAL_TO_GLOBAL_MAP=67,
1391d643ce63SMatthew Knepley                MATOP_SET_VALUES_LOCAL=68,
1392d643ce63SMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=69,
1393d643ce63SMatthew Knepley                MATOP_GET_ROW_MAX=70,
1394d643ce63SMatthew Knepley                MATOP_CONVERT=71,
1395d643ce63SMatthew Knepley                MATOP_SET_COLORING=72,
1396d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1397d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1398d643ce63SMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1399d643ce63SMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1400ba16a8cbSKris Buschelman                MATOP_MULT_CON=77,
1401ba16a8cbSKris Buschelman                MATOP_MULT_TRANSPOSE_CON=78,
1402ba16a8cbSKris Buschelman                MATOP_ILU_FACTOR_SYMBOLIC_CON=79,
1403d643ce63SMatthew Knepley                MATOP_PERMUTE_SPARSIFY=80,
1404d643ce63SMatthew Knepley                MATOP_MULT_MULTIPLE=81,
140541acf15aSKris Buschelman                MATOP_SOLVE_MULTIPLE=82,
140641acf15aSKris Buschelman                MATOP_GET_INERTIA=83,
140741acf15aSKris Buschelman                MATOP_LOAD=84,
140841acf15aSKris Buschelman                MATOP_IS_SYMMETRIC=85,
140941acf15aSKris Buschelman                MATOP_IS_HERMITIAN=86,
141041acf15aSKris Buschelman                MATOP_IS_STRUCTURALLY_SYMMETRIC=87,
141141acf15aSKris Buschelman                MATOP_PB_RELAX=88,
141241acf15aSKris Buschelman                MATOP_GET_VECS=89,
141341acf15aSKris Buschelman                MATOP_MAT_MULT=90,
141441acf15aSKris Buschelman                MATOP_MAT_MULT_SYMBOLIC=91,
141541acf15aSKris Buschelman                MATOP_MAT_MULT_NUMERIC=92,
141641acf15aSKris Buschelman                MATOP_PTAP=93,
141741acf15aSKris Buschelman                MATOP_PTAP_SYMBOLIC=94,
1418bc011b1eSHong Zhang                MATOP_PTAP_NUMERIC=95,
1419bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE=96,
1420bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97,
14217ba1a0bfSKris Buschelman                MATOP_MAT_MULTTRANSPOSE_NUMERIC=98,
14227ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_SEQAIJ=99,
14237ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_SEQAIJ=100,
14247ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_MPIAIJ=101,
142587d4246cSBarry Smith                MATOP_PTAP_NUMERIC_MPIAIJ=102,
1426f5edf698SHong Zhang                MATOP_SET_VALUES_ROW = 105,
1427d5c63f83SSatish Balay                MATOP_GET_ROW_UTRIANGULAR=108,
14282bebee5dSHong Zhang                MATOP_RESTORE_ROW_UTRIANGULAR=109,
142969db28dcSHong Zhang                MATOP_MATSOLVE=110,
143069db28dcSHong Zhang                MATOP_GET_REDUNDANTMATRIX=111
1431fae171e0SBarry Smith              } MatOperation;
1432be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1433be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1434be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1435be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1436112a2221SBarry Smith 
143790ace30eSBarry Smith /*
143890ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
143990ace30eSBarry Smith  stored in a universal format. By changing the format with
1440fb9695e5SSatish Balay  PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will
144190ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
144290ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
144390ace30eSBarry Smith  read into matrices of the same time.
144490ace30eSBarry Smith */
144590ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
144690ace30eSBarry Smith 
1447be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
14491f4e1ec7SBarry Smith 
1450d9274352SBarry Smith /*S
1451d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1452d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1453d9274352SBarry Smith 
1454f7a9e4ceSBarry Smith    Level: advanced
1455d9274352SBarry Smith 
1456d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1457d9274352SBarry Smith 
14586e1639daSBarry Smith   Users manual sections:
14596e1639daSBarry Smith .   sec_singular
14606e1639daSBarry Smith 
1461d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1462d9274352SBarry Smith S*/
146374637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1464d9274352SBarry Smith 
1465be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1466281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*);
1467be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1468be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1469be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
1470be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat);
147174637425SBarry Smith 
1472be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1473be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1474be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
147546129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
14763f1d51d7SBarry Smith 
1477be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1478be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1479be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1480c4f061fbSSatish Balay 
1481be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1482b0a32e0cSBarry Smith 
1483be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
148404f1ad80SBarry Smith 
1485fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
14863ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec);
14873ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1488e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1489e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1490e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace);
1491e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1492e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat);
1493e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal);
1494e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt);
1495e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *);
1496e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat);
1497e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1498e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1499e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal);
1500e884886fSBarry Smith 
1501e884886fSBarry Smith 
1502e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1503e884886fSBarry Smith 
1504e884886fSBarry Smith /*E
1505e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1506e884886fSBarry Smith 
1507e884886fSBarry Smith    Level: beginner
1508e884886fSBarry Smith 
1509e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1510e884886fSBarry Smith E*/
1511e884886fSBarry Smith #define MatMFFDType const char*
1512e884886fSBarry Smith #define MATMFFD_DS  "ds"
1513e884886fSBarry Smith #define MATMFFD_WP  "wp"
1514e884886fSBarry Smith 
1515e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,MatMFFDType);
1516e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1517e884886fSBarry Smith 
1518e884886fSBarry Smith /*MC
1519e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1520e884886fSBarry Smith 
1521e884886fSBarry Smith    Synopsis:
1522e884886fSBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1523e884886fSBarry Smith 
1524e884886fSBarry Smith    Not Collective
1525e884886fSBarry Smith 
1526e884886fSBarry Smith    Input Parameters:
1527e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1528e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1529e884886fSBarry Smith .  name_create - name of routine to create method context
1530e884886fSBarry Smith -  routine_create - routine to create method context
1531e884886fSBarry Smith 
1532e884886fSBarry Smith    Level: developer
1533e884886fSBarry Smith 
1534e884886fSBarry Smith    Notes:
1535e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1536e884886fSBarry Smith 
1537e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1538e884886fSBarry Smith    is ignored.
1539e884886fSBarry Smith 
1540e884886fSBarry Smith    Sample usage:
1541e884886fSBarry Smith .vb
1542e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1543e884886fSBarry Smith                "MyHCreate",MyHCreate);
1544e884886fSBarry Smith .ve
1545e884886fSBarry Smith 
1546e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1547e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1548e884886fSBarry Smith    or at runtime via the option
1549e884886fSBarry Smith $     -snes_mf_type my_h
1550e884886fSBarry Smith 
1551e884886fSBarry Smith .keywords: MatMFFD, register
1552e884886fSBarry Smith 
1553e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1554e884886fSBarry Smith M*/
1555e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1556e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1557e884886fSBarry Smith #else
1558e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1559e884886fSBarry Smith #endif
1560e884886fSBarry Smith 
1561e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]);
1562e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void);
1563e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal);
1564e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth);
1565e884886fSBarry Smith 
1566e884886fSBarry Smith 
1567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1568be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
15697dbadf16SMatthew Knepley 
1570e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
15712eac72dbSBarry Smith #endif
1572