xref: /petsc/include/petscmat.h (revision 1cbb95d3a31c71eb14a99cc255a79c6d2d35ee66)
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 
104e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]);
105c06d978dSMatthew Knepley 
106f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
107a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
108a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
109f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
110f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,MatType);
111be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
112be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
113be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
114be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
11530de9b25SBarry Smith 
116f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
117f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
118f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
119f69a0ea3SMatthew Knepley 
12030de9b25SBarry Smith /*MC
12130de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
12230de9b25SBarry Smith 
12330de9b25SBarry Smith    Synopsis:
124c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
12530de9b25SBarry Smith 
12630de9b25SBarry Smith    Not Collective
12730de9b25SBarry Smith 
12830de9b25SBarry Smith    Input Parameters:
12930de9b25SBarry Smith +  name - name of a new user-defined matrix type
13030de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
13130de9b25SBarry Smith .  name_create - name of routine to create method context
13230de9b25SBarry Smith -  routine_create - routine to create method context
13330de9b25SBarry Smith 
13430de9b25SBarry Smith    Notes:
13530de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
13630de9b25SBarry Smith 
13730de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
13830de9b25SBarry Smith    is ignored.
13930de9b25SBarry Smith 
14030de9b25SBarry Smith    Sample usage:
14130de9b25SBarry Smith .vb
14230de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
14330de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
14430de9b25SBarry Smith .ve
14530de9b25SBarry Smith 
14630de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
14730de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
14830de9b25SBarry Smith    or at runtime via the option
14930de9b25SBarry Smith $     -mat_type my_mat
15030de9b25SBarry Smith 
15130de9b25SBarry Smith    Level: advanced
15230de9b25SBarry Smith 
153ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
15430de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
15530de9b25SBarry Smith 
15630de9b25SBarry Smith .keywords: Mat, register
15730de9b25SBarry Smith 
15830de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
15930de9b25SBarry Smith 
16030de9b25SBarry Smith M*/
161273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
162273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
163273d9f13SBarry Smith #else
164273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
16530de9b25SBarry Smith #endif
16630de9b25SBarry Smith 
1679c666560SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
168273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
169b0a32e0cSBarry Smith extern PetscFList MatList;
17028988994SBarry Smith 
171be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
172be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
173be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
174ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
175ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
176ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
177ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
178ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
179ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
180ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
181be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
182ba966639SSatish 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)
183ba966639SSatish 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)
184ba966639SSatish 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)
185ba966639SSatish 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)
186ba966639SSatish 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))
187ba966639SSatish 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))
188ba966639SSatish 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))
189ba966639SSatish 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)
190ba966639SSatish 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)
191ba966639SSatish 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)
192ba966639SSatish 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)
193ba966639SSatish 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))
194ba966639SSatish 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))
195ba966639SSatish 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))
19663ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
1978d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
1988d7a6e47SBarry Smith 
199be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
200be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*);
201be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*);
202be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
203ba966639SSatish 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)
204ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
205ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
206ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
207ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
208ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
209ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
210be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
211ba966639SSatish 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)
212ba966639SSatish 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)
213ba966639SSatish 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)
214ba966639SSatish 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)
215ba966639SSatish 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))
216ba966639SSatish 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))
217ba966639SSatish 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))
218ba966639SSatish 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)
219ba966639SSatish 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)
220ba966639SSatish 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)
221ba966639SSatish 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)
222ba966639SSatish 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))
223ba966639SSatish 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))
224ba966639SSatish 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))
225be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
226be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
227ba966639SSatish 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)
228ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
229ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
230ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
231ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
232ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
233ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
234be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
235ba966639SSatish 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)
236ba966639SSatish 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)
237ba966639SSatish 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)
238ba966639SSatish 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)
239ba966639SSatish 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))
240ba966639SSatish 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))
241ba966639SSatish 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))
242ba966639SSatish 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)
243ba966639SSatish 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)
244ba966639SSatish 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)
245ba966639SSatish 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)
246ba966639SSatish 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))
247ba966639SSatish 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))
248ba966639SSatish 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))
249be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
250ba966639SSatish 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)
251ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
252be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
253be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
254ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
255ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*);
256284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
2571472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
2581472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
2592a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*);
2602a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter);
2612a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*);
2628cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
263793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat);
264b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat);
265793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2665a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*);
2671472f72bSBarry Smith 
268f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
269be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
27021c89e3eSBarry Smith 
2710c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
27299cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat);
27399cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat);
27499cafbc1SBarry Smith 
2758ed539a5SBarry Smith /* ------------------------------------------------------------*/
276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
277be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
27887d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
279f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
28084cb2905SBarry Smith 
2812ef4de8bSBarry Smith /*S
2822ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
2832ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
2842ef4de8bSBarry Smith 
2852ef4de8bSBarry Smith    Level: beginner
2862ef4de8bSBarry Smith 
2872ef4de8bSBarry Smith   Concepts: matrix; linear operator
2882ef4de8bSBarry Smith 
289d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
2902ef4de8bSBarry Smith S*/
291435da068SBarry Smith typedef struct {
292c1ac3661SBarry Smith   PetscInt k,j,i,c;
293435da068SBarry Smith } MatStencil;
2942ef4de8bSBarry Smith 
295be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
297be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
298435da068SBarry Smith 
299be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
300be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
301be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
3023a7fca6bSBarry Smith 
303d91e6319SBarry Smith /*E
304d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
305d91e6319SBarry Smith      to continue to add values to it
306d91e6319SBarry Smith 
307d91e6319SBarry Smith     Level: beginner
308d91e6319SBarry Smith 
309d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
310d91e6319SBarry Smith E*/
3116d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
312be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
313be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
314be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
3154f9c727eSBarry Smith 
316be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt    MatSetValue_Row;
317be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt    MatSetValue_Column;
318be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscScalar MatSetValue_Value;
3191d73ed98SBarry Smith 
32030de9b25SBarry Smith /*MC
32130de9b25SBarry Smith    MatSetValue - Set a single entry into a matrix.
32230de9b25SBarry Smith 
32330de9b25SBarry Smith    Synopsis:
324c1ac3661SBarry Smith    PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode);
32530de9b25SBarry Smith 
32630de9b25SBarry Smith    Not collective
32730de9b25SBarry Smith 
32830de9b25SBarry Smith    Input Parameters:
32930de9b25SBarry Smith +  m - the matrix
33030de9b25SBarry Smith .  row - the row location of the entry
33130de9b25SBarry Smith .  col - the column location of the entry
33230de9b25SBarry Smith .  value - the value to insert
33330de9b25SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
33430de9b25SBarry Smith 
33530de9b25SBarry Smith    Notes:
33630de9b25SBarry Smith    For efficiency one should use MatSetValues() and set several or many
33730de9b25SBarry Smith    values simultaneously if possible.
33830de9b25SBarry Smith 
33930de9b25SBarry Smith    Level: beginner
34030de9b25SBarry Smith 
34130de9b25SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
34230de9b25SBarry Smith M*/
343b951964fSBarry Smith #define MatSetValue(v,i,j,va,mode) \
344f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \
3451d73ed98SBarry Smith    MatSetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode))
34630de9b25SBarry Smith 
347ea06a074SBarry Smith #define MatGetValue(v,i,j,va) \
348f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,0) || \
3491d73ed98SBarry Smith    MatGetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&va))
35030de9b25SBarry Smith 
351d91e6319SBarry Smith #define MatSetValueLocal(v,i,j,va,mode) \
352f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \
3531d73ed98SBarry Smith    MatSetValuesLocal(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode))
35430de9b25SBarry Smith 
355d91e6319SBarry Smith /*E
356d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
357d91e6319SBarry Smith 
358d91e6319SBarry Smith     Level: beginner
359d91e6319SBarry Smith 
3600a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
361d91e6319SBarry Smith 
362d91e6319SBarry Smith .seealso: MatSetOption()
363d91e6319SBarry Smith E*/
364290bbb0aSBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_COLUMN_ORIENTED,MAT_ROWS_SORTED,
365945da6c4SHong Zhang /*3-4*/       MAT_COLUMNS_SORTED,MAT_NO_NEW_NONZERO_LOCATIONS,
366945da6c4SHong Zhang /*5-6*/       MAT_YES_NEW_NONZERO_LOCATIONS,MAT_SYMMETRIC,
367945da6c4SHong Zhang /*7-8*/       MAT_STRUCTURALLY_SYMMETRIC,MAT_NO_NEW_DIAGONALS,
368945da6c4SHong Zhang /*9-10*/      MAT_YES_NEW_DIAGONALS,MAT_INODE_LIMIT_1,
369945da6c4SHong Zhang /*11-13*/     MAT_INODE_LIMIT_2,MAT_INODE_LIMIT_3,MAT_INODE_LIMIT_4,
370945da6c4SHong Zhang /*14-15*/     MAT_INODE_LIMIT_5,MAT_IGNORE_OFF_PROC_ENTRIES,
371945da6c4SHong Zhang /*16-17*/     MAT_ROWS_UNSORTED,MAT_COLUMNS_UNSORTED,
372945da6c4SHong Zhang /*18*/        MAT_NEW_NONZERO_LOCATION_ERR,
373945da6c4SHong Zhang /*19-20*/     MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
374945da6c4SHong Zhang /*21-22*/     MAT_KEEP_ZEROED_ROWS,MAT_IGNORE_ZERO_ENTRIES,
375945da6c4SHong Zhang               MAT_USE_INODES,MAT_DO_NOT_USE_INODES,
376945da6c4SHong Zhang /*25-26*/     MAT_NOT_SYMMETRIC,MAT_HERMITIAN,
377290bbb0aSBarry Smith               MAT_NOT_STRUCTURALLY_SYMMETRIC,MAT_NOT_HERMITIAN,
378945da6c4SHong Zhang /*29-30*/     MAT_SYMMETRY_ETERNAL,MAT_NOT_SYMMETRY_ETERNAL,
379290bbb0aSBarry Smith               MAT_USE_COMPRESSEDROW,MAT_DO_NOT_USE_COMPRESSEDROW,
380945da6c4SHong Zhang /*33-34*/     MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
3810df259c2SBarry Smith               MAT_GETROW_UPPERTRIANGULAR,MAT_DO_NOT_IGNORE_ZERO_ENTRIES} MatOption;
382290bbb0aSBarry Smith extern const char *MatOptions[];
383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption);
384be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,MatType*);
385ba966639SSatish Balay PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),MatType,t)
38684cb2905SBarry Smith 
387be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
388be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
389be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
390f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat);
391f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat);
392be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
393be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
394be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
395be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
396ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
398be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
399ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
400be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
4017b80b807SBarry Smith 
402be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
403be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
404ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
405be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
406be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
407ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
408e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t)
409*1cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*);
410be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
411ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
412be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
413be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
4142bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat);
415f5edf698SHong Zhang 
416d91e6319SBarry Smith /*E
417d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
418d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
419d91e6319SBarry Smith 
420d91e6319SBarry Smith     Level: beginner
421d91e6319SBarry Smith 
422d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
423d91e6319SBarry Smith 
424d91e6319SBarry Smith .seealso: MatDuplicate()
425d91e6319SBarry Smith E*/
4262e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
4272e8a6d31SBarry Smith 
428f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,MatType,MatReuse,Mat*);
429ba966639SSatish Balay PetscPolymorphicFunction(MatConvert,(Mat A,MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
430be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
431ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
432ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
43394a9d846SBarry Smith 
434d91e6319SBarry Smith /*E
435d91e6319SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
436d91e6319SBarry Smith 
437d91e6319SBarry Smith     Level: beginner
438d91e6319SBarry Smith 
439d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
440d91e6319SBarry Smith 
44194b7f48cSBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
442d91e6319SBarry Smith E*/
443c537a176SHong Zhang typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
444cb5b572fSBarry Smith 
445be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
446be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
447be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
448ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
449e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t)
450be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
451ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
452*1cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*);
453*1cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t)
454be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
455be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
456f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,MatType,Mat*);
457ba966639SSatish Balay PetscPolymorphicFunction(MatLoad,(PetscViewer v,MatType t),(v,t,&a),Mat,a)
4587b80b807SBarry Smith 
4598f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4608f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
4618f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4628f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
463d4fbbf0eSBarry Smith 
464d91e6319SBarry Smith /*S
465d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
466d91e6319SBarry Smith 
467d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
468d91e6319SBarry Smith 
469d91e6319SBarry Smith    Level: intermediate
470d91e6319SBarry Smith 
471d91e6319SBarry Smith   Concepts: matrix^nonzero information
472d91e6319SBarry Smith 
473d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
474d91e6319SBarry Smith S*/
4754e220ebcSLois Curfman McInnes typedef struct {
476b0a32e0cSBarry Smith   PetscLogDouble rows_global,columns_global;         /* number of global rows and columns */
477b0a32e0cSBarry Smith   PetscLogDouble rows_local,columns_local;           /* number of local rows and columns */
478b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
479b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
480b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
481b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
482b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
483b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
484b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4854e220ebcSLois Curfman McInnes } MatInfo;
4864e220ebcSLois Curfman McInnes 
487d9274352SBarry Smith /*E
488d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
489d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
490d9274352SBarry Smith 
491d9274352SBarry Smith     Level: beginner
492d9274352SBarry Smith 
493d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
494d9274352SBarry Smith 
495d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
496d9274352SBarry Smith E*/
4977b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
499be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
501985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]);
502985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]);
503985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
50486697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec);
505be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,Mat*);
506ba966639SSatish Balay PetscPolymorphicFunction(MatTranspose,(Mat A),(A,&t),Mat,t)
507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
508ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
509be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
513ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
515be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
517be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
5187b80b807SBarry Smith 
519be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
520ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
521be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
522f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
523f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
524f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
525f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5267b80b807SBarry Smith 
527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
529be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5305ef9f2a5SBarry Smith 
531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5348ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**);
535f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
5367b80b807SBarry Smith 
537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
538be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *);
540abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *);
541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
542be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
543be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
544be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
547be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscScalar**,Mat*);
55043eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
551cd116e26SSatish Balay #include "petscctable.h"
55243eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
55343eb5e2fSMatthew Knepley #else
55443eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
55543eb5e2fSMatthew Knepley #endif
5568efafbd8SBarry Smith 
557be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5587b80b807SBarry Smith 
559be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
561be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
56222440eb1SKris Buschelman 
563be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
565be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
56622440eb1SKris Buschelman 
567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
568be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
569be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
570bc011b1eSHong Zhang 
571f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
57204aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure);
573be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat);
5741c741599SHong Zhang 
575f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
576f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5777b80b807SBarry Smith 
578be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
579be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
580f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
581f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
582be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
583be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
584052efed2SBarry Smith 
585be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
586be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
58790f02eecSBarry Smith 
588be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
5890c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
590ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
591be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
592be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
59369db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
594bd481603SBarry Smith 
595bd481603SBarry Smith /*MC
5960d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
597bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
598bd481603SBarry Smith 
599bd481603SBarry Smith    Synopsis:
600c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
601bd481603SBarry Smith 
602bd481603SBarry Smith    Collective on MPI_Comm
603bd481603SBarry Smith 
604bd481603SBarry Smith    Input Parameters:
605bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
606859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
607859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
608bd481603SBarry Smith 
609bd481603SBarry Smith    Output Parameters:
610bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
611bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
612bd481603SBarry Smith 
613bd481603SBarry Smith 
614bd481603SBarry Smith    Level: intermediate
615bd481603SBarry Smith 
616bd481603SBarry Smith    Notes:
617bd481603SBarry Smith    See the chapter in the users manual on performance for more details
618bd481603SBarry Smith 
6191d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
620bd481603SBarry Smith 
621bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
622bd481603SBarry Smith 
623bd481603SBarry Smith   Concepts: preallocation^Matrix
624bd481603SBarry Smith 
625bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
626bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
627bd481603SBarry Smith M*/
628c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6297c922b88SBarry Smith { \
630c1ac3661SBarry Smith   PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
631c1ac3661SBarry Smith   _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\
632c1ac3661SBarry Smith   _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\
633c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
634c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp;
6357c922b88SBarry Smith 
636bd481603SBarry Smith /*MC
6370d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
638bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
639bd481603SBarry Smith 
640bd481603SBarry Smith    Synopsis:
641c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
642bd481603SBarry Smith 
643bd481603SBarry Smith    Collective on MPI_Comm
644bd481603SBarry Smith 
645bd481603SBarry Smith    Input Parameters:
646bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
647859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
648859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
649bd481603SBarry Smith 
650bd481603SBarry Smith    Output Parameters:
651bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
652bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
653bd481603SBarry Smith 
654bd481603SBarry Smith 
655bd481603SBarry Smith    Level: intermediate
656bd481603SBarry Smith 
657bd481603SBarry Smith    Notes:
658bd481603SBarry Smith    See the chapter in the users manual on performance for more details
659bd481603SBarry Smith 
6601d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
661bd481603SBarry Smith 
662bd481603SBarry Smith   Concepts: preallocation^Matrix
663bd481603SBarry Smith 
664bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
665bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
666bd481603SBarry Smith M*/
667222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
668222b16d4SBarry Smith { \
669c1ac3661SBarry Smith   PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__end; \
670c1ac3661SBarry Smith   _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\
671c1ac3661SBarry Smith   _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\
672c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
673c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp;
674222b16d4SBarry Smith 
675bd481603SBarry Smith /*MC
676bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
677bd481603SBarry Smith        inserted using a local number of the rows and columns
678bd481603SBarry Smith 
679bd481603SBarry Smith    Synopsis:
680c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
681bd481603SBarry Smith 
682bd481603SBarry Smith    Not Collective
683bd481603SBarry Smith 
684bd481603SBarry Smith    Input Parameters:
685bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
686bd481603SBarry Smith .  nrows - the number of rows indicated
6871d73ed98SBarry Smith .  rows - the indices of the rows
688bd481603SBarry Smith .  ncols - the number of columns in the matrix
689bd481603SBarry Smith .  cols - the columns indicated
690bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
691bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
692bd481603SBarry Smith 
693bd481603SBarry Smith 
694bd481603SBarry Smith    Level: intermediate
695bd481603SBarry Smith 
696bd481603SBarry Smith    Notes:
697bd481603SBarry Smith    See the chapter in the users manual on performance for more details
698bd481603SBarry Smith 
6991d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
700bd481603SBarry Smith 
701bd481603SBarry Smith   Concepts: preallocation^Matrix
702bd481603SBarry Smith 
703bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
704bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
705bd481603SBarry Smith M*/
706c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
707c4f061fbSSatish Balay {\
708c1ac3661SBarry Smith   PetscInt __l;\
709ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
710ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
711c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
712ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
713c4f061fbSSatish Balay   }\
714c4f061fbSSatish Balay }
715c4f061fbSSatish Balay 
716bd481603SBarry Smith /*MC
717bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
718bd481603SBarry Smith        inserted using a local number of the rows and columns
719bd481603SBarry Smith 
720bd481603SBarry Smith    Synopsis:
721c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
722bd481603SBarry Smith 
723bd481603SBarry Smith    Not Collective
724bd481603SBarry Smith 
725bd481603SBarry Smith    Input Parameters:
726bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
727bd481603SBarry Smith .  nrows - the number of rows indicated
7281d73ed98SBarry Smith .  rows - the indices of the rows
729bd481603SBarry Smith .  ncols - the number of columns in the matrix
730bd481603SBarry Smith .  cols - the columns indicated
731bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
732bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
733bd481603SBarry Smith 
734bd481603SBarry Smith 
735bd481603SBarry Smith    Level: intermediate
736bd481603SBarry Smith 
737bd481603SBarry Smith    Notes:
738bd481603SBarry Smith    See the chapter in the users manual on performance for more details
739bd481603SBarry Smith 
740bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
741bd481603SBarry Smith 
742bd481603SBarry Smith   Concepts: preallocation^Matrix
743bd481603SBarry Smith 
744bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
745bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
746bd481603SBarry Smith M*/
747d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
748d3d32019SBarry Smith {\
749c1ac3661SBarry Smith   PetscInt __l;\
750d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
751d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
752d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
753d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
754d3d32019SBarry Smith   }\
755d3d32019SBarry Smith }
756d3d32019SBarry Smith 
757bd481603SBarry Smith /*MC
758bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
759bd481603SBarry Smith        inserted using a local number of the rows and columns
760bd481603SBarry Smith 
761bd481603SBarry Smith    Synopsis:
762c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
763bd481603SBarry Smith 
764bd481603SBarry Smith    Not Collective
765bd481603SBarry Smith 
766bd481603SBarry Smith    Input Parameters:
76764075487SBarry Smith +  row - the row
768bd481603SBarry Smith .  ncols - the number of columns in the matrix
769a50f8bf6SHong Zhang -  cols - the columns indicated
770a50f8bf6SHong Zhang 
771a50f8bf6SHong Zhang    Output Parameters:
772a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
773bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
774bd481603SBarry Smith 
775bd481603SBarry Smith 
776bd481603SBarry Smith    Level: intermediate
777bd481603SBarry Smith 
778bd481603SBarry Smith    Notes:
779bd481603SBarry Smith    See the chapter in the users manual on performance for more details
780bd481603SBarry Smith 
781bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
782bd481603SBarry Smith 
783bd481603SBarry Smith   Concepts: preallocation^Matrix
784bd481603SBarry Smith 
785bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
786bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
787bd481603SBarry Smith M*/
788c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
789c1ac3661SBarry Smith { PetscInt __i; \
7908ee2e534SBarry Smith   if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
7918ee2e534SBarry 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);\
7927c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
79364075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
7947cc688d7SBarry Smith     else dnz[row - __rstart]++;\
7957c922b88SBarry Smith   }\
7967c922b88SBarry Smith }
7977c922b88SBarry Smith 
798bd481603SBarry Smith /*MC
799bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
800bd481603SBarry Smith        inserted using a local number of the rows and columns
801bd481603SBarry Smith 
802bd481603SBarry Smith    Synopsis:
803c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
804bd481603SBarry Smith 
805bd481603SBarry Smith    Not Collective
806bd481603SBarry Smith 
807bd481603SBarry Smith    Input Parameters:
808bd481603SBarry Smith +  nrows - the number of rows indicated
8091d73ed98SBarry Smith .  rows - the indices of the rows
810bd481603SBarry Smith .  ncols - the number of columns in the matrix
811bd481603SBarry Smith .  cols - the columns indicated
812bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
813bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
814bd481603SBarry Smith 
815bd481603SBarry Smith 
816bd481603SBarry Smith    Level: intermediate
817bd481603SBarry Smith 
818bd481603SBarry Smith    Notes:
819bd481603SBarry Smith    See the chapter in the users manual on performance for more details
820bd481603SBarry Smith 
821bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
822bd481603SBarry Smith 
823bd481603SBarry Smith   Concepts: preallocation^Matrix
824bd481603SBarry Smith 
825bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
826bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
827bd481603SBarry Smith M*/
828d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
829c1ac3661SBarry Smith { PetscInt __i; \
830d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
831d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
832d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
833d3d32019SBarry Smith   }\
834d3d32019SBarry Smith }
835d3d32019SBarry Smith 
836bd481603SBarry Smith /*MC
8370d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
838bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
839bd481603SBarry Smith 
840bd481603SBarry Smith    Synopsis:
841c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
842bd481603SBarry Smith 
843bd481603SBarry Smith    Collective on MPI_Comm
844bd481603SBarry Smith 
845bd481603SBarry Smith    Input Parameters:
846bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
847bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
848bd481603SBarry Smith 
849bd481603SBarry Smith 
850bd481603SBarry Smith    Level: intermediate
851bd481603SBarry Smith 
852bd481603SBarry Smith    Notes:
853bd481603SBarry Smith    See the chapter in the users manual on performance for more details
854bd481603SBarry Smith 
855bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
856bd481603SBarry Smith 
857bd481603SBarry Smith   Concepts: preallocation^Matrix
858bd481603SBarry Smith 
859bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
860bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
861bd481603SBarry Smith M*/
862ef66eb69SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree(dnz);CHKERRQ(_4_ierr);}
8637c922b88SBarry Smith 
864bd481603SBarry Smith 
865bd481603SBarry Smith 
8667b80b807SBarry Smith /* Routines unique to particular data structures */
867be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
868ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
869be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBDiagGetData(Mat,PetscInt*,PetscInt*,PetscInt*[],PetscInt*[],PetscScalar***);
870698d4c6aSKris Buschelman 
871be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
872be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
873698d4c6aSKris Buschelman 
874be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
875be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
876be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
877c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
878c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
8797b80b807SBarry Smith 
880a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
881a96a251dSBarry Smith 
882be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
883ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
884be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
885ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
886be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
887ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
888be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDensePreallocation(Mat,PetscScalar[]);
889be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]);
890273d9f13SBarry Smith 
891be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
892ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
893be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
894be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
89553803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
896be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
897be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
898be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]);
899be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
900be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
901284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
902be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
903be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
904be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
905be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
906273d9f13SBarry Smith 
907be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
908ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*);
9091b807ce4Svictorle 
910be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
911be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
9122e8a6d31SBarry Smith 
913be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
9143a7fca6bSBarry Smith 
9157b80b807SBarry Smith /*
9167b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
91794b7f48cSBarry Smith   done through the KSP and PC interfaces.
9187b80b807SBarry Smith */
9197b80b807SBarry Smith 
920d9274352SBarry Smith /*E
921d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
922d9274352SBarry Smith        with an optional dynamic library name, for example
923d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
924d9274352SBarry Smith 
925d9274352SBarry Smith    Level: beginner
926d9274352SBarry Smith 
927e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
928e5a9bf91SBarry Smith 
929d9274352SBarry Smith .seealso: MatGetOrdering()
930d9274352SBarry Smith E*/
93149773a63SBarry Smith #define MatOrderingType char*
932b12f92e5SBarry Smith #define MATORDERING_NATURAL   "natural"
933b12f92e5SBarry Smith #define MATORDERING_ND        "nd"
934b12f92e5SBarry Smith #define MATORDERING_1WD       "1wd"
935b12f92e5SBarry Smith #define MATORDERING_RCM       "rcm"
936b12f92e5SBarry Smith #define MATORDERING_QMD       "qmd"
937b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength"
93862152c8bSBarry Smith #define MATORDERING_DSC_ND    "dsc_nd"
93962152c8bSBarry Smith #define MATORDERING_DSC_MMD   "dsc_mmd"
94062152c8bSBarry Smith #define MATORDERING_DSC_MDF   "dsc_mdf"
941c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
942c06d978dSMatthew Knepley #define MATORDERING_IDENTITY  "identity"
943c06d978dSMatthew Knepley #define MATORDERING_REVERSE   "reverse"
944b12f92e5SBarry Smith 
945be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
946be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
947be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
94830de9b25SBarry Smith 
94930de9b25SBarry Smith /*MC
95030de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
95130de9b25SBarry Smith                                matrix package.
95230de9b25SBarry Smith 
95330de9b25SBarry Smith    Synopsis:
954c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
95530de9b25SBarry Smith 
95630de9b25SBarry Smith    Not Collective
95730de9b25SBarry Smith 
95830de9b25SBarry Smith    Input Parameters:
95930de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
96030de9b25SBarry Smith .  path - location of library where creation routine is
96130de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
96230de9b25SBarry Smith -  function - function pointer that creates the ordering
96330de9b25SBarry Smith 
96430de9b25SBarry Smith    Level: developer
96530de9b25SBarry Smith 
96630de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
96730de9b25SBarry Smith    is ignored.
96830de9b25SBarry Smith 
96930de9b25SBarry Smith    Sample usage:
97030de9b25SBarry Smith .vb
97130de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
97230de9b25SBarry Smith                "MyOrder",MyOrder);
97330de9b25SBarry Smith .ve
97430de9b25SBarry Smith 
97530de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
97630de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
97730de9b25SBarry Smith    or at runtime via the option
9782401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
97930de9b25SBarry Smith 
980ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
98130de9b25SBarry Smith 
98230de9b25SBarry Smith .keywords: matrix, ordering, register
98330de9b25SBarry Smith 
98430de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
98530de9b25SBarry Smith M*/
986aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
987f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
988b12f92e5SBarry Smith #else
989f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
990b12f92e5SBarry Smith #endif
99130de9b25SBarry Smith 
992be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
993be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
9942bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
995b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
996d4fbbf0eSBarry Smith 
997be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
998a2ce50c7SBarry Smith 
999d91e6319SBarry Smith /*S
10002401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1001b00f7748SHong Zhang 
100215e8a5b3SHong Zhang    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE
1003b00f7748SHong Zhang 
100415e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1005b00f7748SHong Zhang 
1006b00f7748SHong Zhang    Level: developer
1007b00f7748SHong Zhang 
1008d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1009d7d82daaSBarry Smith           MatFactorInfoInitialize()
1010b00f7748SHong Zhang 
1011b00f7748SHong Zhang S*/
1012b00f7748SHong Zhang typedef struct {
10130a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
1014fbf22428SSatish Balay   PetscReal     shiftpd;        /* if true, shift until positive pivots */
10152cea7109SBarry Smith   PetscReal     shift_fraction; /* record shift fraction taken */
101615e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
101715e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1018b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
101915e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
102067e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1021348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1022bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1023bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
102415e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
102515e8a5b3SHong Zhang } MatFactorInfo;
1026ffa6d0a5SLois Curfman McInnes 
1027be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
1028be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,MatFactorInfo*);
1029be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1030be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,MatFactorInfo*,Mat*);
1031be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,MatFactorInfo*);
1032be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,MatFactorInfo*);
1033be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1034be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1035be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1036be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,MatFactorInfo*);
1037be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,MatFactorInfo*,Mat*);
1038be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,MatFactorInfo*,Mat *);
1039be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1040be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1041be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1042be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1043be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1044be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1045be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1046be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
10478ed539a5SBarry Smith 
1048be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1049bb5a7306SBarry Smith 
1050d91e6319SBarry Smith /*E
1051d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1052bb1eb677SSatish Balay 
1053d91e6319SBarry Smith     Level: beginner
1054d91e6319SBarry Smith 
1055d9274352SBarry Smith    May be bitwise ORd together
1056d9274352SBarry Smith 
1057d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1058d91e6319SBarry Smith 
10594e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10604e7234bfSBarry Smith 
1061a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax()
1062d91e6319SBarry Smith E*/
1063ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1064ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1065ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
106684cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1067be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1068be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10698ed539a5SBarry Smith 
1070d4fbbf0eSBarry Smith /*
1071639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1072639f9d9dSBarry Smith */
1073b12f92e5SBarry Smith 
1074d9274352SBarry Smith /*E
1075d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1076d9274352SBarry Smith        with an optional dynamic library name, for example
1077d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1078d9274352SBarry Smith 
1079d9274352SBarry Smith    Level: beginner
1080d9274352SBarry Smith 
1081d9274352SBarry Smith .seealso: MatGetColoring()
1082d9274352SBarry Smith E*/
1083e5a9bf91SBarry Smith #define MatColoringType const char*
1084b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1085b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1086b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1087b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1088b12f92e5SBarry Smith 
10892e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,MatColoringType,ISColoring*);
10902e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
109130de9b25SBarry Smith 
109230de9b25SBarry Smith /*MC
109330de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
109430de9b25SBarry Smith                                matrix package.
109530de9b25SBarry Smith 
109630de9b25SBarry Smith    Synopsis:
1097c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
109830de9b25SBarry Smith 
109930de9b25SBarry Smith    Not Collective
110030de9b25SBarry Smith 
110130de9b25SBarry Smith    Input Parameters:
110230de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
110330de9b25SBarry Smith .  path - location of library where creation routine is
110430de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
110530de9b25SBarry Smith -  function - function pointer that creates the coloring
110630de9b25SBarry Smith 
110730de9b25SBarry Smith    Level: developer
110830de9b25SBarry Smith 
110930de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
111030de9b25SBarry Smith    is ignored.
111130de9b25SBarry Smith 
111230de9b25SBarry Smith    Sample usage:
111330de9b25SBarry Smith .vb
111430de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
111530de9b25SBarry Smith                "MyColor",MyColor);
111630de9b25SBarry Smith .ve
111730de9b25SBarry Smith 
111830de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
111930de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
112030de9b25SBarry Smith    or at runtime via the option
112130de9b25SBarry Smith $     -mat_coloring_type my_color
112230de9b25SBarry Smith 
1123ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
112430de9b25SBarry Smith 
112530de9b25SBarry Smith .keywords: matrix, Coloring, register
112630de9b25SBarry Smith 
112730de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
112830de9b25SBarry Smith M*/
1129aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1130f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1131b12f92e5SBarry Smith #else
1132f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1133b12f92e5SBarry Smith #endif
113430de9b25SBarry Smith 
11352bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1136f1144a30SSatish Balay 
1137be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1138be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1139be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1140639f9d9dSBarry Smith 
1141d9274352SBarry Smith /*S
1142d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1143d9274352SBarry Smith         and coloring
1144639f9d9dSBarry Smith 
1145d9274352SBarry Smith    Level: beginner
1146d9274352SBarry Smith 
1147d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1148d9274352SBarry Smith 
1149d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1150d9274352SBarry Smith S*/
1151e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1152639f9d9dSBarry Smith 
1153be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1154be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1155be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1156be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
11574a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1158be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1159be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring,PetscInt);
1160be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring,PetscInt*);
1161be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1162be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1163be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1164be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring);
1165be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1166be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1167639f9d9dSBarry Smith /*
11680752156aSBarry Smith     These routines are for partitioning matrices: currently used only
11693eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
11700752156aSBarry Smith */
1171ca161407SBarry Smith 
1172d9274352SBarry Smith /*S
1173d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1174d9274352SBarry Smith 
1175d9274352SBarry Smith    Level: beginner
1176d9274352SBarry Smith 
1177d9274352SBarry Smith   Concepts: partitioning
1178d9274352SBarry Smith 
1179743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1180d9274352SBarry Smith S*/
118191e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1182d9274352SBarry Smith 
1183d9274352SBarry Smith /*E
11845bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1185d9274352SBarry Smith        with an optional dynamic library name, for example
1186d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1187d9274352SBarry Smith 
1188d9274352SBarry Smith    Level: beginner
1189d9274352SBarry Smith 
1190b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1191d9274352SBarry Smith E*/
1192e5a9bf91SBarry Smith #define MatPartitioningType const char*
11938ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1194c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
11958ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
119671306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
119771306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
119871306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
119971306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
120071306c60Sjroman 
1201ca161407SBarry Smith 
1202be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
12032e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1204be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1205be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1206be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1207be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1208be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1209be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
12102aabb6bbSBarry Smith 
1211be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
121230de9b25SBarry Smith 
121330de9b25SBarry Smith /*MC
121430de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
121530de9b25SBarry Smith    matrix package.
121630de9b25SBarry Smith 
121730de9b25SBarry Smith    Synopsis:
1218c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
121930de9b25SBarry Smith 
122030de9b25SBarry Smith    Not Collective
122130de9b25SBarry Smith 
122230de9b25SBarry Smith    Input Parameters:
122330de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
122430de9b25SBarry Smith .  path - location of library where creation routine is
122530de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
122630de9b25SBarry Smith -  function - function pointer that creates the partitioning type
122730de9b25SBarry Smith 
122830de9b25SBarry Smith    Level: developer
122930de9b25SBarry Smith 
123030de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
123130de9b25SBarry Smith    is ignored.
123230de9b25SBarry Smith 
123330de9b25SBarry Smith    Sample usage:
123430de9b25SBarry Smith .vb
123530de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
123630de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
123730de9b25SBarry Smith .ve
123830de9b25SBarry Smith 
123930de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
124030de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
124130de9b25SBarry Smith    or at runtime via the option
124230de9b25SBarry Smith $     -mat_partitioning_type my_part
124330de9b25SBarry Smith 
1244ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
124530de9b25SBarry Smith 
124630de9b25SBarry Smith .keywords: matrix, partitioning, register
124730de9b25SBarry Smith 
124830de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
124930de9b25SBarry Smith M*/
1250aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1251f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
12522aabb6bbSBarry Smith #else
1253f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
12542aabb6bbSBarry Smith #endif
12552aabb6bbSBarry Smith 
12562bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1257f1144a30SSatish Balay 
1258be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1259be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
12602bad1931SBarry Smith 
1261be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1262be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1263be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1264ca161407SBarry Smith 
1265be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
12660752156aSBarry Smith 
1267be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
126971306c60Sjroman 
1270954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
1271be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
127271306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1273be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1274be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
127571306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1277be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1278be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
127971306c60Sjroman 
128071306c60Sjroman #define MP_PARTY_OPT "opt"
128171306c60Sjroman #define MP_PARTY_LIN "lin"
128271306c60Sjroman #define MP_PARTY_SCA "sca"
128371306c60Sjroman #define MP_PARTY_RAN "ran"
128471306c60Sjroman #define MP_PARTY_GBF "gbf"
128571306c60Sjroman #define MP_PARTY_GCF "gcf"
128671306c60Sjroman #define MP_PARTY_BUB "bub"
128771306c60Sjroman #define MP_PARTY_DEF "def"
1288be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
128971306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
129071306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
129171306c60Sjroman #define MP_PARTY_NONE "no"
1292be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1294be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1295be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
129671306c60Sjroman 
129771306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1298be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1299be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1300be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1301be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1302be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
130371306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1304be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1305be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1306be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
130771306c60Sjroman 
13080752156aSBarry Smith /*
13090a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1310d4fbbf0eSBarry Smith */
13111c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13121c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13131c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13141c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13151c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13167c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13177c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13181c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13191c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13207c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13217c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13221c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13231c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
13241c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
13251c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13261c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13271c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13281c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13291c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13301c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13311c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13321c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
13331c1c02c0SLois Curfman McInnes                MATOP_COMPRESS=22,
13341c1c02c0SLois Curfman McInnes                MATOP_SET_OPTION=23,
13351c1c02c0SLois Curfman McInnes                MATOP_ZERO_ENTRIES=24,
13361c1c02c0SLois Curfman McInnes                MATOP_ZERO_ROWS=25,
13371c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_SYMBOLIC=26,
13381c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_NUMERIC=27,
13391c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_SYMBOLIC=28,
13401c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_NUMERIC=29,
1341d643ce63SMatthew Knepley                MATOP_SETUP_PREALLOCATION=30,
1342d643ce63SMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=31,
1343d643ce63SMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=32,
1344d643ce63SMatthew Knepley                MATOP_GET_ARRAY=33,
1345d643ce63SMatthew Knepley                MATOP_RESTORE_ARRAY=34,
1346e26f1c3cSBarry Smith                MATOP_DUPLICATE=35,
1347d643ce63SMatthew Knepley                MATOP_FORWARD_SOLVE=36,
1348d643ce63SMatthew Knepley                MATOP_BACKWARD_SOLVE=37,
1349d643ce63SMatthew Knepley                MATOP_ILUFACTOR=38,
1350d643ce63SMatthew Knepley                MATOP_ICCFACTOR=39,
1351d643ce63SMatthew Knepley                MATOP_AXPY=40,
1352d643ce63SMatthew Knepley                MATOP_GET_SUBMATRICES=41,
1353d643ce63SMatthew Knepley                MATOP_INCREASE_OVERLAP=42,
1354d643ce63SMatthew Knepley                MATOP_GET_VALUES=43,
1355d643ce63SMatthew Knepley                MATOP_COPY=44,
1356d643ce63SMatthew Knepley                MATOP_PRINT_HELP=45,
1357d643ce63SMatthew Knepley                MATOP_SCALE=46,
1358d643ce63SMatthew Knepley                MATOP_SHIFT=47,
1359d643ce63SMatthew Knepley                MATOP_DIAGONAL_SHIFT=48,
1360d643ce63SMatthew Knepley                MATOP_ILUDT_FACTOR=49,
1361d643ce63SMatthew Knepley                MATOP_GET_BLOCK_SIZE=50,
1362d643ce63SMatthew Knepley                MATOP_GET_ROW_IJ=51,
1363d643ce63SMatthew Knepley                MATOP_RESTORE_ROW_IJ=52,
1364d643ce63SMatthew Knepley                MATOP_GET_COLUMN_IJ=53,
1365d643ce63SMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=54,
1366d643ce63SMatthew Knepley                MATOP_FDCOLORING_CREATE=55,
1367d643ce63SMatthew Knepley                MATOP_COLORING_PATCH=56,
1368d643ce63SMatthew Knepley                MATOP_SET_UNFACTORED=57,
1369d643ce63SMatthew Knepley                MATOP_PERMUTE=58,
1370d643ce63SMatthew Knepley                MATOP_SET_VALUES_BLOCKED=59,
1371d643ce63SMatthew Knepley                MATOP_GET_SUBMATRIX=60,
1372d643ce63SMatthew Knepley                MATOP_DESTROY=61,
1373d643ce63SMatthew Knepley                MATOP_VIEW=62,
1374d643ce63SMatthew Knepley                MATOP_GET_MAPS=63,
1375d643ce63SMatthew Knepley                MATOP_USE_SCALED_FORM=64,
1376d643ce63SMatthew Knepley                MATOP_SCALE_SYSTEM=65,
1377d643ce63SMatthew Knepley                MATOP_UNSCALE_SYSTEM=66,
1378ba16a8cbSKris Buschelman                MATOP_SET_LOCAL_TO_GLOBAL_MAP=67,
1379d643ce63SMatthew Knepley                MATOP_SET_VALUES_LOCAL=68,
1380d643ce63SMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=69,
1381d643ce63SMatthew Knepley                MATOP_GET_ROW_MAX=70,
1382d643ce63SMatthew Knepley                MATOP_CONVERT=71,
1383d643ce63SMatthew Knepley                MATOP_SET_COLORING=72,
1384d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1385d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1386d643ce63SMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1387d643ce63SMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1388ba16a8cbSKris Buschelman                MATOP_MULT_CON=77,
1389ba16a8cbSKris Buschelman                MATOP_MULT_TRANSPOSE_CON=78,
1390ba16a8cbSKris Buschelman                MATOP_ILU_FACTOR_SYMBOLIC_CON=79,
1391d643ce63SMatthew Knepley                MATOP_PERMUTE_SPARSIFY=80,
1392d643ce63SMatthew Knepley                MATOP_MULT_MULTIPLE=81,
139341acf15aSKris Buschelman                MATOP_SOLVE_MULTIPLE=82,
139441acf15aSKris Buschelman                MATOP_GET_INERTIA=83,
139541acf15aSKris Buschelman                MATOP_LOAD=84,
139641acf15aSKris Buschelman                MATOP_IS_SYMMETRIC=85,
139741acf15aSKris Buschelman                MATOP_IS_HERMITIAN=86,
139841acf15aSKris Buschelman                MATOP_IS_STRUCTURALLY_SYMMETRIC=87,
139941acf15aSKris Buschelman                MATOP_PB_RELAX=88,
140041acf15aSKris Buschelman                MATOP_GET_VECS=89,
140141acf15aSKris Buschelman                MATOP_MAT_MULT=90,
140241acf15aSKris Buschelman                MATOP_MAT_MULT_SYMBOLIC=91,
140341acf15aSKris Buschelman                MATOP_MAT_MULT_NUMERIC=92,
140441acf15aSKris Buschelman                MATOP_PTAP=93,
140541acf15aSKris Buschelman                MATOP_PTAP_SYMBOLIC=94,
1406bc011b1eSHong Zhang                MATOP_PTAP_NUMERIC=95,
1407bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE=96,
1408bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97,
14097ba1a0bfSKris Buschelman                MATOP_MAT_MULTTRANSPOSE_NUMERIC=98,
14107ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_SEQAIJ=99,
14117ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_SEQAIJ=100,
14127ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_MPIAIJ=101,
141387d4246cSBarry Smith                MATOP_PTAP_NUMERIC_MPIAIJ=102,
1414f5edf698SHong Zhang                MATOP_SET_VALUES_ROW = 105,
1415d5c63f83SSatish Balay                MATOP_GET_ROW_UTRIANGULAR=108,
14162bebee5dSHong Zhang                MATOP_RESTORE_ROW_UTRIANGULAR=109,
141769db28dcSHong Zhang                MATOP_MATSOLVE=110,
141869db28dcSHong Zhang                MATOP_GET_REDUNDANTMATRIX=111
1419fae171e0SBarry Smith              } MatOperation;
1420be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1421be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1422be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1423be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1424112a2221SBarry Smith 
142590ace30eSBarry Smith /*
142690ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
142790ace30eSBarry Smith  stored in a universal format. By changing the format with
1428fb9695e5SSatish Balay  PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will
142990ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
143090ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
143190ace30eSBarry Smith  read into matrices of the same time.
143290ace30eSBarry Smith */
143390ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
143490ace30eSBarry Smith 
1435be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1436be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
14371f4e1ec7SBarry Smith 
1438d9274352SBarry Smith /*S
1439d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1440d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1441d9274352SBarry Smith 
1442f7a9e4ceSBarry Smith    Level: advanced
1443d9274352SBarry Smith 
1444d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1445d9274352SBarry Smith 
14466e1639daSBarry Smith   Users manual sections:
14476e1639daSBarry Smith .   sec_singular
14486e1639daSBarry Smith 
1449d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1450d9274352SBarry Smith S*/
145174637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1452d9274352SBarry Smith 
1453be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1454281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*);
1455be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1456be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1457be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
1458be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat);
145974637425SBarry Smith 
1460be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1461be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1462be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
146346129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
14643f1d51d7SBarry Smith 
1465be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1466be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1467be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1468c4f061fbSSatish Balay 
1469be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1470b0a32e0cSBarry Smith 
1471be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
147204f1ad80SBarry Smith 
1473fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
14743ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec);
14753ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1476e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1477e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1478e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace);
1479e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1480e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat);
1481e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal);
1482e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt);
1483e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *);
1484e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat);
1485e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1486e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1487e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal);
1488e884886fSBarry Smith 
1489e884886fSBarry Smith 
1490e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1491e884886fSBarry Smith 
1492e884886fSBarry Smith /*E
1493e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1494e884886fSBarry Smith 
1495e884886fSBarry Smith    Level: beginner
1496e884886fSBarry Smith 
1497e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1498e884886fSBarry Smith E*/
1499e884886fSBarry Smith #define MatMFFDType const char*
1500e884886fSBarry Smith #define MATMFFD_DS  "ds"
1501e884886fSBarry Smith #define MATMFFD_WP  "wp"
1502e884886fSBarry Smith 
1503e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,MatMFFDType);
1504e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1505e884886fSBarry Smith 
1506e884886fSBarry Smith /*MC
1507e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1508e884886fSBarry Smith 
1509e884886fSBarry Smith    Synopsis:
1510e884886fSBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1511e884886fSBarry Smith 
1512e884886fSBarry Smith    Not Collective
1513e884886fSBarry Smith 
1514e884886fSBarry Smith    Input Parameters:
1515e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1516e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1517e884886fSBarry Smith .  name_create - name of routine to create method context
1518e884886fSBarry Smith -  routine_create - routine to create method context
1519e884886fSBarry Smith 
1520e884886fSBarry Smith    Level: developer
1521e884886fSBarry Smith 
1522e884886fSBarry Smith    Notes:
1523e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1524e884886fSBarry Smith 
1525e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1526e884886fSBarry Smith    is ignored.
1527e884886fSBarry Smith 
1528e884886fSBarry Smith    Sample usage:
1529e884886fSBarry Smith .vb
1530e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1531e884886fSBarry Smith                "MyHCreate",MyHCreate);
1532e884886fSBarry Smith .ve
1533e884886fSBarry Smith 
1534e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1535e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1536e884886fSBarry Smith    or at runtime via the option
1537e884886fSBarry Smith $     -snes_mf_type my_h
1538e884886fSBarry Smith 
1539e884886fSBarry Smith .keywords: MatMFFD, register
1540e884886fSBarry Smith 
1541e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1542e884886fSBarry Smith M*/
1543e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1544e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1545e884886fSBarry Smith #else
1546e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1547e884886fSBarry Smith #endif
1548e884886fSBarry Smith 
1549e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]);
1550e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void);
1551e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal);
1552e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth);
1553e884886fSBarry Smith 
1554e884886fSBarry Smith 
1555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
15577dbadf16SMatthew Knepley 
1558e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
15592eac72dbSBarry Smith #endif
1560