xref: /petsc/include/petscmat.h (revision 5a7f1df34c47719bcb2902509fa26888bd3a7c48)
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*/
29273d9f13SBarry Smith #define MATSAME            "same"
30273d9f13SBarry Smith #define MATSEQMAIJ         "seqmaij"
31273d9f13SBarry Smith #define MATMPIMAIJ         "mpimaij"
32209238afSKris Buschelman #define MATMAIJ            "maij"
33273d9f13SBarry Smith #define MATIS              "is"
34273d9f13SBarry Smith #define MATMPIROWBS        "mpirowbs"
35273d9f13SBarry Smith #define MATSEQAIJ          "seqaij"
36273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
37209238afSKris Buschelman #define MATAIJ             "aij"
38273d9f13SBarry Smith #define MATSHELL           "shell"
39273d9f13SBarry Smith #define MATSEQBDIAG        "seqbdiag"
40273d9f13SBarry Smith #define MATMPIBDIAG        "mpibdiag"
41209238afSKris Buschelman #define MATBDIAG           "bdiag"
42209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
43273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
44209238afSKris Buschelman #define MATDENSE           "dense"
45273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
46273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
47209238afSKris Buschelman #define MATBAIJ            "baij"
48273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
49273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
50273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
51209238afSKris Buschelman #define MATSBAIJ           "sbaij"
52cebc7f6cSBarry Smith #define MATDAAD            "daad"
53cebc7f6cSBarry Smith #define MATMFFD            "mffd"
54c8a8475eSBarry Smith #define MATNORMAL          "normal"
55ab92ecdeSBarry Smith #define MATLRC             "lrc"
56b3a1e11cSKris Buschelman #define MATSEQAIJSPOOLES   "seqaijspooles"
57d10c748bSKris Buschelman #define MATMPIAIJSPOOLES   "mpiaijspooles"
589abb65ffSKris Buschelman #define MATSEQSBAIJSPOOLES "seqsbaijspooles"
5922191285SKris Buschelman #define MATMPISBAIJSPOOLES "mpisbaijspooles"
60bb4d4166SHong Zhang #define MATAIJSPOOLES      "aijspooles"
61bb4d4166SHong Zhang #define MATSBAIJSPOOLES    "sbaijspooles"
6214b396bbSKris Buschelman #define MATSUPERLU         "superlu"
631d96aa28SKris Buschelman #define MATSUPERLU_DIST    "superlu_dist"
641677d0b8SKris Buschelman #define MATUMFPACK         "umfpack"
65e8aa55a4SKris Buschelman #define MATESSL            "essl"
664eb8e494SKris Buschelman #define MATLUSOL           "lusol"
67397b6df1SKris Buschelman #define MATAIJMUMPS        "aijmumps"
68397b6df1SKris Buschelman #define MATSBAIJMUMPS      "sbaijmumps"
698da957c5SKris Buschelman #define MATDSCPACK         "dscpack"
707065e2aaSKris Buschelman #define MATMATLAB          "matlab"
714b8d99adSRichard Tran Mills #define MATSEQCSRPERM      "seqcsrperm"
7218def556SRichard Tran Mills #define MATMPICSRPERM      "mpicsrperm"
73814655b8SBarry Smith #define MATCSRPERM         "csrperm"
7481824310SBarry Smith #define MATSEQCRL          "seqcrl"
75c02b24c5SBarry Smith #define MATMPICRL          "mpicrl"
76c02b24c5SBarry Smith #define MATCRL             "crl"
77c0aa2d19SHong Zhang #define MATPLAPACK         "plapack"
782a6744ebSBarry Smith #define MATSCATTER         "scatter"
79421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
80793850ffSBarry Smith #define MATCOMPOSITE       "composite"
81*5a7f1df3SHong Zhang #define MATSEQFFTW         "seqfftw"
82f69a0ea3SMatthew Knepley #define MatType const char*
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;
90c06d978dSMatthew Knepley 
91ceb03754SKris Buschelman /*E
92ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
93ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
94ceb03754SKris Buschelman 
95ceb03754SKris Buschelman     Level: beginner
96ceb03754SKris Buschelman 
97ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
98ceb03754SKris Buschelman 
990c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
100ceb03754SKris Buschelman E*/
101ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse;
102ceb03754SKris Buschelman 
103be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(char *);
104c06d978dSMatthew Knepley 
105f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
106a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
107a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
108f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
109f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,MatType);
110be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
111be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
112be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
113be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
11430de9b25SBarry Smith 
115f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
116f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
117f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
118f69a0ea3SMatthew Knepley 
11930de9b25SBarry Smith /*MC
12030de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
12130de9b25SBarry Smith 
12230de9b25SBarry Smith    Synopsis:
123c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
12430de9b25SBarry Smith 
12530de9b25SBarry Smith    Not Collective
12630de9b25SBarry Smith 
12730de9b25SBarry Smith    Input Parameters:
12830de9b25SBarry Smith +  name - name of a new user-defined matrix type
12930de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
13030de9b25SBarry Smith .  name_create - name of routine to create method context
13130de9b25SBarry Smith -  routine_create - routine to create method context
13230de9b25SBarry Smith 
13330de9b25SBarry Smith    Notes:
13430de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
13530de9b25SBarry Smith 
13630de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
13730de9b25SBarry Smith    is ignored.
13830de9b25SBarry Smith 
13930de9b25SBarry Smith    Sample usage:
14030de9b25SBarry Smith .vb
14130de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
14230de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
14330de9b25SBarry Smith .ve
14430de9b25SBarry Smith 
14530de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
14630de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
14730de9b25SBarry Smith    or at runtime via the option
14830de9b25SBarry Smith $     -mat_type my_mat
14930de9b25SBarry Smith 
15030de9b25SBarry Smith    Level: advanced
15130de9b25SBarry Smith 
152ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
15330de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
15430de9b25SBarry Smith 
15530de9b25SBarry Smith .keywords: Mat, register
15630de9b25SBarry Smith 
15730de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
15830de9b25SBarry Smith 
15930de9b25SBarry Smith M*/
160273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
161273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
162273d9f13SBarry Smith #else
163273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
16430de9b25SBarry Smith #endif
16530de9b25SBarry Smith 
166273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
167b0a32e0cSBarry Smith extern PetscFList MatList;
16828988994SBarry Smith 
169be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
170be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
171be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
172ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
173ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
174ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
175ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
176ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
177ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
178ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
179be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
180ba966639SSatish 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)
181ba966639SSatish 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)
182ba966639SSatish 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)
183ba966639SSatish 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)
184ba966639SSatish 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))
185ba966639SSatish 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))
186ba966639SSatish 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))
187ba966639SSatish 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)
188ba966639SSatish 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)
189ba966639SSatish 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)
190ba966639SSatish 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)
191ba966639SSatish 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))
192ba966639SSatish 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))
193ba966639SSatish 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))
1942fb0ec9aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
1958d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
1968d7a6e47SBarry Smith 
197be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
198be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*);
199be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*);
200be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
201ba966639SSatish 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)
202ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
203ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
204ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
205ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
206ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
207ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
208be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
209ba966639SSatish 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)
210ba966639SSatish 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)
211ba966639SSatish 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)
212ba966639SSatish 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)
213ba966639SSatish 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))
214ba966639SSatish 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))
215ba966639SSatish 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))
216ba966639SSatish 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)
217ba966639SSatish 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)
218ba966639SSatish 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)
219ba966639SSatish 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)
220ba966639SSatish 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))
221ba966639SSatish 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))
222ba966639SSatish 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))
223be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
224be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
225ba966639SSatish 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)
226ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
227ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
228ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
229ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
230ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
231ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
232be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
233ba966639SSatish 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)
234ba966639SSatish 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)
235ba966639SSatish 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)
236ba966639SSatish 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)
237ba966639SSatish 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))
238ba966639SSatish 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))
239ba966639SSatish 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))
240ba966639SSatish 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)
241ba966639SSatish 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)
242ba966639SSatish 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)
243ba966639SSatish 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)
244ba966639SSatish 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))
245ba966639SSatish 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))
246ba966639SSatish 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))
247be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
248ba966639SSatish 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)
249ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
250be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
251be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
252ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
253ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*);
254284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
2551472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
2561472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
2572a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*);
2582a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter);
2592a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*);
2608cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
261793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat);
262793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
263*5a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*);
2641472f72bSBarry Smith 
265f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
266be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
26721c89e3eSBarry Smith 
2680c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
26999cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat);
27099cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat);
27199cafbc1SBarry Smith 
2728ed539a5SBarry Smith /* ------------------------------------------------------------*/
273be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
274be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
27587d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
276f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
27784cb2905SBarry Smith 
2782ef4de8bSBarry Smith /*S
2792ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
2802ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
2812ef4de8bSBarry Smith 
2822ef4de8bSBarry Smith    Level: beginner
2832ef4de8bSBarry Smith 
2842ef4de8bSBarry Smith   Concepts: matrix; linear operator
2852ef4de8bSBarry Smith 
286d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
2872ef4de8bSBarry Smith S*/
288435da068SBarry Smith typedef struct {
289c1ac3661SBarry Smith   PetscInt k,j,i,c;
290435da068SBarry Smith } MatStencil;
2912ef4de8bSBarry Smith 
292be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
294be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
295435da068SBarry Smith 
296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
297be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
298be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
2993a7fca6bSBarry Smith 
300d91e6319SBarry Smith /*E
301d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
302d91e6319SBarry Smith      to continue to add values to it
303d91e6319SBarry Smith 
304d91e6319SBarry Smith     Level: beginner
305d91e6319SBarry Smith 
306d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
307d91e6319SBarry Smith E*/
3086d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
309be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
310be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
311be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
3124f9c727eSBarry Smith 
313be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt    MatSetValue_Row;
314be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt    MatSetValue_Column;
315be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscScalar MatSetValue_Value;
3161d73ed98SBarry Smith 
31730de9b25SBarry Smith /*MC
31830de9b25SBarry Smith    MatSetValue - Set a single entry into a matrix.
31930de9b25SBarry Smith 
32030de9b25SBarry Smith    Synopsis:
321c1ac3661SBarry Smith    PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode);
32230de9b25SBarry Smith 
32330de9b25SBarry Smith    Not collective
32430de9b25SBarry Smith 
32530de9b25SBarry Smith    Input Parameters:
32630de9b25SBarry Smith +  m - the matrix
32730de9b25SBarry Smith .  row - the row location of the entry
32830de9b25SBarry Smith .  col - the column location of the entry
32930de9b25SBarry Smith .  value - the value to insert
33030de9b25SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
33130de9b25SBarry Smith 
33230de9b25SBarry Smith    Notes:
33330de9b25SBarry Smith    For efficiency one should use MatSetValues() and set several or many
33430de9b25SBarry Smith    values simultaneously if possible.
33530de9b25SBarry Smith 
33630de9b25SBarry Smith    Level: beginner
33730de9b25SBarry Smith 
33830de9b25SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
33930de9b25SBarry Smith M*/
340b951964fSBarry Smith #define MatSetValue(v,i,j,va,mode) \
341f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \
3421d73ed98SBarry Smith    MatSetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode))
34330de9b25SBarry Smith 
344ea06a074SBarry Smith #define MatGetValue(v,i,j,va) \
345f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,0) || \
3461d73ed98SBarry Smith    MatGetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&va))
34730de9b25SBarry Smith 
348d91e6319SBarry Smith #define MatSetValueLocal(v,i,j,va,mode) \
349f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \
3501d73ed98SBarry Smith    MatSetValuesLocal(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode))
35130de9b25SBarry Smith 
352d91e6319SBarry Smith /*E
353d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
354d91e6319SBarry Smith 
355d91e6319SBarry Smith     Level: beginner
356d91e6319SBarry Smith 
3570a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
358d91e6319SBarry Smith 
359d91e6319SBarry Smith .seealso: MatSetOption()
360d91e6319SBarry Smith E*/
361290bbb0aSBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_COLUMN_ORIENTED,MAT_ROWS_SORTED,
362290bbb0aSBarry Smith               MAT_COLUMNS_SORTED,MAT_NO_NEW_NONZERO_LOCATIONS,
363290bbb0aSBarry Smith               MAT_YES_NEW_NONZERO_LOCATIONS,MAT_SYMMETRIC,
364290bbb0aSBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,MAT_NO_NEW_DIAGONALS,
365290bbb0aSBarry Smith               MAT_YES_NEW_DIAGONALS,MAT_INODE_LIMIT_1,MAT_INODE_LIMIT_2,
366290bbb0aSBarry Smith               MAT_INODE_LIMIT_3,MAT_INODE_LIMIT_4,MAT_INODE_LIMIT_5,
367290bbb0aSBarry Smith               MAT_IGNORE_OFF_PROC_ENTRIES,MAT_ROWS_UNSORTED,
368290bbb0aSBarry Smith               MAT_COLUMNS_UNSORTED,MAT_NEW_NONZERO_LOCATION_ERR,
369290bbb0aSBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
370290bbb0aSBarry Smith               MAT_KEEP_ZEROED_ROWS,MAT_IGNORE_ZERO_ENTRIES,MAT_USE_INODES,
371290bbb0aSBarry Smith               MAT_DO_NOT_USE_INODES,MAT_NOT_SYMMETRIC,MAT_HERMITIAN,
372290bbb0aSBarry Smith               MAT_NOT_STRUCTURALLY_SYMMETRIC,MAT_NOT_HERMITIAN,
373290bbb0aSBarry Smith               MAT_SYMMETRY_ETERNAL,MAT_NOT_SYMMETRY_ETERNAL,
374290bbb0aSBarry Smith               MAT_USE_COMPRESSEDROW,MAT_DO_NOT_USE_COMPRESSEDROW,
375290bbb0aSBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,MAT_GETROW_UPPERTRIANGULAR} MatOption;
376290bbb0aSBarry Smith extern const char *MatOptions[];
377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption);
378be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,MatType*);
379ba966639SSatish Balay PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),MatType,t)
38084cb2905SBarry Smith 
381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
382be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
384f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat);
385f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat);
386be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
387be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
388be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
389be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
390ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
391be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
392be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
393ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
394be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
3957b80b807SBarry Smith 
396be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
398ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
400be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
401ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
402e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t)
403be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
404ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
405be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
406be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
4072bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat);
408f5edf698SHong Zhang 
409d91e6319SBarry Smith /*E
410d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
411d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
412d91e6319SBarry Smith 
413d91e6319SBarry Smith     Level: beginner
414d91e6319SBarry Smith 
415d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
416d91e6319SBarry Smith 
417d91e6319SBarry Smith .seealso: MatDuplicate()
418d91e6319SBarry Smith E*/
4192e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
4202e8a6d31SBarry Smith 
421f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,MatType,MatReuse,Mat*);
422ba966639SSatish Balay PetscPolymorphicFunction(MatConvert,(Mat A,MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
423be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
424ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
425ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
42694a9d846SBarry Smith 
427d91e6319SBarry Smith /*E
428d91e6319SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
429d91e6319SBarry Smith 
430d91e6319SBarry Smith     Level: beginner
431d91e6319SBarry Smith 
432d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
433d91e6319SBarry Smith 
43494b7f48cSBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
435d91e6319SBarry Smith E*/
436c537a176SHong Zhang typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
437cb5b572fSBarry Smith 
438be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
439be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
440be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
441ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
442e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t)
443be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
444ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
445be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscTruth*);
446ba966639SSatish Balay PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,&t),PetscTruth,t)
447be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
449f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,MatType,Mat*);
450ba966639SSatish Balay PetscPolymorphicFunction(MatLoad,(PetscViewer v,MatType t),(v,t,&a),Mat,a)
4517b80b807SBarry Smith 
452be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
453be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
454be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
455be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
456d4fbbf0eSBarry Smith 
457d91e6319SBarry Smith /*S
458d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
459d91e6319SBarry Smith 
460d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
461d91e6319SBarry Smith 
462d91e6319SBarry Smith    Level: intermediate
463d91e6319SBarry Smith 
464d91e6319SBarry Smith   Concepts: matrix^nonzero information
465d91e6319SBarry Smith 
466d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
467d91e6319SBarry Smith S*/
4684e220ebcSLois Curfman McInnes typedef struct {
469b0a32e0cSBarry Smith   PetscLogDouble rows_global,columns_global;         /* number of global rows and columns */
470b0a32e0cSBarry Smith   PetscLogDouble rows_local,columns_local;           /* number of local rows and columns */
471b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
472b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
473b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
474b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
475b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
476b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
477b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4784e220ebcSLois Curfman McInnes } MatInfo;
4794e220ebcSLois Curfman McInnes 
480d9274352SBarry Smith /*E
481d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
482d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
483d9274352SBarry Smith 
484d9274352SBarry Smith     Level: beginner
485d9274352SBarry Smith 
486d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
487d9274352SBarry Smith 
488d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
489d9274352SBarry Smith E*/
4907b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
491be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
492be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
493be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
494be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec);
495be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,Mat*);
496ba966639SSatish Balay PetscPolymorphicFunction(MatTranspose,(Mat A),(A,&t),Mat,t)
497be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
498ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
499be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
501be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
502be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
503ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
504be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
505be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
506be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
5087b80b807SBarry Smith 
509be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
510ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
512f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
513f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
514f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
515f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5167b80b807SBarry Smith 
517be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
519be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5205ef9f2a5SBarry Smith 
521be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
522be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
523be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5248ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**);
5257b80b807SBarry Smith 
526be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *);
529abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *);
530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
535be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
536be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
538be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscScalar**,Mat*);
53943eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
540cd116e26SSatish Balay #include "petscctable.h"
54143eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
54243eb5e2fSMatthew Knepley #else
54343eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
54443eb5e2fSMatthew Knepley #endif
5458efafbd8SBarry Smith 
546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5477b80b807SBarry Smith 
548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
55122440eb1SKris Buschelman 
552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
553be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
55522440eb1SKris Buschelman 
556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
557be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
558be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
559bc011b1eSHong Zhang 
560f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
56104aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure);
562be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat);
5631c741599SHong Zhang 
564f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
565f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5667b80b807SBarry Smith 
567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
568be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
569f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
570f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
571be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
572be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
573052efed2SBarry Smith 
574be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
575be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
57690f02eecSBarry Smith 
577be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
5780c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
579ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
580be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
581be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
582bd481603SBarry Smith 
583bd481603SBarry Smith /*MC
5840d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
585bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
586bd481603SBarry Smith 
587bd481603SBarry Smith    Synopsis:
588c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
589bd481603SBarry Smith 
590bd481603SBarry Smith    Collective on MPI_Comm
591bd481603SBarry Smith 
592bd481603SBarry Smith    Input Parameters:
593bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
594859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
595859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
596bd481603SBarry Smith 
597bd481603SBarry Smith    Output Parameters:
598bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
599bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
600bd481603SBarry Smith 
601bd481603SBarry Smith 
602bd481603SBarry Smith    Level: intermediate
603bd481603SBarry Smith 
604bd481603SBarry Smith    Notes:
605bd481603SBarry Smith    See the chapter in the users manual on performance for more details
606bd481603SBarry Smith 
6071d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
608bd481603SBarry Smith 
609bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
610bd481603SBarry Smith 
611bd481603SBarry Smith   Concepts: preallocation^Matrix
612bd481603SBarry Smith 
613bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
614bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
615bd481603SBarry Smith M*/
616c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6177c922b88SBarry Smith { \
618c1ac3661SBarry Smith   PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
619c1ac3661SBarry Smith   _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\
620c1ac3661SBarry Smith   _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\
621c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
622c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp;
6237c922b88SBarry Smith 
624bd481603SBarry Smith /*MC
6250d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
626bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
627bd481603SBarry Smith 
628bd481603SBarry Smith    Synopsis:
629c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
630bd481603SBarry Smith 
631bd481603SBarry Smith    Collective on MPI_Comm
632bd481603SBarry Smith 
633bd481603SBarry Smith    Input Parameters:
634bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
635859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
636859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
637bd481603SBarry Smith 
638bd481603SBarry Smith    Output Parameters:
639bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
640bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
641bd481603SBarry Smith 
642bd481603SBarry Smith 
643bd481603SBarry Smith    Level: intermediate
644bd481603SBarry Smith 
645bd481603SBarry Smith    Notes:
646bd481603SBarry Smith    See the chapter in the users manual on performance for more details
647bd481603SBarry Smith 
6481d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
649bd481603SBarry Smith 
650bd481603SBarry Smith   Concepts: preallocation^Matrix
651bd481603SBarry Smith 
652bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
653bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
654bd481603SBarry Smith M*/
655222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
656222b16d4SBarry Smith { \
657c1ac3661SBarry Smith   PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__end; \
658c1ac3661SBarry Smith   _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\
659c1ac3661SBarry Smith   _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\
660c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
661c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp;
662222b16d4SBarry Smith 
663bd481603SBarry Smith /*MC
664bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
665bd481603SBarry Smith        inserted using a local number of the rows and columns
666bd481603SBarry Smith 
667bd481603SBarry Smith    Synopsis:
668c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
669bd481603SBarry Smith 
670bd481603SBarry Smith    Not Collective
671bd481603SBarry Smith 
672bd481603SBarry Smith    Input Parameters:
673bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
674bd481603SBarry Smith .  nrows - the number of rows indicated
6751d73ed98SBarry Smith .  rows - the indices of the rows
676bd481603SBarry Smith .  ncols - the number of columns in the matrix
677bd481603SBarry Smith .  cols - the columns indicated
678bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
679bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
680bd481603SBarry Smith 
681bd481603SBarry Smith 
682bd481603SBarry Smith    Level: intermediate
683bd481603SBarry Smith 
684bd481603SBarry Smith    Notes:
685bd481603SBarry Smith    See the chapter in the users manual on performance for more details
686bd481603SBarry Smith 
6871d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
688bd481603SBarry Smith 
689bd481603SBarry Smith   Concepts: preallocation^Matrix
690bd481603SBarry Smith 
691bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
692bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
693bd481603SBarry Smith M*/
694c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
695c4f061fbSSatish Balay {\
696c1ac3661SBarry Smith   PetscInt __l;\
697ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
698ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
699c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
700ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
701c4f061fbSSatish Balay   }\
702c4f061fbSSatish Balay }
703c4f061fbSSatish Balay 
704bd481603SBarry Smith /*MC
705bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
706bd481603SBarry Smith        inserted using a local number of the rows and columns
707bd481603SBarry Smith 
708bd481603SBarry Smith    Synopsis:
709c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
710bd481603SBarry Smith 
711bd481603SBarry Smith    Not Collective
712bd481603SBarry Smith 
713bd481603SBarry Smith    Input Parameters:
714bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
715bd481603SBarry Smith .  nrows - the number of rows indicated
7161d73ed98SBarry Smith .  rows - the indices of the rows
717bd481603SBarry Smith .  ncols - the number of columns in the matrix
718bd481603SBarry Smith .  cols - the columns indicated
719bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
720bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
721bd481603SBarry Smith 
722bd481603SBarry Smith 
723bd481603SBarry Smith    Level: intermediate
724bd481603SBarry Smith 
725bd481603SBarry Smith    Notes:
726bd481603SBarry Smith    See the chapter in the users manual on performance for more details
727bd481603SBarry Smith 
728bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
729bd481603SBarry Smith 
730bd481603SBarry Smith   Concepts: preallocation^Matrix
731bd481603SBarry Smith 
732bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
733bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
734bd481603SBarry Smith M*/
735d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
736d3d32019SBarry Smith {\
737c1ac3661SBarry Smith   PetscInt __l;\
738d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
739d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
740d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
741d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
742d3d32019SBarry Smith   }\
743d3d32019SBarry Smith }
744d3d32019SBarry Smith 
745bd481603SBarry Smith /*MC
746bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
747bd481603SBarry Smith        inserted using a local number of the rows and columns
748bd481603SBarry Smith 
749bd481603SBarry Smith    Synopsis:
750c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
751bd481603SBarry Smith 
752bd481603SBarry Smith    Not Collective
753bd481603SBarry Smith 
754bd481603SBarry Smith    Input Parameters:
75564075487SBarry Smith +  row - the row
756bd481603SBarry Smith .  ncols - the number of columns in the matrix
757a50f8bf6SHong Zhang -  cols - the columns indicated
758a50f8bf6SHong Zhang 
759a50f8bf6SHong Zhang    Output Parameters:
760a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
761bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
762bd481603SBarry Smith 
763bd481603SBarry Smith 
764bd481603SBarry Smith    Level: intermediate
765bd481603SBarry Smith 
766bd481603SBarry Smith    Notes:
767bd481603SBarry Smith    See the chapter in the users manual on performance for more details
768bd481603SBarry Smith 
769bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
770bd481603SBarry Smith 
771bd481603SBarry Smith   Concepts: preallocation^Matrix
772bd481603SBarry Smith 
773bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
774bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
775bd481603SBarry Smith M*/
776c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
777c1ac3661SBarry Smith { PetscInt __i; \
7788ee2e534SBarry Smith   if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
7798ee2e534SBarry 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);\
7807c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
78164075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
7827cc688d7SBarry Smith     else dnz[row - __rstart]++;\
7837c922b88SBarry Smith   }\
7847c922b88SBarry Smith }
7857c922b88SBarry Smith 
786bd481603SBarry Smith /*MC
787bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
788bd481603SBarry Smith        inserted using a local number of the rows and columns
789bd481603SBarry Smith 
790bd481603SBarry Smith    Synopsis:
791c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
792bd481603SBarry Smith 
793bd481603SBarry Smith    Not Collective
794bd481603SBarry Smith 
795bd481603SBarry Smith    Input Parameters:
796bd481603SBarry Smith +  nrows - the number of rows indicated
7971d73ed98SBarry Smith .  rows - the indices of the rows
798bd481603SBarry Smith .  ncols - the number of columns in the matrix
799bd481603SBarry Smith .  cols - the columns indicated
800bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
801bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
802bd481603SBarry Smith 
803bd481603SBarry Smith 
804bd481603SBarry Smith    Level: intermediate
805bd481603SBarry Smith 
806bd481603SBarry Smith    Notes:
807bd481603SBarry Smith    See the chapter in the users manual on performance for more details
808bd481603SBarry Smith 
809bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
810bd481603SBarry Smith 
811bd481603SBarry Smith   Concepts: preallocation^Matrix
812bd481603SBarry Smith 
813bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
814bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
815bd481603SBarry Smith M*/
816d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
817c1ac3661SBarry Smith { PetscInt __i; \
818d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
819d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
820d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
821d3d32019SBarry Smith   }\
822d3d32019SBarry Smith }
823d3d32019SBarry Smith 
824bd481603SBarry Smith /*MC
8250d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
826bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
827bd481603SBarry Smith 
828bd481603SBarry Smith    Synopsis:
829c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
830bd481603SBarry Smith 
831bd481603SBarry Smith    Collective on MPI_Comm
832bd481603SBarry Smith 
833bd481603SBarry Smith    Input Parameters:
834bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
835bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
836bd481603SBarry Smith 
837bd481603SBarry Smith 
838bd481603SBarry Smith    Level: intermediate
839bd481603SBarry Smith 
840bd481603SBarry Smith    Notes:
841bd481603SBarry Smith    See the chapter in the users manual on performance for more details
842bd481603SBarry Smith 
843bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
844bd481603SBarry Smith 
845bd481603SBarry Smith   Concepts: preallocation^Matrix
846bd481603SBarry Smith 
847bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
848bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
849bd481603SBarry Smith M*/
850ef66eb69SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree(dnz);CHKERRQ(_4_ierr);}
8517c922b88SBarry Smith 
852bd481603SBarry Smith 
853bd481603SBarry Smith 
8547b80b807SBarry Smith /* Routines unique to particular data structures */
855be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
856ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
857be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBDiagGetData(Mat,PetscInt*,PetscInt*,PetscInt*[],PetscInt*[],PetscScalar***);
858698d4c6aSKris Buschelman 
859be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
860be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
861698d4c6aSKris Buschelman 
862be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
863be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
864be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
865c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
866c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
8677b80b807SBarry Smith 
868a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
869a96a251dSBarry Smith 
870be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
871ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
872be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
873ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
874be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
875ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
876be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDensePreallocation(Mat,PetscScalar[]);
877be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]);
878273d9f13SBarry Smith 
879be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
880ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
881be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
882be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
88353803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
884be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
885be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
886be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDensePreallocation(Mat,PetscScalar[]);
887be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]);
888be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
889be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
890284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
891be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
892be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
893be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
894be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
895273d9f13SBarry Smith 
896be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
897ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*);
8981b807ce4Svictorle 
899be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
900be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
9012e8a6d31SBarry Smith 
902be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
9033a7fca6bSBarry Smith 
9047b80b807SBarry Smith /*
9057b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
90694b7f48cSBarry Smith   done through the KSP and PC interfaces.
9077b80b807SBarry Smith */
9087b80b807SBarry Smith 
909d9274352SBarry Smith /*E
910d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
911d9274352SBarry Smith        with an optional dynamic library name, for example
912d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
913d9274352SBarry Smith 
914d9274352SBarry Smith    Level: beginner
915d9274352SBarry Smith 
916d9274352SBarry Smith .seealso: MatGetOrdering()
917d9274352SBarry Smith E*/
91849773a63SBarry Smith #define MatOrderingType char*
919b12f92e5SBarry Smith #define MATORDERING_NATURAL   "natural"
920b12f92e5SBarry Smith #define MATORDERING_ND        "nd"
921b12f92e5SBarry Smith #define MATORDERING_1WD       "1wd"
922b12f92e5SBarry Smith #define MATORDERING_RCM       "rcm"
923b12f92e5SBarry Smith #define MATORDERING_QMD       "qmd"
924b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength"
92562152c8bSBarry Smith #define MATORDERING_DSC_ND    "dsc_nd"
92662152c8bSBarry Smith #define MATORDERING_DSC_MMD   "dsc_mmd"
92762152c8bSBarry Smith #define MATORDERING_DSC_MDF   "dsc_mdf"
928c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
929c06d978dSMatthew Knepley #define MATORDERING_IDENTITY  "identity"
930c06d978dSMatthew Knepley #define MATORDERING_REVERSE   "reverse"
931b12f92e5SBarry Smith 
932be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
933be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
934be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
93530de9b25SBarry Smith 
93630de9b25SBarry Smith /*MC
93730de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
93830de9b25SBarry Smith                                matrix package.
93930de9b25SBarry Smith 
94030de9b25SBarry Smith    Synopsis:
941c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
94230de9b25SBarry Smith 
94330de9b25SBarry Smith    Not Collective
94430de9b25SBarry Smith 
94530de9b25SBarry Smith    Input Parameters:
94630de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
94730de9b25SBarry Smith .  path - location of library where creation routine is
94830de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
94930de9b25SBarry Smith -  function - function pointer that creates the ordering
95030de9b25SBarry Smith 
95130de9b25SBarry Smith    Level: developer
95230de9b25SBarry Smith 
95330de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
95430de9b25SBarry Smith    is ignored.
95530de9b25SBarry Smith 
95630de9b25SBarry Smith    Sample usage:
95730de9b25SBarry Smith .vb
95830de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
95930de9b25SBarry Smith                "MyOrder",MyOrder);
96030de9b25SBarry Smith .ve
96130de9b25SBarry Smith 
96230de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
96330de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
96430de9b25SBarry Smith    or at runtime via the option
9652401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
96630de9b25SBarry Smith 
967ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
96830de9b25SBarry Smith 
96930de9b25SBarry Smith .keywords: matrix, ordering, register
97030de9b25SBarry Smith 
97130de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
97230de9b25SBarry Smith M*/
973aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
974f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
975b12f92e5SBarry Smith #else
976f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
977b12f92e5SBarry Smith #endif
97830de9b25SBarry Smith 
979be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
980be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
9812bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
982b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
983d4fbbf0eSBarry Smith 
984be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
985a2ce50c7SBarry Smith 
986d91e6319SBarry Smith /*S
9872401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
988b00f7748SHong Zhang 
98915e8a5b3SHong Zhang    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE
990b00f7748SHong Zhang 
99115e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
992b00f7748SHong Zhang 
993b00f7748SHong Zhang    Level: developer
994b00f7748SHong Zhang 
995d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
996d7d82daaSBarry Smith           MatFactorInfoInitialize()
997b00f7748SHong Zhang 
998b00f7748SHong Zhang S*/
999b00f7748SHong Zhang typedef struct {
10000a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
1001fbf22428SSatish Balay   PetscReal     shiftpd;        /* if true, shift until positive pivots */
10022cea7109SBarry Smith   PetscReal     shift_fraction; /* record shift fraction taken */
100315e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
100415e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1005b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
100615e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
1007f6275e2eSBarry Smith   PetscReal     fill;           /* expected fill; nonzeros in factored matrix/nonzeros in original matrix*/
1008348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1009bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1010bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
101115e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
101215e8a5b3SHong Zhang } MatFactorInfo;
1013ffa6d0a5SLois Curfman McInnes 
1014be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
1015be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,MatFactorInfo*);
1016be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1017be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,MatFactorInfo*,Mat*);
1018be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,MatFactorInfo*);
1019be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,MatFactorInfo*);
1020be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1021be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1022be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1023be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,MatFactorInfo*);
1024be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,MatFactorInfo*,Mat*);
1025be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,MatFactorInfo*,Mat *);
1026be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1027be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1028be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1029be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1030be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1031be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1032be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1033be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
10348ed539a5SBarry Smith 
1035be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1036bb5a7306SBarry Smith 
1037d91e6319SBarry Smith /*E
1038d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1039bb1eb677SSatish Balay 
1040d91e6319SBarry Smith     Level: beginner
1041d91e6319SBarry Smith 
1042d9274352SBarry Smith    May be bitwise ORd together
1043d9274352SBarry Smith 
1044d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1045d91e6319SBarry Smith 
10464e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10474e7234bfSBarry Smith 
1048a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax()
1049d91e6319SBarry Smith E*/
1050ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1051ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1052ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
105384cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1054be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1055be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10568ed539a5SBarry Smith 
1057d4fbbf0eSBarry Smith /*
1058639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1059639f9d9dSBarry Smith */
1060b12f92e5SBarry Smith 
1061d9274352SBarry Smith /*E
1062d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1063d9274352SBarry Smith        with an optional dynamic library name, for example
1064d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1065d9274352SBarry Smith 
1066d9274352SBarry Smith    Level: beginner
1067d9274352SBarry Smith 
1068d9274352SBarry Smith .seealso: MatGetColoring()
1069d9274352SBarry Smith E*/
107049773a63SBarry Smith #define MatColoringType char*
1071b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1072b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1073b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1074b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1075b12f92e5SBarry Smith 
1076be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*);
1077be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatColoringType,ISColoring *));
107830de9b25SBarry Smith 
107930de9b25SBarry Smith /*MC
108030de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
108130de9b25SBarry Smith                                matrix package.
108230de9b25SBarry Smith 
108330de9b25SBarry Smith    Synopsis:
1084c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
108530de9b25SBarry Smith 
108630de9b25SBarry Smith    Not Collective
108730de9b25SBarry Smith 
108830de9b25SBarry Smith    Input Parameters:
108930de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
109030de9b25SBarry Smith .  path - location of library where creation routine is
109130de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
109230de9b25SBarry Smith -  function - function pointer that creates the coloring
109330de9b25SBarry Smith 
109430de9b25SBarry Smith    Level: developer
109530de9b25SBarry Smith 
109630de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
109730de9b25SBarry Smith    is ignored.
109830de9b25SBarry Smith 
109930de9b25SBarry Smith    Sample usage:
110030de9b25SBarry Smith .vb
110130de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
110230de9b25SBarry Smith                "MyColor",MyColor);
110330de9b25SBarry Smith .ve
110430de9b25SBarry Smith 
110530de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
110630de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
110730de9b25SBarry Smith    or at runtime via the option
110830de9b25SBarry Smith $     -mat_coloring_type my_color
110930de9b25SBarry Smith 
1110ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
111130de9b25SBarry Smith 
111230de9b25SBarry Smith .keywords: matrix, Coloring, register
111330de9b25SBarry Smith 
111430de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
111530de9b25SBarry Smith M*/
1116aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1117f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1118b12f92e5SBarry Smith #else
1119f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1120b12f92e5SBarry Smith #endif
112130de9b25SBarry Smith 
11222bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1123f1144a30SSatish Balay 
1124be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1125be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1126be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1127639f9d9dSBarry Smith 
1128d9274352SBarry Smith /*S
1129d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1130d9274352SBarry Smith         and coloring
1131639f9d9dSBarry Smith 
1132d9274352SBarry Smith    Level: beginner
1133d9274352SBarry Smith 
1134d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1135d9274352SBarry Smith 
1136d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1137d9274352SBarry Smith S*/
1138e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1139639f9d9dSBarry Smith 
1140be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1141be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1142be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1143be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
11444a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1145be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1146be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring,PetscInt);
1147be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring,PetscInt*);
1148be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1149be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1150be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1151be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring);
1152be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1153be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1154639f9d9dSBarry Smith /*
11550752156aSBarry Smith     These routines are for partitioning matrices: currently used only
11563eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
11570752156aSBarry Smith */
1158ca161407SBarry Smith 
1159d9274352SBarry Smith /*S
1160d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1161d9274352SBarry Smith 
1162d9274352SBarry Smith    Level: beginner
1163d9274352SBarry Smith 
1164d9274352SBarry Smith   Concepts: partitioning
1165d9274352SBarry Smith 
1166743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1167d9274352SBarry Smith S*/
116891e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1169d9274352SBarry Smith 
1170d9274352SBarry Smith /*E
11715bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1172d9274352SBarry Smith        with an optional dynamic library name, for example
1173d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1174d9274352SBarry Smith 
1175d9274352SBarry Smith    Level: beginner
1176d9274352SBarry Smith 
1177b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1178d9274352SBarry Smith E*/
117949773a63SBarry Smith #define MatPartitioningType char*
11808ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1181c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
11828ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
118371306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
118471306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
118571306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
118671306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
118771306c60Sjroman 
1188ca161407SBarry Smith 
1189be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1190be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1191be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1192be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1193be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1194be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1195be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1196be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
11972aabb6bbSBarry Smith 
1198be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
119930de9b25SBarry Smith 
120030de9b25SBarry Smith /*MC
120130de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
120230de9b25SBarry Smith    matrix package.
120330de9b25SBarry Smith 
120430de9b25SBarry Smith    Synopsis:
1205c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
120630de9b25SBarry Smith 
120730de9b25SBarry Smith    Not Collective
120830de9b25SBarry Smith 
120930de9b25SBarry Smith    Input Parameters:
121030de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
121130de9b25SBarry Smith .  path - location of library where creation routine is
121230de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
121330de9b25SBarry Smith -  function - function pointer that creates the partitioning type
121430de9b25SBarry Smith 
121530de9b25SBarry Smith    Level: developer
121630de9b25SBarry Smith 
121730de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
121830de9b25SBarry Smith    is ignored.
121930de9b25SBarry Smith 
122030de9b25SBarry Smith    Sample usage:
122130de9b25SBarry Smith .vb
122230de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
122330de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
122430de9b25SBarry Smith .ve
122530de9b25SBarry Smith 
122630de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
122730de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
122830de9b25SBarry Smith    or at runtime via the option
122930de9b25SBarry Smith $     -mat_partitioning_type my_part
123030de9b25SBarry Smith 
1231ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
123230de9b25SBarry Smith 
123330de9b25SBarry Smith .keywords: matrix, partitioning, register
123430de9b25SBarry Smith 
123530de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
123630de9b25SBarry Smith M*/
1237aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1238f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
12392aabb6bbSBarry Smith #else
1240f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
12412aabb6bbSBarry Smith #endif
12422aabb6bbSBarry Smith 
12432bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1244f1144a30SSatish Balay 
1245be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1246be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
12472bad1931SBarry Smith 
1248be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1249be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1250be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1251ca161407SBarry Smith 
1252be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
12530752156aSBarry Smith 
1254be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1255be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
125671306c60Sjroman 
1257954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
1258be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
125971306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1260be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1261be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
126271306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1263be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1264be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1265be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
126671306c60Sjroman 
126771306c60Sjroman #define MP_PARTY_OPT "opt"
126871306c60Sjroman #define MP_PARTY_LIN "lin"
126971306c60Sjroman #define MP_PARTY_SCA "sca"
127071306c60Sjroman #define MP_PARTY_RAN "ran"
127171306c60Sjroman #define MP_PARTY_GBF "gbf"
127271306c60Sjroman #define MP_PARTY_GCF "gcf"
127371306c60Sjroman #define MP_PARTY_BUB "bub"
127471306c60Sjroman #define MP_PARTY_DEF "def"
1275be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
127671306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
127771306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
127871306c60Sjroman #define MP_PARTY_NONE "no"
1279be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1280be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1281be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1282be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
128371306c60Sjroman 
128471306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1285be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1286be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1287be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1288be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1289be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
129071306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1291be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1292be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
129471306c60Sjroman 
12950752156aSBarry Smith /*
12960a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1297d4fbbf0eSBarry Smith */
12981c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
12991c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13001c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13011c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13021c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13037c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13047c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13051c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13061c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13077c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13087c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13091c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13101c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
13111c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
13121c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13131c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13141c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13151c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13161c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13171c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13181c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13191c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
13201c1c02c0SLois Curfman McInnes                MATOP_COMPRESS=22,
13211c1c02c0SLois Curfman McInnes                MATOP_SET_OPTION=23,
13221c1c02c0SLois Curfman McInnes                MATOP_ZERO_ENTRIES=24,
13231c1c02c0SLois Curfman McInnes                MATOP_ZERO_ROWS=25,
13241c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_SYMBOLIC=26,
13251c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_NUMERIC=27,
13261c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_SYMBOLIC=28,
13271c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_NUMERIC=29,
1328d643ce63SMatthew Knepley                MATOP_SETUP_PREALLOCATION=30,
1329d643ce63SMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=31,
1330d643ce63SMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=32,
1331d643ce63SMatthew Knepley                MATOP_GET_ARRAY=33,
1332d643ce63SMatthew Knepley                MATOP_RESTORE_ARRAY=34,
1333e26f1c3cSBarry Smith                MATOP_DUPLICATE=35,
1334d643ce63SMatthew Knepley                MATOP_FORWARD_SOLVE=36,
1335d643ce63SMatthew Knepley                MATOP_BACKWARD_SOLVE=37,
1336d643ce63SMatthew Knepley                MATOP_ILUFACTOR=38,
1337d643ce63SMatthew Knepley                MATOP_ICCFACTOR=39,
1338d643ce63SMatthew Knepley                MATOP_AXPY=40,
1339d643ce63SMatthew Knepley                MATOP_GET_SUBMATRICES=41,
1340d643ce63SMatthew Knepley                MATOP_INCREASE_OVERLAP=42,
1341d643ce63SMatthew Knepley                MATOP_GET_VALUES=43,
1342d643ce63SMatthew Knepley                MATOP_COPY=44,
1343d643ce63SMatthew Knepley                MATOP_PRINT_HELP=45,
1344d643ce63SMatthew Knepley                MATOP_SCALE=46,
1345d643ce63SMatthew Knepley                MATOP_SHIFT=47,
1346d643ce63SMatthew Knepley                MATOP_DIAGONAL_SHIFT=48,
1347d643ce63SMatthew Knepley                MATOP_ILUDT_FACTOR=49,
1348d643ce63SMatthew Knepley                MATOP_GET_BLOCK_SIZE=50,
1349d643ce63SMatthew Knepley                MATOP_GET_ROW_IJ=51,
1350d643ce63SMatthew Knepley                MATOP_RESTORE_ROW_IJ=52,
1351d643ce63SMatthew Knepley                MATOP_GET_COLUMN_IJ=53,
1352d643ce63SMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=54,
1353d643ce63SMatthew Knepley                MATOP_FDCOLORING_CREATE=55,
1354d643ce63SMatthew Knepley                MATOP_COLORING_PATCH=56,
1355d643ce63SMatthew Knepley                MATOP_SET_UNFACTORED=57,
1356d643ce63SMatthew Knepley                MATOP_PERMUTE=58,
1357d643ce63SMatthew Knepley                MATOP_SET_VALUES_BLOCKED=59,
1358d643ce63SMatthew Knepley                MATOP_GET_SUBMATRIX=60,
1359d643ce63SMatthew Knepley                MATOP_DESTROY=61,
1360d643ce63SMatthew Knepley                MATOP_VIEW=62,
1361d643ce63SMatthew Knepley                MATOP_GET_MAPS=63,
1362d643ce63SMatthew Knepley                MATOP_USE_SCALED_FORM=64,
1363d643ce63SMatthew Knepley                MATOP_SCALE_SYSTEM=65,
1364d643ce63SMatthew Knepley                MATOP_UNSCALE_SYSTEM=66,
1365ba16a8cbSKris Buschelman                MATOP_SET_LOCAL_TO_GLOBAL_MAP=67,
1366d643ce63SMatthew Knepley                MATOP_SET_VALUES_LOCAL=68,
1367d643ce63SMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=69,
1368d643ce63SMatthew Knepley                MATOP_GET_ROW_MAX=70,
1369d643ce63SMatthew Knepley                MATOP_CONVERT=71,
1370d643ce63SMatthew Knepley                MATOP_SET_COLORING=72,
1371d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1372d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1373d643ce63SMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1374d643ce63SMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1375ba16a8cbSKris Buschelman                MATOP_MULT_CON=77,
1376ba16a8cbSKris Buschelman                MATOP_MULT_TRANSPOSE_CON=78,
1377ba16a8cbSKris Buschelman                MATOP_ILU_FACTOR_SYMBOLIC_CON=79,
1378d643ce63SMatthew Knepley                MATOP_PERMUTE_SPARSIFY=80,
1379d643ce63SMatthew Knepley                MATOP_MULT_MULTIPLE=81,
138041acf15aSKris Buschelman                MATOP_SOLVE_MULTIPLE=82,
138141acf15aSKris Buschelman                MATOP_GET_INERTIA=83,
138241acf15aSKris Buschelman                MATOP_LOAD=84,
138341acf15aSKris Buschelman                MATOP_IS_SYMMETRIC=85,
138441acf15aSKris Buschelman                MATOP_IS_HERMITIAN=86,
138541acf15aSKris Buschelman                MATOP_IS_STRUCTURALLY_SYMMETRIC=87,
138641acf15aSKris Buschelman                MATOP_PB_RELAX=88,
138741acf15aSKris Buschelman                MATOP_GET_VECS=89,
138841acf15aSKris Buschelman                MATOP_MAT_MULT=90,
138941acf15aSKris Buschelman                MATOP_MAT_MULT_SYMBOLIC=91,
139041acf15aSKris Buschelman                MATOP_MAT_MULT_NUMERIC=92,
139141acf15aSKris Buschelman                MATOP_PTAP=93,
139241acf15aSKris Buschelman                MATOP_PTAP_SYMBOLIC=94,
1393bc011b1eSHong Zhang                MATOP_PTAP_NUMERIC=95,
1394bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE=96,
1395bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97,
13967ba1a0bfSKris Buschelman                MATOP_MAT_MULTTRANSPOSE_NUMERIC=98,
13977ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_SEQAIJ=99,
13987ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_SEQAIJ=100,
13997ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_MPIAIJ=101,
140087d4246cSBarry Smith                MATOP_PTAP_NUMERIC_MPIAIJ=102,
1401f5edf698SHong Zhang                MATOP_SET_VALUES_ROW = 105,
1402d5c63f83SSatish Balay                MATOP_GET_ROW_UTRIANGULAR=108,
14032bebee5dSHong Zhang                MATOP_RESTORE_ROW_UTRIANGULAR=109,
14042bebee5dSHong Zhang                MATOP_MATSOLVE=110
1405fae171e0SBarry Smith              } MatOperation;
1406be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1407be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1408be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1409be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1410112a2221SBarry Smith 
141190ace30eSBarry Smith /*
141290ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
141390ace30eSBarry Smith  stored in a universal format. By changing the format with
1414fb9695e5SSatish Balay  PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will
141590ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
141690ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
141790ace30eSBarry Smith  read into matrices of the same time.
141890ace30eSBarry Smith */
141990ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
142090ace30eSBarry Smith 
1421be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1422be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
14231f4e1ec7SBarry Smith 
1424d9274352SBarry Smith /*S
1425d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1426d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1427d9274352SBarry Smith 
1428f7a9e4ceSBarry Smith    Level: advanced
1429d9274352SBarry Smith 
1430d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1431d9274352SBarry Smith 
14326e1639daSBarry Smith   Users manual sections:
14336e1639daSBarry Smith .   sec_singular
14346e1639daSBarry Smith 
1435d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1436d9274352SBarry Smith S*/
143774637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1438d9274352SBarry Smith 
1439be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1440281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*);
1441be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1442be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1443be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
1444be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat);
144574637425SBarry Smith 
1446be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1447be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
144946129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
14503f1d51d7SBarry Smith 
1451be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1452be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1453be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1454c4f061fbSSatish Balay 
1455be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1456b0a32e0cSBarry Smith 
1457be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
145804f1ad80SBarry Smith 
1459be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1460be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
14617dbadf16SMatthew Knepley 
1462e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
14632eac72dbSBarry Smith #endif
1464