xref: /petsc/include/petscmat.h (revision 814655b831001cc2c5390a0269bf1b7bc476caca)
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"
55b3a1e11cSKris Buschelman #define MATSEQAIJSPOOLES   "seqaijspooles"
56d10c748bSKris Buschelman #define MATMPIAIJSPOOLES   "mpiaijspooles"
579abb65ffSKris Buschelman #define MATSEQSBAIJSPOOLES "seqsbaijspooles"
5822191285SKris Buschelman #define MATMPISBAIJSPOOLES "mpisbaijspooles"
59bb4d4166SHong Zhang #define MATAIJSPOOLES      "aijspooles"
60bb4d4166SHong Zhang #define MATSBAIJSPOOLES    "sbaijspooles"
6114b396bbSKris Buschelman #define MATSUPERLU         "superlu"
621d96aa28SKris Buschelman #define MATSUPERLU_DIST    "superlu_dist"
631677d0b8SKris Buschelman #define MATUMFPACK         "umfpack"
64e8aa55a4SKris Buschelman #define MATESSL            "essl"
654eb8e494SKris Buschelman #define MATLUSOL           "lusol"
66397b6df1SKris Buschelman #define MATAIJMUMPS        "aijmumps"
67397b6df1SKris Buschelman #define MATSBAIJMUMPS      "sbaijmumps"
688da957c5SKris Buschelman #define MATDSCPACK         "dscpack"
697065e2aaSKris Buschelman #define MATMATLAB          "matlab"
704b8d99adSRichard Tran Mills #define MATSEQCSRPERM      "seqcsrperm"
7118def556SRichard Tran Mills #define MATMPICSRPERM      "mpicsrperm"
72*814655b8SBarry Smith #define MATCSRPERM         "csrperm"
73f69a0ea3SMatthew Knepley #define MatType const char*
74d91e6319SBarry Smith 
75c06d978dSMatthew Knepley /* Logging support */
76552e946dSBarry Smith #define    MAT_FILE_COOKIE 1211216    /* used to indicate matrices in binary files */
77be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE;
78be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
79be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE;
80be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE;
816849ba73SBarry Smith extern PetscEvent  MAT_Mult, MAT_MultMatrixFree, MAT_Mults, MAT_MultConstrained, MAT_MultAdd, MAT_MultTranspose;
826849ba73SBarry Smith extern PetscEvent  MAT_MultTransposeConstrained, MAT_MultTransposeAdd, MAT_Solve, MAT_Solves, MAT_SolveAdd, MAT_SolveTranspose;
836849ba73SBarry Smith extern PetscEvent  MAT_SolveTransposeAdd, MAT_Relax, MAT_ForwardSolve, MAT_BackwardSolve, MAT_LUFactor, MAT_LUFactorSymbolic;
846849ba73SBarry Smith extern PetscEvent  MAT_LUFactorNumeric, MAT_CholeskyFactor, MAT_CholeskyFactorSymbolic, MAT_CholeskyFactorNumeric, MAT_ILUFactor;
856849ba73SBarry Smith extern PetscEvent  MAT_ILUFactorSymbolic, MAT_ICCFactorSymbolic, MAT_Copy, MAT_Convert, MAT_Scale, MAT_AssemblyBegin;
866849ba73SBarry Smith extern PetscEvent  MAT_AssemblyEnd, MAT_SetValues, MAT_GetValues, MAT_GetRow, MAT_GetSubMatrices, MAT_GetColoring, MAT_GetOrdering;
876849ba73SBarry Smith extern PetscEvent  MAT_IncreaseOverlap, MAT_Partitioning, MAT_ZeroEntries, MAT_Load, MAT_View, MAT_AXPY, MAT_FDColoringCreate;
886849ba73SBarry Smith extern PetscEvent  MAT_FDColoringApply, MAT_Transpose, MAT_FDColoringFunction;
89cb00f407SKris Buschelman extern PetscEvent  MAT_MatMult, MAT_MatMultSymbolic, MAT_MatMultNumeric;
90534c1384SKris Buschelman extern PetscEvent  MAT_PtAP, MAT_PtAPSymbolic, MAT_PtAPNumeric;
91bc011b1eSHong Zhang extern PetscEvent  MAT_MatMultTranspose, MAT_MatMultTransposeSymbolic, MAT_MatMultTransposeNumeric;
92c06d978dSMatthew Knepley 
93ceb03754SKris Buschelman /*E
94ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
95ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
96ceb03754SKris Buschelman 
97ceb03754SKris Buschelman     Level: beginner
98ceb03754SKris Buschelman 
99ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
100ceb03754SKris Buschelman 
1010c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
102ceb03754SKris Buschelman E*/
103ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse;
104ceb03754SKris Buschelman 
105be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(char *);
106c06d978dSMatthew Knepley 
107f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
108ba966639SSatish Balay PetscPolymorphicFunction(MatCreate,(MPI_Comm comm,void *ctx),(comm,&A),Mat,A)
109ba966639SSatish Balay PetscPolymorphicFunction(MatCreate,(void *ctx),(PETSC_COMM_WORLD,&A),Mat,A)
110f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
111f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,MatType);
112be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
113be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
114be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
115be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
11630de9b25SBarry Smith 
117f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
118f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
119f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
120f69a0ea3SMatthew Knepley 
12130de9b25SBarry Smith /*MC
12230de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
12330de9b25SBarry Smith 
12430de9b25SBarry Smith    Synopsis:
125c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
12630de9b25SBarry Smith 
12730de9b25SBarry Smith    Not Collective
12830de9b25SBarry Smith 
12930de9b25SBarry Smith    Input Parameters:
13030de9b25SBarry Smith +  name - name of a new user-defined matrix type
13130de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
13230de9b25SBarry Smith .  name_create - name of routine to create method context
13330de9b25SBarry Smith -  routine_create - routine to create method context
13430de9b25SBarry Smith 
13530de9b25SBarry Smith    Notes:
13630de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
13730de9b25SBarry Smith 
13830de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
13930de9b25SBarry Smith    is ignored.
14030de9b25SBarry Smith 
14130de9b25SBarry Smith    Sample usage:
14230de9b25SBarry Smith .vb
14330de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
14430de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
14530de9b25SBarry Smith .ve
14630de9b25SBarry Smith 
14730de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
14830de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
14930de9b25SBarry Smith    or at runtime via the option
15030de9b25SBarry Smith $     -mat_type my_mat
15130de9b25SBarry Smith 
15230de9b25SBarry Smith    Level: advanced
15330de9b25SBarry Smith 
154ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
15530de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
15630de9b25SBarry Smith 
15730de9b25SBarry Smith .keywords: Mat, register
15830de9b25SBarry Smith 
15930de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
16030de9b25SBarry Smith 
16130de9b25SBarry Smith M*/
162273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
163273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
164273d9f13SBarry Smith #else
165273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
16630de9b25SBarry Smith #endif
16730de9b25SBarry Smith 
168273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
169b0a32e0cSBarry Smith extern PetscFList MatList;
17028988994SBarry Smith 
171be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
172be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
173be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
174ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
175ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
176ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
177ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
178ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
179ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
180ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
181be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
182ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(comm,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
183ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(comm,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
184ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(comm,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
185ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(comm,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
186ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(comm,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
187ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(comm,m,n,M,N,0,nnz,0,onz,A))
188ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(comm,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
189ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(PETSC_COMM_WORLD,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
190ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(PETSC_COMM_WORLD,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
191ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(PETSC_COMM_WORLD,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
192ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N),(PETSC_COMM_WORLD,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
193ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(PETSC_COMM_WORLD,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
194ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(PETSC_COMM_WORLD,m,n,M,N,0,nnz,0,onz,A))
195ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(PETSC_COMM_WORLD,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
196be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
197be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*);
198be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*);
199be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
200ba966639SSatish 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)
201ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
202ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
203ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
204ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
205ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
206ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
207be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
208ba966639SSatish 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)
209ba966639SSatish 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)
210ba966639SSatish 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)
211ba966639SSatish 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)
212ba966639SSatish 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))
213ba966639SSatish 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))
214ba966639SSatish 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))
215ba966639SSatish 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)
216ba966639SSatish 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)
217ba966639SSatish 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)
218ba966639SSatish 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)
219ba966639SSatish 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))
220ba966639SSatish 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))
221ba966639SSatish 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))
222be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
223be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
224ba966639SSatish 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)
225ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
226ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
227ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
228ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
229ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
230ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
231be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
232ba966639SSatish 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)
233ba966639SSatish 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)
234ba966639SSatish 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)
235ba966639SSatish 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)
236ba966639SSatish 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))
237ba966639SSatish 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))
238ba966639SSatish 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))
239ba966639SSatish 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)
240ba966639SSatish 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)
241ba966639SSatish 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)
242ba966639SSatish 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)
243ba966639SSatish 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))
244ba966639SSatish 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))
245ba966639SSatish 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))
246be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
247ba966639SSatish 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)
248ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
249be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
250be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
251ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
252f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
253be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
25421c89e3eSBarry Smith 
255be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPrintHelp(Mat);
256be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetPetscMaps(Mat,PetscMap*,PetscMap*);
257ec0117caSBarry Smith 
2580c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
2598ed539a5SBarry Smith /* ------------------------------------------------------------*/
260be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
261be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
26284cb2905SBarry Smith 
2632ef4de8bSBarry Smith /*S
2642ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
2652ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
2662ef4de8bSBarry Smith 
2672ef4de8bSBarry Smith    Level: beginner
2682ef4de8bSBarry Smith 
2692ef4de8bSBarry Smith   Concepts: matrix; linear operator
2702ef4de8bSBarry Smith 
271d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
2722ef4de8bSBarry Smith S*/
273435da068SBarry Smith typedef struct {
274c1ac3661SBarry Smith   PetscInt k,j,i,c;
275435da068SBarry Smith } MatStencil;
2762ef4de8bSBarry Smith 
277be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
278be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
279be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
280435da068SBarry Smith 
281be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
282be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
283be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
2843a7fca6bSBarry Smith 
285d91e6319SBarry Smith /*E
286d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
287d91e6319SBarry Smith      to continue to add values to it
288d91e6319SBarry Smith 
289d91e6319SBarry Smith     Level: beginner
290d91e6319SBarry Smith 
291d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
292d91e6319SBarry Smith E*/
2936d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
294be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
295be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
2974f9c727eSBarry Smith 
298be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt    MatSetValue_Row;
299be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt    MatSetValue_Column;
300be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscScalar MatSetValue_Value;
3011d73ed98SBarry Smith 
30230de9b25SBarry Smith /*MC
30330de9b25SBarry Smith    MatSetValue - Set a single entry into a matrix.
30430de9b25SBarry Smith 
30530de9b25SBarry Smith    Synopsis:
306c1ac3661SBarry Smith    PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode);
30730de9b25SBarry Smith 
30830de9b25SBarry Smith    Not collective
30930de9b25SBarry Smith 
31030de9b25SBarry Smith    Input Parameters:
31130de9b25SBarry Smith +  m - the matrix
31230de9b25SBarry Smith .  row - the row location of the entry
31330de9b25SBarry Smith .  col - the column location of the entry
31430de9b25SBarry Smith .  value - the value to insert
31530de9b25SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
31630de9b25SBarry Smith 
31730de9b25SBarry Smith    Notes:
31830de9b25SBarry Smith    For efficiency one should use MatSetValues() and set several or many
31930de9b25SBarry Smith    values simultaneously if possible.
32030de9b25SBarry Smith 
32130de9b25SBarry Smith    Level: beginner
32230de9b25SBarry Smith 
32330de9b25SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
32430de9b25SBarry Smith M*/
325b951964fSBarry Smith #define MatSetValue(v,i,j,va,mode) \
326f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \
3271d73ed98SBarry Smith    MatSetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode))
32830de9b25SBarry Smith 
329ea06a074SBarry Smith #define MatGetValue(v,i,j,va) \
330f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,0) || \
3311d73ed98SBarry Smith    MatGetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&va))
33230de9b25SBarry Smith 
333d91e6319SBarry Smith #define MatSetValueLocal(v,i,j,va,mode) \
334f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \
3351d73ed98SBarry Smith    MatSetValuesLocal(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode))
33630de9b25SBarry Smith 
337d91e6319SBarry Smith /*E
338d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
339d91e6319SBarry Smith 
340d91e6319SBarry Smith     Level: beginner
341d91e6319SBarry Smith 
3420a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
343d91e6319SBarry Smith 
344d91e6319SBarry Smith .seealso: MatSetOption()
345d91e6319SBarry Smith E*/
3466d4a8577SBarry Smith typedef enum {MAT_ROW_ORIENTED=1,MAT_COLUMN_ORIENTED=2,MAT_ROWS_SORTED=4,
3476d4a8577SBarry Smith               MAT_COLUMNS_SORTED=8,MAT_NO_NEW_NONZERO_LOCATIONS=16,
3486d4a8577SBarry Smith               MAT_YES_NEW_NONZERO_LOCATIONS=32,MAT_SYMMETRIC=64,
3496ca9ecd3SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC=65,MAT_NO_NEW_DIAGONALS=66,
3506ca9ecd3SBarry Smith               MAT_YES_NEW_DIAGONALS=67,MAT_INODE_LIMIT_1=68,MAT_INODE_LIMIT_2=69,
3516ca9ecd3SBarry Smith               MAT_INODE_LIMIT_3=70,MAT_INODE_LIMIT_4=71,MAT_INODE_LIMIT_5=72,
3526ca9ecd3SBarry Smith               MAT_IGNORE_OFF_PROC_ENTRIES=73,MAT_ROWS_UNSORTED=74,
3534787f768SSatish Balay               MAT_COLUMNS_UNSORTED=75,MAT_NEW_NONZERO_LOCATION_ERR=76,
3547c922b88SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR=77,MAT_USE_HASH_TABLE=78,
3552bad1931SBarry Smith               MAT_KEEP_ZEROED_ROWS=79,MAT_IGNORE_ZERO_ENTRIES=80,MAT_USE_INODES=81,
356efcf0fc3SBarry Smith               MAT_DO_NOT_USE_INODES=82,MAT_NOT_SYMMETRIC=83,MAT_HERMITIAN=84,
3579a4540c5SBarry Smith               MAT_NOT_STRUCTURALLY_SYMMETRIC=85,MAT_NOT_HERMITIAN=86,
358d487561eSHong Zhang               MAT_SYMMETRY_ETERNAL=87,MAT_NOT_SYMMETRY_ETERNAL=88,
359941593c8SHong Zhang               MAT_USE_COMPRESSEDROW=89,MAT_DO_NOT_USE_COMPRESSEDROW=90,
360941593c8SHong Zhang               MAT_IGNORE_LOWER_TRIANGULAR=91,MAT_ERROR_LOWER_TRIANGULAR=92} MatOption;
361be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption);
362be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,MatType*);
363ba966639SSatish Balay PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),MatType,t)
36484cb2905SBarry Smith 
365be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
366be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
367be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
368be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
370be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
372ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
373be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
374be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
375ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
376be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
3777b80b807SBarry Smith 
378be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
379be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
380ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
382be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
383ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
384ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0.0,&t),PetscTruth,t)
385be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
386ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
387be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
388be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
3892eac72dbSBarry Smith 
390d91e6319SBarry Smith /*E
391d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
392d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
393d91e6319SBarry Smith 
394d91e6319SBarry Smith     Level: beginner
395d91e6319SBarry Smith 
396d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
397d91e6319SBarry Smith 
398d91e6319SBarry Smith .seealso: MatDuplicate()
399d91e6319SBarry Smith E*/
4002e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
4012e8a6d31SBarry Smith 
402be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvertRegister(const char[],const char[],const char[],PetscErrorCode (*)(Mat,MatType,MatReuse,Mat*));
403273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
404273d9f13SBarry Smith #define MatConvertRegisterDynamic(a,b,c,d) MatConvertRegister(a,b,c,0)
405273d9f13SBarry Smith #else
406273d9f13SBarry Smith #define MatConvertRegisterDynamic(a,b,c,d) MatConvertRegister(a,b,c,d)
407273d9f13SBarry Smith #endif
408be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvertRegisterAll(const char[]);
409be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvertRegisterDestroy(void);
410273d9f13SBarry Smith extern PetscTruth MatConvertRegisterAllCalled;
411b0a32e0cSBarry Smith extern PetscFList MatConvertList;
412f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,MatType,MatReuse,Mat*);
413ba966639SSatish Balay PetscPolymorphicFunction(MatConvert,(Mat A,MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
414be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
415ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
416ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
41794a9d846SBarry Smith 
418d91e6319SBarry Smith /*E
419d91e6319SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
420d91e6319SBarry Smith 
421d91e6319SBarry Smith     Level: beginner
422d91e6319SBarry Smith 
423d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
424d91e6319SBarry Smith 
42594b7f48cSBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
426d91e6319SBarry Smith E*/
427c537a176SHong Zhang typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
428cb5b572fSBarry Smith 
429be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
430be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
431be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
432ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
433ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0.0,&t),PetscTruth,t)
434be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
435ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
436be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscTruth*);
437ba966639SSatish Balay PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,&t),PetscTruth,t)
438be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
439be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
440f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,MatType,Mat*);
441ba966639SSatish Balay PetscPolymorphicFunction(MatLoad,(PetscViewer v,MatType t),(v,t,&a),Mat,a)
4427b80b807SBarry Smith 
443be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
444be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
445be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
446be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
447d4fbbf0eSBarry Smith 
448d91e6319SBarry Smith /*S
449d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
450d91e6319SBarry Smith 
451d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
452d91e6319SBarry Smith 
453d91e6319SBarry Smith    Level: intermediate
454d91e6319SBarry Smith 
455d91e6319SBarry Smith   Concepts: matrix^nonzero information
456d91e6319SBarry Smith 
457d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
458d91e6319SBarry Smith S*/
4594e220ebcSLois Curfman McInnes typedef struct {
460b0a32e0cSBarry Smith   PetscLogDouble rows_global,columns_global;         /* number of global rows and columns */
461b0a32e0cSBarry Smith   PetscLogDouble rows_local,columns_local;           /* number of local rows and columns */
462b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
463b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
464b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
465b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
466b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
467b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
468b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4694e220ebcSLois Curfman McInnes } MatInfo;
4704e220ebcSLois Curfman McInnes 
471d9274352SBarry Smith /*E
472d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
473d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
474d9274352SBarry Smith 
475d9274352SBarry Smith     Level: beginner
476d9274352SBarry Smith 
477d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
478d9274352SBarry Smith 
479d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
480d9274352SBarry Smith E*/
4817b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
482be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
483be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
484be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
485be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec);
486be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,Mat*);
487ba966639SSatish Balay PetscPolymorphicFunction(MatTranspose,(Mat A),(A,&t),Mat,t)
488be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
489ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
490be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
491be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
492be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
493be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
494ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
495be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
496be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
497be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
4997b80b807SBarry Smith 
500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
501ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
502be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
503f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
504f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
505f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
506f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5077b80b807SBarry Smith 
508be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
509be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5115ef9f2a5SBarry Smith 
512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
513be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5157b80b807SBarry Smith 
516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
517be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *);
519be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
520be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
521be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
522be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
523be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
524be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
525be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
526be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscScalar**,Mat*);
5288efafbd8SBarry Smith 
529be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5307b80b807SBarry Smith 
531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
53422440eb1SKris Buschelman 
535be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
536be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
53822440eb1SKris Buschelman 
539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
542bc011b1eSHong Zhang 
543f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
544f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat);
545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat);
5461c741599SHong Zhang 
547f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
548f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5497b80b807SBarry Smith 
550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
552f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
553f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
556052efed2SBarry Smith 
557be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
558be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
55990f02eecSBarry Smith 
560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
5610c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
562ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
563be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
565bd481603SBarry Smith 
566bd481603SBarry Smith /*MC
567bd481603SBarry Smith    MatPreallocInitialize - Begins the block of code that will count the number of nonzeros per
568bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
569bd481603SBarry Smith 
570bd481603SBarry Smith    Synopsis:
571c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
572bd481603SBarry Smith 
573bd481603SBarry Smith    Collective on MPI_Comm
574bd481603SBarry Smith 
575bd481603SBarry Smith    Input Parameters:
576bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
577859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
578859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
579bd481603SBarry Smith 
580bd481603SBarry Smith    Output Parameters:
581bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
582bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
583bd481603SBarry Smith 
584bd481603SBarry Smith 
585bd481603SBarry Smith    Level: intermediate
586bd481603SBarry Smith 
587bd481603SBarry Smith    Notes:
588bd481603SBarry Smith    See the chapter in the users manual on performance for more details
589bd481603SBarry Smith 
5901d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
591bd481603SBarry Smith 
592bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
593bd481603SBarry Smith 
594bd481603SBarry Smith   Concepts: preallocation^Matrix
595bd481603SBarry Smith 
596bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
597bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
598bd481603SBarry Smith M*/
599c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6007c922b88SBarry Smith { \
601c1ac3661SBarry Smith   PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
602c1ac3661SBarry Smith   _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\
603c1ac3661SBarry Smith   _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\
604c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
605c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp;
6067c922b88SBarry Smith 
607bd481603SBarry Smith /*MC
608bd481603SBarry Smith    MatPreallocSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
609bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
610bd481603SBarry Smith 
611bd481603SBarry Smith    Synopsis:
612c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
613bd481603SBarry Smith 
614bd481603SBarry Smith    Collective on MPI_Comm
615bd481603SBarry Smith 
616bd481603SBarry Smith    Input Parameters:
617bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
618859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
619859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
620bd481603SBarry Smith 
621bd481603SBarry Smith    Output Parameters:
622bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
623bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
624bd481603SBarry Smith 
625bd481603SBarry Smith 
626bd481603SBarry Smith    Level: intermediate
627bd481603SBarry Smith 
628bd481603SBarry Smith    Notes:
629bd481603SBarry Smith    See the chapter in the users manual on performance for more details
630bd481603SBarry Smith 
6311d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
632bd481603SBarry Smith 
633bd481603SBarry Smith   Concepts: preallocation^Matrix
634bd481603SBarry Smith 
635bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
636bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
637bd481603SBarry Smith M*/
638222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
639222b16d4SBarry Smith { \
640c1ac3661SBarry Smith   PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__end; \
641c1ac3661SBarry Smith   _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\
642c1ac3661SBarry Smith   _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\
643c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
644c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp;
645222b16d4SBarry Smith 
646bd481603SBarry Smith /*MC
647bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
648bd481603SBarry Smith        inserted using a local number of the rows and columns
649bd481603SBarry Smith 
650bd481603SBarry Smith    Synopsis:
651c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
652bd481603SBarry Smith 
653bd481603SBarry Smith    Not Collective
654bd481603SBarry Smith 
655bd481603SBarry Smith    Input Parameters:
656bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
657bd481603SBarry Smith .  nrows - the number of rows indicated
6581d73ed98SBarry Smith .  rows - the indices of the rows
659bd481603SBarry Smith .  ncols - the number of columns in the matrix
660bd481603SBarry Smith .  cols - the columns indicated
661bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
662bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
663bd481603SBarry Smith 
664bd481603SBarry Smith 
665bd481603SBarry Smith    Level: intermediate
666bd481603SBarry Smith 
667bd481603SBarry Smith    Notes:
668bd481603SBarry Smith    See the chapter in the users manual on performance for more details
669bd481603SBarry Smith 
6701d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
671bd481603SBarry Smith 
672bd481603SBarry Smith   Concepts: preallocation^Matrix
673bd481603SBarry Smith 
674bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
675bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
676bd481603SBarry Smith M*/
677c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
678c4f061fbSSatish Balay {\
679c1ac3661SBarry Smith   PetscInt __l;\
680ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
681ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
682c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
683ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
684c4f061fbSSatish Balay   }\
685c4f061fbSSatish Balay }
686c4f061fbSSatish Balay 
687bd481603SBarry Smith /*MC
688bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
689bd481603SBarry Smith        inserted using a local number of the rows and columns
690bd481603SBarry Smith 
691bd481603SBarry Smith    Synopsis:
692c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
693bd481603SBarry Smith 
694bd481603SBarry Smith    Not Collective
695bd481603SBarry Smith 
696bd481603SBarry Smith    Input Parameters:
697bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
698bd481603SBarry Smith .  nrows - the number of rows indicated
6991d73ed98SBarry Smith .  rows - the indices of the rows
700bd481603SBarry Smith .  ncols - the number of columns in the matrix
701bd481603SBarry Smith .  cols - the columns indicated
702bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
703bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
704bd481603SBarry Smith 
705bd481603SBarry Smith 
706bd481603SBarry Smith    Level: intermediate
707bd481603SBarry Smith 
708bd481603SBarry Smith    Notes:
709bd481603SBarry Smith    See the chapter in the users manual on performance for more details
710bd481603SBarry Smith 
711bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
712bd481603SBarry Smith 
713bd481603SBarry Smith   Concepts: preallocation^Matrix
714bd481603SBarry Smith 
715bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
716bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
717bd481603SBarry Smith M*/
718d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
719d3d32019SBarry Smith {\
720c1ac3661SBarry Smith   PetscInt __l;\
721d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
722d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
723d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
724d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
725d3d32019SBarry Smith   }\
726d3d32019SBarry Smith }
727d3d32019SBarry Smith 
728bd481603SBarry Smith /*MC
729bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
730bd481603SBarry Smith        inserted using a local number of the rows and columns
731bd481603SBarry Smith 
732bd481603SBarry Smith    Synopsis:
733c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
734bd481603SBarry Smith 
735bd481603SBarry Smith    Not Collective
736bd481603SBarry Smith 
737bd481603SBarry Smith    Input Parameters:
738bd481603SBarry Smith +  nrows - the number of rows indicated
7391d73ed98SBarry Smith .  rows - the indices of the rows
740bd481603SBarry Smith .  ncols - the number of columns in the matrix
741a50f8bf6SHong Zhang -  cols - the columns indicated
742a50f8bf6SHong Zhang 
743a50f8bf6SHong Zhang    Output Parameters:
744a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
745bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
746bd481603SBarry Smith 
747bd481603SBarry Smith 
748bd481603SBarry Smith    Level: intermediate
749bd481603SBarry Smith 
750bd481603SBarry Smith    Notes:
751bd481603SBarry Smith    See the chapter in the users manual on performance for more details
752bd481603SBarry Smith 
753bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
754bd481603SBarry Smith 
755bd481603SBarry Smith   Concepts: preallocation^Matrix
756bd481603SBarry Smith 
757bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
758bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
759bd481603SBarry Smith M*/
760c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
761c1ac3661SBarry Smith { PetscInt __i; \
7627c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
7637c922b88SBarry Smith     if (cols[__i] < __start || cols[__i] >= __end) onz[row - __rstart]++; \
7647c922b88SBarry Smith   }\
7657c922b88SBarry Smith   dnz[row - __rstart] = nc - onz[row - __rstart];\
7667c922b88SBarry Smith }
7677c922b88SBarry Smith 
768bd481603SBarry Smith /*MC
769bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
770bd481603SBarry Smith        inserted using a local number of the rows and columns
771bd481603SBarry Smith 
772bd481603SBarry Smith    Synopsis:
773c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
774bd481603SBarry Smith 
775bd481603SBarry Smith    Not Collective
776bd481603SBarry Smith 
777bd481603SBarry Smith    Input Parameters:
778bd481603SBarry Smith +  nrows - the number of rows indicated
7791d73ed98SBarry Smith .  rows - the indices of the rows
780bd481603SBarry Smith .  ncols - the number of columns in the matrix
781bd481603SBarry Smith .  cols - the columns indicated
782bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
783bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
784bd481603SBarry Smith 
785bd481603SBarry Smith 
786bd481603SBarry Smith    Level: intermediate
787bd481603SBarry Smith 
788bd481603SBarry Smith    Notes:
789bd481603SBarry Smith    See the chapter in the users manual on performance for more details
790bd481603SBarry Smith 
791bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
792bd481603SBarry Smith 
793bd481603SBarry Smith   Concepts: preallocation^Matrix
794bd481603SBarry Smith 
795bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
796bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
797bd481603SBarry Smith M*/
798d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
799c1ac3661SBarry Smith { PetscInt __i; \
800d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
801d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
802d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
803d3d32019SBarry Smith   }\
804d3d32019SBarry Smith }
805d3d32019SBarry Smith 
806bd481603SBarry Smith /*MC
807bd481603SBarry Smith    MatPreallocFinalize - Ends the block of code that will count the number of nonzeros per
808bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
809bd481603SBarry Smith 
810bd481603SBarry Smith    Synopsis:
811c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
812bd481603SBarry Smith 
813bd481603SBarry Smith    Collective on MPI_Comm
814bd481603SBarry Smith 
815bd481603SBarry Smith    Input Parameters:
816bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
817bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
818bd481603SBarry Smith 
819bd481603SBarry Smith 
820bd481603SBarry Smith    Level: intermediate
821bd481603SBarry Smith 
822bd481603SBarry Smith    Notes:
823bd481603SBarry Smith    See the chapter in the users manual on performance for more details
824bd481603SBarry Smith 
825bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
826bd481603SBarry Smith 
827bd481603SBarry Smith   Concepts: preallocation^Matrix
828bd481603SBarry Smith 
829bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
830bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
831bd481603SBarry Smith M*/
832ef66eb69SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree(dnz);CHKERRQ(_4_ierr);}
8337c922b88SBarry Smith 
834bd481603SBarry Smith 
835bd481603SBarry Smith 
8367b80b807SBarry Smith /* Routines unique to particular data structures */
837be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
838ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
839be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBDiagGetData(Mat,PetscInt*,PetscInt*,PetscInt*[],PetscInt*[],PetscScalar***);
840698d4c6aSKris Buschelman 
841be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
842be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
843698d4c6aSKris Buschelman 
844be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
845be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
846be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
8477b80b807SBarry Smith 
848a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
849a96a251dSBarry Smith 
850be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
851ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
852be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
853ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
854be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
855ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
856be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDensePreallocation(Mat,PetscScalar[]);
857be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]);
858be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
859273d9f13SBarry Smith 
860be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
861ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
862be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
863be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
864be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
865be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
866be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDensePreallocation(Mat,PetscScalar[]);
867be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]);
868be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
869be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
870be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
871be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
872be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
873be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
874273d9f13SBarry Smith 
875be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
8761b807ce4Svictorle 
877be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
878be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
8792e8a6d31SBarry Smith 
880be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
8813a7fca6bSBarry Smith 
8827b80b807SBarry Smith /*
8837b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
88494b7f48cSBarry Smith   done through the KSP and PC interfaces.
8857b80b807SBarry Smith */
8867b80b807SBarry Smith 
887d9274352SBarry Smith /*E
888d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
889d9274352SBarry Smith        with an optional dynamic library name, for example
890d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
891d9274352SBarry Smith 
892d9274352SBarry Smith    Level: beginner
893d9274352SBarry Smith 
894d9274352SBarry Smith .seealso: MatGetOrdering()
895d9274352SBarry Smith E*/
89649773a63SBarry Smith #define MatOrderingType char*
897b12f92e5SBarry Smith #define MATORDERING_NATURAL   "natural"
898b12f92e5SBarry Smith #define MATORDERING_ND        "nd"
899b12f92e5SBarry Smith #define MATORDERING_1WD       "1wd"
900b12f92e5SBarry Smith #define MATORDERING_RCM       "rcm"
901b12f92e5SBarry Smith #define MATORDERING_QMD       "qmd"
902b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength"
90362152c8bSBarry Smith #define MATORDERING_DSC_ND    "dsc_nd"
90462152c8bSBarry Smith #define MATORDERING_DSC_MMD   "dsc_mmd"
90562152c8bSBarry Smith #define MATORDERING_DSC_MDF   "dsc_mdf"
906c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
907c06d978dSMatthew Knepley #define MATORDERING_IDENTITY  "identity"
908c06d978dSMatthew Knepley #define MATORDERING_REVERSE   "reverse"
909b12f92e5SBarry Smith 
910be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
911be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
912be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
91330de9b25SBarry Smith 
91430de9b25SBarry Smith /*MC
91530de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
91630de9b25SBarry Smith                                matrix package.
91730de9b25SBarry Smith 
91830de9b25SBarry Smith    Synopsis:
919c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
92030de9b25SBarry Smith 
92130de9b25SBarry Smith    Not Collective
92230de9b25SBarry Smith 
92330de9b25SBarry Smith    Input Parameters:
92430de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
92530de9b25SBarry Smith .  path - location of library where creation routine is
92630de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
92730de9b25SBarry Smith -  function - function pointer that creates the ordering
92830de9b25SBarry Smith 
92930de9b25SBarry Smith    Level: developer
93030de9b25SBarry Smith 
93130de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
93230de9b25SBarry Smith    is ignored.
93330de9b25SBarry Smith 
93430de9b25SBarry Smith    Sample usage:
93530de9b25SBarry Smith .vb
93630de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
93730de9b25SBarry Smith                "MyOrder",MyOrder);
93830de9b25SBarry Smith .ve
93930de9b25SBarry Smith 
94030de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
94130de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
94230de9b25SBarry Smith    or at runtime via the option
94330de9b25SBarry Smith $     -pc_ilu_mat_ordering_type my_order
94430de9b25SBarry Smith $     -pc_lu_mat_ordering_type my_order
94530de9b25SBarry Smith 
946ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
94730de9b25SBarry Smith 
94830de9b25SBarry Smith .keywords: matrix, ordering, register
94930de9b25SBarry Smith 
95030de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
95130de9b25SBarry Smith M*/
952aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
953f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
954b12f92e5SBarry Smith #else
955f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
956b12f92e5SBarry Smith #endif
95730de9b25SBarry Smith 
958be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
959be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
9602bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
961b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
962d4fbbf0eSBarry Smith 
963be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
964a2ce50c7SBarry Smith 
965d91e6319SBarry Smith /*S
96615e8a5b3SHong Zhang    MatFactorInfo - Data based into the matrix factorization routines
967b00f7748SHong Zhang 
96815e8a5b3SHong Zhang    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE
969b00f7748SHong Zhang 
97015e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
971b00f7748SHong Zhang 
972b00f7748SHong Zhang    Level: developer
973b00f7748SHong Zhang 
974d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
975d7d82daaSBarry Smith           MatFactorInfoInitialize()
976b00f7748SHong Zhang 
977b00f7748SHong Zhang S*/
978b00f7748SHong Zhang typedef struct {
9790a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
9800a29876aSHong Zhang   PetscTruth    shiftpd;         /* if true, shift until positive pivots */
9812cea7109SBarry Smith   PetscReal     shift_fraction; /* record shift fraction taken */
98215e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
98315e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
984b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
98515e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
986f6275e2eSBarry Smith   PetscReal     fill;           /* expected fill; nonzeros in factored matrix/nonzeros in original matrix*/
987348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
988bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
989bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
99015e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
99115e8a5b3SHong Zhang } MatFactorInfo;
992ffa6d0a5SLois Curfman McInnes 
993be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
994be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,MatFactorInfo*);
995be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
996be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,MatFactorInfo*,Mat*);
997be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,MatFactorInfo*);
998be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,MatFactorInfo*);
999be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1000be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1001be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1002be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,MatFactorInfo*);
1003be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,MatFactorInfo*,Mat*);
1004be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,MatFactorInfo*,Mat *);
1005be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1006be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1007be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1008be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1009be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1010be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1011be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1012be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
10138ed539a5SBarry Smith 
1014be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1015bb5a7306SBarry Smith 
1016d91e6319SBarry Smith /*E
1017d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1018bb1eb677SSatish Balay 
1019d91e6319SBarry Smith     Level: beginner
1020d91e6319SBarry Smith 
1021d9274352SBarry Smith    May be bitwise ORd together
1022d9274352SBarry Smith 
1023d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1024d91e6319SBarry Smith 
10254e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10264e7234bfSBarry Smith 
1027a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax()
1028d91e6319SBarry Smith E*/
1029ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1030ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1031ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
103284cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1033be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1034be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10358ed539a5SBarry Smith 
1036d4fbbf0eSBarry Smith /*
1037639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1038639f9d9dSBarry Smith */
1039b12f92e5SBarry Smith 
1040d9274352SBarry Smith /*E
1041d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1042d9274352SBarry Smith        with an optional dynamic library name, for example
1043d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1044d9274352SBarry Smith 
1045d9274352SBarry Smith    Level: beginner
1046d9274352SBarry Smith 
1047d9274352SBarry Smith .seealso: MatGetColoring()
1048d9274352SBarry Smith E*/
104949773a63SBarry Smith #define MatColoringType char*
1050b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1051b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1052b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1053b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1054b12f92e5SBarry Smith 
1055be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*);
1056be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatColoringType,ISColoring *));
105730de9b25SBarry Smith 
105830de9b25SBarry Smith /*MC
105930de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
106030de9b25SBarry Smith                                matrix package.
106130de9b25SBarry Smith 
106230de9b25SBarry Smith    Synopsis:
1063c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
106430de9b25SBarry Smith 
106530de9b25SBarry Smith    Not Collective
106630de9b25SBarry Smith 
106730de9b25SBarry Smith    Input Parameters:
106830de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
106930de9b25SBarry Smith .  path - location of library where creation routine is
107030de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
107130de9b25SBarry Smith -  function - function pointer that creates the coloring
107230de9b25SBarry Smith 
107330de9b25SBarry Smith    Level: developer
107430de9b25SBarry Smith 
107530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
107630de9b25SBarry Smith    is ignored.
107730de9b25SBarry Smith 
107830de9b25SBarry Smith    Sample usage:
107930de9b25SBarry Smith .vb
108030de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
108130de9b25SBarry Smith                "MyColor",MyColor);
108230de9b25SBarry Smith .ve
108330de9b25SBarry Smith 
108430de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
108530de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
108630de9b25SBarry Smith    or at runtime via the option
108730de9b25SBarry Smith $     -mat_coloring_type my_color
108830de9b25SBarry Smith 
1089ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
109030de9b25SBarry Smith 
109130de9b25SBarry Smith .keywords: matrix, Coloring, register
109230de9b25SBarry Smith 
109330de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
109430de9b25SBarry Smith M*/
1095aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1096f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1097b12f92e5SBarry Smith #else
1098f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1099b12f92e5SBarry Smith #endif
110030de9b25SBarry Smith 
11012bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1102f1144a30SSatish Balay 
1103be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1104be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1105be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1106639f9d9dSBarry Smith 
1107d9274352SBarry Smith /*S
1108d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1109d9274352SBarry Smith         and coloring
1110639f9d9dSBarry Smith 
1111d9274352SBarry Smith    Level: beginner
1112d9274352SBarry Smith 
1113d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1114d9274352SBarry Smith 
1115d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1116d9274352SBarry Smith S*/
1117e2a1c21fSSatish Balay typedef struct _p_MatFDColoring *MatFDColoring;
1118639f9d9dSBarry Smith 
1119be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1120be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1121be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1122be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1123be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1124be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring,PetscInt);
1125be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring,PetscInt*);
1126be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1127be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1128be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1129be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring);
1130be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1131be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1132639f9d9dSBarry Smith /*
11330752156aSBarry Smith     These routines are for partitioning matrices: currently used only
11343eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
11350752156aSBarry Smith */
1136ca161407SBarry Smith 
1137d9274352SBarry Smith /*S
1138d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1139d9274352SBarry Smith 
1140d9274352SBarry Smith    Level: beginner
1141d9274352SBarry Smith 
1142d9274352SBarry Smith   Concepts: partitioning
1143d9274352SBarry Smith 
1144743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1145d9274352SBarry Smith S*/
114691e9ee9fSBarry Smith typedef struct _p_MatPartitioning *MatPartitioning;
1147d9274352SBarry Smith 
1148d9274352SBarry Smith /*E
11495bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1150d9274352SBarry Smith        with an optional dynamic library name, for example
1151d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1152d9274352SBarry Smith 
1153d9274352SBarry Smith    Level: beginner
1154d9274352SBarry Smith 
1155b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1156d9274352SBarry Smith E*/
115749773a63SBarry Smith #define MatPartitioningType char*
11588ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1159c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
11608ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
116171306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
116271306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
116371306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
116471306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
116571306c60Sjroman 
1166ca161407SBarry Smith 
1167be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1168be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1169be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1170be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1171be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1172be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1173be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1174be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
11752aabb6bbSBarry Smith 
1176be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
117730de9b25SBarry Smith 
117830de9b25SBarry Smith /*MC
117930de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
118030de9b25SBarry Smith    matrix package.
118130de9b25SBarry Smith 
118230de9b25SBarry Smith    Synopsis:
1183c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
118430de9b25SBarry Smith 
118530de9b25SBarry Smith    Not Collective
118630de9b25SBarry Smith 
118730de9b25SBarry Smith    Input Parameters:
118830de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
118930de9b25SBarry Smith .  path - location of library where creation routine is
119030de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
119130de9b25SBarry Smith -  function - function pointer that creates the partitioning type
119230de9b25SBarry Smith 
119330de9b25SBarry Smith    Level: developer
119430de9b25SBarry Smith 
119530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
119630de9b25SBarry Smith    is ignored.
119730de9b25SBarry Smith 
119830de9b25SBarry Smith    Sample usage:
119930de9b25SBarry Smith .vb
120030de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
120130de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
120230de9b25SBarry Smith .ve
120330de9b25SBarry Smith 
120430de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
120530de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
120630de9b25SBarry Smith    or at runtime via the option
120730de9b25SBarry Smith $     -mat_partitioning_type my_part
120830de9b25SBarry Smith 
1209ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
121030de9b25SBarry Smith 
121130de9b25SBarry Smith .keywords: matrix, partitioning, register
121230de9b25SBarry Smith 
121330de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
121430de9b25SBarry Smith M*/
1215aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1216f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
12172aabb6bbSBarry Smith #else
1218f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
12192aabb6bbSBarry Smith #endif
12202aabb6bbSBarry Smith 
12212bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1222f1144a30SSatish Balay 
1223be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1224be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
12252bad1931SBarry Smith 
1226be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1227be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1228be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1229ca161407SBarry Smith 
1230be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
12310752156aSBarry Smith 
1232be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1233be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
123471306c60Sjroman 
123571306c60Sjroman typedef enum { MP_CHACO_MULTILEVEL_KL, MP_CHACO_SPECTRAL, MP_CHACO_LINEAR,
123671306c60Sjroman     MP_CHACO_RANDOM, MP_CHACO_SCATTERED } MPChacoGlobalType;
1237be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
123871306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1239be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1240be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
124171306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1242be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1243be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1244be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
124571306c60Sjroman 
124671306c60Sjroman #define MP_PARTY_OPT "opt"
124771306c60Sjroman #define MP_PARTY_LIN "lin"
124871306c60Sjroman #define MP_PARTY_SCA "sca"
124971306c60Sjroman #define MP_PARTY_RAN "ran"
125071306c60Sjroman #define MP_PARTY_GBF "gbf"
125171306c60Sjroman #define MP_PARTY_GCF "gcf"
125271306c60Sjroman #define MP_PARTY_BUB "bub"
125371306c60Sjroman #define MP_PARTY_DEF "def"
1254be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
125571306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
125671306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
125771306c60Sjroman #define MP_PARTY_NONE "no"
1258be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1259be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1260be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1261be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
126271306c60Sjroman 
126371306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1264be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1265be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1266be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1267be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
126971306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1270be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1271be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1272be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
127371306c60Sjroman 
12740752156aSBarry Smith /*
12750a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1276d4fbbf0eSBarry Smith */
12771c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
12781c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
12791c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
12801c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
12811c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
12827c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
12837c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
12841c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
12851c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
12867c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
12877c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
12881c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
12891c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
12901c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
12911c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
12921c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
12931c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
12941c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
12951c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
12961c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
12971c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
12981c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
12991c1c02c0SLois Curfman McInnes                MATOP_COMPRESS=22,
13001c1c02c0SLois Curfman McInnes                MATOP_SET_OPTION=23,
13011c1c02c0SLois Curfman McInnes                MATOP_ZERO_ENTRIES=24,
13021c1c02c0SLois Curfman McInnes                MATOP_ZERO_ROWS=25,
13031c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_SYMBOLIC=26,
13041c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_NUMERIC=27,
13051c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_SYMBOLIC=28,
13061c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_NUMERIC=29,
1307d643ce63SMatthew Knepley                MATOP_SETUP_PREALLOCATION=30,
1308d643ce63SMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=31,
1309d643ce63SMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=32,
1310d643ce63SMatthew Knepley                MATOP_GET_ARRAY=33,
1311d643ce63SMatthew Knepley                MATOP_RESTORE_ARRAY=34,
1312d643ce63SMatthew Knepley                MATOP_DUPLCIATE=35,
1313d643ce63SMatthew Knepley                MATOP_FORWARD_SOLVE=36,
1314d643ce63SMatthew Knepley                MATOP_BACKWARD_SOLVE=37,
1315d643ce63SMatthew Knepley                MATOP_ILUFACTOR=38,
1316d643ce63SMatthew Knepley                MATOP_ICCFACTOR=39,
1317d643ce63SMatthew Knepley                MATOP_AXPY=40,
1318d643ce63SMatthew Knepley                MATOP_GET_SUBMATRICES=41,
1319d643ce63SMatthew Knepley                MATOP_INCREASE_OVERLAP=42,
1320d643ce63SMatthew Knepley                MATOP_GET_VALUES=43,
1321d643ce63SMatthew Knepley                MATOP_COPY=44,
1322d643ce63SMatthew Knepley                MATOP_PRINT_HELP=45,
1323d643ce63SMatthew Knepley                MATOP_SCALE=46,
1324d643ce63SMatthew Knepley                MATOP_SHIFT=47,
1325d643ce63SMatthew Knepley                MATOP_DIAGONAL_SHIFT=48,
1326d643ce63SMatthew Knepley                MATOP_ILUDT_FACTOR=49,
1327d643ce63SMatthew Knepley                MATOP_GET_BLOCK_SIZE=50,
1328d643ce63SMatthew Knepley                MATOP_GET_ROW_IJ=51,
1329d643ce63SMatthew Knepley                MATOP_RESTORE_ROW_IJ=52,
1330d643ce63SMatthew Knepley                MATOP_GET_COLUMN_IJ=53,
1331d643ce63SMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=54,
1332d643ce63SMatthew Knepley                MATOP_FDCOLORING_CREATE=55,
1333d643ce63SMatthew Knepley                MATOP_COLORING_PATCH=56,
1334d643ce63SMatthew Knepley                MATOP_SET_UNFACTORED=57,
1335d643ce63SMatthew Knepley                MATOP_PERMUTE=58,
1336d643ce63SMatthew Knepley                MATOP_SET_VALUES_BLOCKED=59,
1337d643ce63SMatthew Knepley                MATOP_GET_SUBMATRIX=60,
1338d643ce63SMatthew Knepley                MATOP_DESTROY=61,
1339d643ce63SMatthew Knepley                MATOP_VIEW=62,
1340d643ce63SMatthew Knepley                MATOP_GET_MAPS=63,
1341d643ce63SMatthew Knepley                MATOP_USE_SCALED_FORM=64,
1342d643ce63SMatthew Knepley                MATOP_SCALE_SYSTEM=65,
1343d643ce63SMatthew Knepley                MATOP_UNSCALE_SYSTEM=66,
1344ba16a8cbSKris Buschelman                MATOP_SET_LOCAL_TO_GLOBAL_MAP=67,
1345d643ce63SMatthew Knepley                MATOP_SET_VALUES_LOCAL=68,
1346d643ce63SMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=69,
1347d643ce63SMatthew Knepley                MATOP_GET_ROW_MAX=70,
1348d643ce63SMatthew Knepley                MATOP_CONVERT=71,
1349d643ce63SMatthew Knepley                MATOP_SET_COLORING=72,
1350d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1351d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1352d643ce63SMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1353d643ce63SMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1354ba16a8cbSKris Buschelman                MATOP_MULT_CON=77,
1355ba16a8cbSKris Buschelman                MATOP_MULT_TRANSPOSE_CON=78,
1356ba16a8cbSKris Buschelman                MATOP_ILU_FACTOR_SYMBOLIC_CON=79,
1357d643ce63SMatthew Knepley                MATOP_PERMUTE_SPARSIFY=80,
1358d643ce63SMatthew Knepley                MATOP_MULT_MULTIPLE=81,
135941acf15aSKris Buschelman                MATOP_SOLVE_MULTIPLE=82,
136041acf15aSKris Buschelman                MATOP_GET_INERTIA=83,
136141acf15aSKris Buschelman                MATOP_LOAD=84,
136241acf15aSKris Buschelman                MATOP_IS_SYMMETRIC=85,
136341acf15aSKris Buschelman                MATOP_IS_HERMITIAN=86,
136441acf15aSKris Buschelman                MATOP_IS_STRUCTURALLY_SYMMETRIC=87,
136541acf15aSKris Buschelman                MATOP_PB_RELAX=88,
136641acf15aSKris Buschelman                MATOP_GET_VECS=89,
136741acf15aSKris Buschelman                MATOP_MAT_MULT=90,
136841acf15aSKris Buschelman                MATOP_MAT_MULT_SYMBOLIC=91,
136941acf15aSKris Buschelman                MATOP_MAT_MULT_NUMERIC=92,
137041acf15aSKris Buschelman                MATOP_PTAP=93,
137141acf15aSKris Buschelman                MATOP_PTAP_SYMBOLIC=94,
1372bc011b1eSHong Zhang                MATOP_PTAP_NUMERIC=95,
1373bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE=96,
1374bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97,
13757ba1a0bfSKris Buschelman                MATOP_MAT_MULTTRANSPOSE_NUMERIC=98,
13767ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_SEQAIJ=99,
13777ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_SEQAIJ=100,
13787ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_MPIAIJ=101,
13797ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_MPIAIJ=102
1380fae171e0SBarry Smith              } MatOperation;
1381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1382be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1384be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1385112a2221SBarry Smith 
138690ace30eSBarry Smith /*
138790ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
138890ace30eSBarry Smith  stored in a universal format. By changing the format with
1389fb9695e5SSatish Balay  PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will
139090ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
139190ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
139290ace30eSBarry Smith  read into matrices of the same time.
139390ace30eSBarry Smith */
139490ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
139590ace30eSBarry Smith 
1396be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsGetColor(Mat,ISColoring *);
1398860d1616SSatish Balay 
1399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
14001f4e1ec7SBarry Smith 
1401d9274352SBarry Smith /*S
1402d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1403d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1404d9274352SBarry Smith 
1405f7a9e4ceSBarry Smith    Level: advanced
1406d9274352SBarry Smith 
1407d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1408d9274352SBarry Smith 
14096e1639daSBarry Smith   Users manual sections:
14106e1639daSBarry Smith .   sec_singular
14116e1639daSBarry Smith 
1412d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1413d9274352SBarry Smith S*/
141474637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1415d9274352SBarry Smith 
1416be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1417be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1418be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1419be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
1420be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat);
142174637425SBarry Smith 
1422be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1423be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1424be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
142546129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
14263f1d51d7SBarry Smith 
1427be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1428be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1429be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1430c4f061fbSSatish Balay 
1431be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1432b0a32e0cSBarry Smith 
1433be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
143404f1ad80SBarry Smith 
1435be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1436be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
14377dbadf16SMatthew Knepley 
1438e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
14392eac72dbSBarry Smith #endif
1440