xref: /petsc/include/petscmat.h (revision 4b8d99ad717e76011c9a5e9ca1a9b87bebd9c096)
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"
70*4b8d99adSRichard Tran Mills #define MATSEQCSRPERM      "seqcsrperm"
71f69a0ea3SMatthew Knepley #define MatType const char*
72d91e6319SBarry Smith 
73c06d978dSMatthew Knepley /* Logging support */
74552e946dSBarry Smith #define    MAT_FILE_COOKIE 1211216    /* used to indicate matrices in binary files */
75be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE;
76be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
77be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE;
78be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE;
796849ba73SBarry Smith extern PetscEvent  MAT_Mult, MAT_MultMatrixFree, MAT_Mults, MAT_MultConstrained, MAT_MultAdd, MAT_MultTranspose;
806849ba73SBarry Smith extern PetscEvent  MAT_MultTransposeConstrained, MAT_MultTransposeAdd, MAT_Solve, MAT_Solves, MAT_SolveAdd, MAT_SolveTranspose;
816849ba73SBarry Smith extern PetscEvent  MAT_SolveTransposeAdd, MAT_Relax, MAT_ForwardSolve, MAT_BackwardSolve, MAT_LUFactor, MAT_LUFactorSymbolic;
826849ba73SBarry Smith extern PetscEvent  MAT_LUFactorNumeric, MAT_CholeskyFactor, MAT_CholeskyFactorSymbolic, MAT_CholeskyFactorNumeric, MAT_ILUFactor;
836849ba73SBarry Smith extern PetscEvent  MAT_ILUFactorSymbolic, MAT_ICCFactorSymbolic, MAT_Copy, MAT_Convert, MAT_Scale, MAT_AssemblyBegin;
846849ba73SBarry Smith extern PetscEvent  MAT_AssemblyEnd, MAT_SetValues, MAT_GetValues, MAT_GetRow, MAT_GetSubMatrices, MAT_GetColoring, MAT_GetOrdering;
856849ba73SBarry Smith extern PetscEvent  MAT_IncreaseOverlap, MAT_Partitioning, MAT_ZeroEntries, MAT_Load, MAT_View, MAT_AXPY, MAT_FDColoringCreate;
866849ba73SBarry Smith extern PetscEvent  MAT_FDColoringApply, MAT_Transpose, MAT_FDColoringFunction;
87cb00f407SKris Buschelman extern PetscEvent  MAT_MatMult, MAT_MatMultSymbolic, MAT_MatMultNumeric;
88534c1384SKris Buschelman extern PetscEvent  MAT_PtAP, MAT_PtAPSymbolic, MAT_PtAPNumeric;
89bc011b1eSHong Zhang extern PetscEvent  MAT_MatMultTranspose, MAT_MatMultTransposeSymbolic, MAT_MatMultTransposeNumeric;
90c06d978dSMatthew Knepley 
91ceb03754SKris Buschelman /*E
92ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
93ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
94ceb03754SKris Buschelman 
95ceb03754SKris Buschelman     Level: beginner
96ceb03754SKris Buschelman 
97ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
98ceb03754SKris Buschelman 
990c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
100ceb03754SKris Buschelman E*/
101ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse;
102ceb03754SKris Buschelman 
103be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(char *);
104c06d978dSMatthew Knepley 
105f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
106f69a0ea3SMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm,void *ctx),(comm,&A),Mat,A);
107f69a0ea3SMatthew Knepley PetscPolymorphicFunction(MatCreate,(void *ctx),(PETSC_COMM_WORLD,&A),Mat,A);
108f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
109f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,MatType);
110be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
111be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
112be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
113be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
11430de9b25SBarry Smith 
115f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
116f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
117f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
118f69a0ea3SMatthew Knepley 
11930de9b25SBarry Smith /*MC
12030de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
12130de9b25SBarry Smith 
12230de9b25SBarry Smith    Synopsis:
123c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
12430de9b25SBarry Smith 
12530de9b25SBarry Smith    Not Collective
12630de9b25SBarry Smith 
12730de9b25SBarry Smith    Input Parameters:
12830de9b25SBarry Smith +  name - name of a new user-defined matrix type
12930de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
13030de9b25SBarry Smith .  name_create - name of routine to create method context
13130de9b25SBarry Smith -  routine_create - routine to create method context
13230de9b25SBarry Smith 
13330de9b25SBarry Smith    Notes:
13430de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
13530de9b25SBarry Smith 
13630de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
13730de9b25SBarry Smith    is ignored.
13830de9b25SBarry Smith 
13930de9b25SBarry Smith    Sample usage:
14030de9b25SBarry Smith .vb
14130de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
14230de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
14330de9b25SBarry Smith .ve
14430de9b25SBarry Smith 
14530de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
14630de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
14730de9b25SBarry Smith    or at runtime via the option
14830de9b25SBarry Smith $     -mat_type my_mat
14930de9b25SBarry Smith 
15030de9b25SBarry Smith    Level: advanced
15130de9b25SBarry Smith 
152ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
15330de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
15430de9b25SBarry Smith 
15530de9b25SBarry Smith .keywords: Mat, register
15630de9b25SBarry Smith 
15730de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
15830de9b25SBarry Smith 
15930de9b25SBarry Smith M*/
160273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
161273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
162273d9f13SBarry Smith #else
163273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
16430de9b25SBarry Smith #endif
16530de9b25SBarry Smith 
166273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
167b0a32e0cSBarry Smith extern PetscFList MatList;
16828988994SBarry Smith 
169be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
170be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
171be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1720c451bc4SBarry Smith PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A);
1730c451bc4SBarry Smith PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A);
1740c451bc4SBarry Smith PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A);
1750c451bc4SBarry Smith PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A);
1760c451bc4SBarry Smith PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A));
1770c451bc4SBarry Smith PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A));
1780c451bc4SBarry Smith PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A));
179be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1800c451bc4SBarry Smith 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);
1810c451bc4SBarry Smith 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);
1820c451bc4SBarry Smith 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);
1830c451bc4SBarry Smith 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);
1840c451bc4SBarry Smith 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));
1850c451bc4SBarry Smith 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));
1860c451bc4SBarry Smith 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));
1870c451bc4SBarry Smith 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);
1880c451bc4SBarry Smith 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);
1890c451bc4SBarry Smith 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);
1900c451bc4SBarry Smith 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);
1910c451bc4SBarry Smith 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));
1920c451bc4SBarry Smith 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));
1930c451bc4SBarry Smith 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));
194be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
195be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*);
196be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*);
197be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1980c451bc4SBarry Smith PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,nz,nnz,&A),Mat,A);
1990c451bc4SBarry Smith PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A);
2000c451bc4SBarry Smith PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A);
2010c451bc4SBarry Smith PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A);
2020c451bc4SBarry Smith PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A));
2030c451bc4SBarry Smith PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A));
2040c451bc4SBarry Smith PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A));
205be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
2060c451bc4SBarry Smith 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);
2070c451bc4SBarry Smith 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);
2080c451bc4SBarry Smith 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);
2090c451bc4SBarry Smith 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);
2100c451bc4SBarry Smith 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));
2110c451bc4SBarry Smith 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));
2120c451bc4SBarry Smith 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));
2130c451bc4SBarry Smith 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);
2140c451bc4SBarry Smith 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);
2150c451bc4SBarry Smith 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);
2160c451bc4SBarry Smith 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);
2170c451bc4SBarry Smith 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));
2180c451bc4SBarry Smith 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));
2190c451bc4SBarry Smith 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));
220be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
221be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
2220c451bc4SBarry Smith PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,nz,nnz,&A),Mat,A);
2230c451bc4SBarry Smith PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A);
2240c451bc4SBarry Smith PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A);
2250c451bc4SBarry Smith PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A);
2260c451bc4SBarry Smith PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A));
2270c451bc4SBarry Smith PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A));
2280c451bc4SBarry Smith PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A));
229be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
2300c451bc4SBarry Smith 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);
2310c451bc4SBarry Smith 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);
2320c451bc4SBarry Smith 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);
2330c451bc4SBarry Smith 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);
2340c451bc4SBarry Smith 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));
2350c451bc4SBarry Smith 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));
2360c451bc4SBarry Smith 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));
2370c451bc4SBarry Smith 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);
2380c451bc4SBarry Smith 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);
2390c451bc4SBarry Smith 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);
2400c451bc4SBarry Smith 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);
2410c451bc4SBarry Smith 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));
2420c451bc4SBarry Smith 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));
2430c451bc4SBarry Smith 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));
244be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
2450c451bc4SBarry Smith PetscPolymorphicFunction(MatCreateShell,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(comm,m,n,M,N,ctx,&A),Mat,A);
2460c451bc4SBarry Smith PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A);
247be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
248be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
2490c451bc4SBarry Smith PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A);
250f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
251be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
25221c89e3eSBarry Smith 
253be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPrintHelp(Mat);
254be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetPetscMaps(Mat,PetscMap*,PetscMap*);
255ec0117caSBarry Smith 
2560c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
2578ed539a5SBarry Smith /* ------------------------------------------------------------*/
258be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
259be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
26084cb2905SBarry Smith 
2612ef4de8bSBarry Smith /*S
2622ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
2632ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
2642ef4de8bSBarry Smith 
2652ef4de8bSBarry Smith    Level: beginner
2662ef4de8bSBarry Smith 
2672ef4de8bSBarry Smith   Concepts: matrix; linear operator
2682ef4de8bSBarry Smith 
269d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
2702ef4de8bSBarry Smith S*/
271435da068SBarry Smith typedef struct {
272c1ac3661SBarry Smith   PetscInt k,j,i,c;
273435da068SBarry Smith } MatStencil;
2742ef4de8bSBarry Smith 
275be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
277be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
278435da068SBarry Smith 
279be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
280be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
281be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
2823a7fca6bSBarry Smith 
283d91e6319SBarry Smith /*E
284d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
285d91e6319SBarry Smith      to continue to add values to it
286d91e6319SBarry Smith 
287d91e6319SBarry Smith     Level: beginner
288d91e6319SBarry Smith 
289d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
290d91e6319SBarry Smith E*/
2916d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
292be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
294be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
2954f9c727eSBarry Smith 
296be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt    MatSetValue_Row;
297be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt    MatSetValue_Column;
298be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscScalar MatSetValue_Value;
2991d73ed98SBarry Smith 
30030de9b25SBarry Smith /*MC
30130de9b25SBarry Smith    MatSetValue - Set a single entry into a matrix.
30230de9b25SBarry Smith 
30330de9b25SBarry Smith    Synopsis:
304c1ac3661SBarry Smith    PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode);
30530de9b25SBarry Smith 
30630de9b25SBarry Smith    Not collective
30730de9b25SBarry Smith 
30830de9b25SBarry Smith    Input Parameters:
30930de9b25SBarry Smith +  m - the matrix
31030de9b25SBarry Smith .  row - the row location of the entry
31130de9b25SBarry Smith .  col - the column location of the entry
31230de9b25SBarry Smith .  value - the value to insert
31330de9b25SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
31430de9b25SBarry Smith 
31530de9b25SBarry Smith    Notes:
31630de9b25SBarry Smith    For efficiency one should use MatSetValues() and set several or many
31730de9b25SBarry Smith    values simultaneously if possible.
31830de9b25SBarry Smith 
31930de9b25SBarry Smith    Level: beginner
32030de9b25SBarry Smith 
32130de9b25SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
32230de9b25SBarry Smith M*/
323b951964fSBarry Smith #define MatSetValue(v,i,j,va,mode) \
324f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \
3251d73ed98SBarry Smith    MatSetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode))
32630de9b25SBarry Smith 
327ea06a074SBarry Smith #define MatGetValue(v,i,j,va) \
328f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,0) || \
3291d73ed98SBarry Smith    MatGetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&va))
33030de9b25SBarry Smith 
331d91e6319SBarry Smith #define MatSetValueLocal(v,i,j,va,mode) \
332f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \
3331d73ed98SBarry Smith    MatSetValuesLocal(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode))
33430de9b25SBarry Smith 
335d91e6319SBarry Smith /*E
336d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
337d91e6319SBarry Smith 
338d91e6319SBarry Smith     Level: beginner
339d91e6319SBarry Smith 
3400a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
341d91e6319SBarry Smith 
342d91e6319SBarry Smith .seealso: MatSetOption()
343d91e6319SBarry Smith E*/
3446d4a8577SBarry Smith typedef enum {MAT_ROW_ORIENTED=1,MAT_COLUMN_ORIENTED=2,MAT_ROWS_SORTED=4,
3456d4a8577SBarry Smith               MAT_COLUMNS_SORTED=8,MAT_NO_NEW_NONZERO_LOCATIONS=16,
3466d4a8577SBarry Smith               MAT_YES_NEW_NONZERO_LOCATIONS=32,MAT_SYMMETRIC=64,
3476ca9ecd3SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC=65,MAT_NO_NEW_DIAGONALS=66,
3486ca9ecd3SBarry Smith               MAT_YES_NEW_DIAGONALS=67,MAT_INODE_LIMIT_1=68,MAT_INODE_LIMIT_2=69,
3496ca9ecd3SBarry Smith               MAT_INODE_LIMIT_3=70,MAT_INODE_LIMIT_4=71,MAT_INODE_LIMIT_5=72,
3506ca9ecd3SBarry Smith               MAT_IGNORE_OFF_PROC_ENTRIES=73,MAT_ROWS_UNSORTED=74,
3514787f768SSatish Balay               MAT_COLUMNS_UNSORTED=75,MAT_NEW_NONZERO_LOCATION_ERR=76,
3527c922b88SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR=77,MAT_USE_HASH_TABLE=78,
3532bad1931SBarry Smith               MAT_KEEP_ZEROED_ROWS=79,MAT_IGNORE_ZERO_ENTRIES=80,MAT_USE_INODES=81,
354efcf0fc3SBarry Smith               MAT_DO_NOT_USE_INODES=82,MAT_NOT_SYMMETRIC=83,MAT_HERMITIAN=84,
3559a4540c5SBarry Smith               MAT_NOT_STRUCTURALLY_SYMMETRIC=85,MAT_NOT_HERMITIAN=86,
356d487561eSHong Zhang               MAT_SYMMETRY_ETERNAL=87,MAT_NOT_SYMMETRY_ETERNAL=88,
357941593c8SHong Zhang               MAT_USE_COMPRESSEDROW=89,MAT_DO_NOT_USE_COMPRESSEDROW=90,
358941593c8SHong Zhang               MAT_IGNORE_LOWER_TRIANGULAR=91,MAT_ERROR_LOWER_TRIANGULAR=92} MatOption;
359be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption);
360be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,MatType*);
3610c451bc4SBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),MatType,t);
36284cb2905SBarry Smith 
363be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
364be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
365be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
366be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
367be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
368be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
3700c451bc4SBarry Smith PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a);
371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
3730c451bc4SBarry Smith PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a);
374be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
3757b80b807SBarry Smith 
376be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
3780c451bc4SBarry Smith PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y));
379be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
3810c451bc4SBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t);
3820c451bc4SBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0.0,&t),PetscTruth,t);
383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
3840c451bc4SBarry Smith PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y));
385be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
386be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
3872eac72dbSBarry Smith 
388d91e6319SBarry Smith /*E
389d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
390d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
391d91e6319SBarry Smith 
392d91e6319SBarry Smith     Level: beginner
393d91e6319SBarry Smith 
394d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
395d91e6319SBarry Smith 
396d91e6319SBarry Smith .seealso: MatDuplicate()
397d91e6319SBarry Smith E*/
3982e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
3992e8a6d31SBarry Smith 
400be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvertRegister(const char[],const char[],const char[],PetscErrorCode (*)(Mat,MatType,MatReuse,Mat*));
401273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
402273d9f13SBarry Smith #define MatConvertRegisterDynamic(a,b,c,d) MatConvertRegister(a,b,c,0)
403273d9f13SBarry Smith #else
404273d9f13SBarry Smith #define MatConvertRegisterDynamic(a,b,c,d) MatConvertRegister(a,b,c,d)
405273d9f13SBarry Smith #endif
406be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvertRegisterAll(const char[]);
407be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvertRegisterDestroy(void);
408273d9f13SBarry Smith extern PetscTruth MatConvertRegisterAllCalled;
409b0a32e0cSBarry Smith extern PetscFList MatConvertList;
410f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,MatType,MatReuse,Mat*);
411f69a0ea3SMatthew Knepley PetscPolymorphicFunction(MatConvert,(Mat A,MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a);
412be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
4130c451bc4SBarry Smith PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a);
4140c451bc4SBarry Smith PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a);
41594a9d846SBarry Smith 
416d91e6319SBarry Smith /*E
417d91e6319SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
418d91e6319SBarry Smith 
419d91e6319SBarry Smith     Level: beginner
420d91e6319SBarry Smith 
421d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
422d91e6319SBarry Smith 
42394b7f48cSBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
424d91e6319SBarry Smith E*/
425c537a176SHong Zhang typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
426cb5b572fSBarry Smith 
427be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
428be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
429be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
4300c451bc4SBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t);
4310c451bc4SBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0.0,&t),PetscTruth,t);
432be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
4330c451bc4SBarry Smith PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t);
434be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscTruth*);
4350c451bc4SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,&t),PetscTruth,t);
436be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
437be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
438f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,MatType,Mat*);
439f69a0ea3SMatthew Knepley PetscPolymorphicFunction(MatLoad,(PetscViewer v,MatType t),(v,t,&a),Mat,a);
4407b80b807SBarry Smith 
441be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
442be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
443be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
444be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
445d4fbbf0eSBarry Smith 
446d91e6319SBarry Smith /*S
447d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
448d91e6319SBarry Smith 
449d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
450d91e6319SBarry Smith 
451d91e6319SBarry Smith    Level: intermediate
452d91e6319SBarry Smith 
453d91e6319SBarry Smith   Concepts: matrix^nonzero information
454d91e6319SBarry Smith 
455d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
456d91e6319SBarry Smith S*/
4574e220ebcSLois Curfman McInnes typedef struct {
458b0a32e0cSBarry Smith   PetscLogDouble rows_global,columns_global;         /* number of global rows and columns */
459b0a32e0cSBarry Smith   PetscLogDouble rows_local,columns_local;           /* number of local rows and columns */
460b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
461b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
462b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
463b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
464b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
465b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
466b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4674e220ebcSLois Curfman McInnes } MatInfo;
4684e220ebcSLois Curfman McInnes 
469d9274352SBarry Smith /*E
470d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
471d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
472d9274352SBarry Smith 
473d9274352SBarry Smith     Level: beginner
474d9274352SBarry Smith 
475d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
476d9274352SBarry Smith 
477d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
478d9274352SBarry Smith E*/
4797b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
480be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
481be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
482be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
483be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec);
484be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,Mat*);
4850c451bc4SBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,&t),Mat,t);
486be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
4870c451bc4SBarry Smith PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t);
488be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
489be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
490be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
491be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
4920c451bc4SBarry Smith PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t);
493be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
494be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
495be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
496be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
4977b80b807SBarry Smith 
498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
4990c451bc4SBarry Smith PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n);
500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
501f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
502f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
503f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
504f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5057b80b807SBarry Smith 
506be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
508be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5095ef9f2a5SBarry Smith 
510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5137b80b807SBarry Smith 
514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
515be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *);
517be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
519be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
520be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
521be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
522be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
523be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
524be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
525be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscScalar**,Mat*);
5268efafbd8SBarry Smith 
527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5287b80b807SBarry Smith 
529be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
53222440eb1SKris Buschelman 
533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
535be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
53622440eb1SKris Buschelman 
537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
538be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
540bc011b1eSHong Zhang 
541f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
542f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat);
543be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat);
5441c741599SHong Zhang 
545f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
546f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5477b80b807SBarry Smith 
548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
550f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
551f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
553be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
554052efed2SBarry Smith 
555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
55790f02eecSBarry Smith 
558be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
5590c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
5600c451bc4SBarry Smith PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y));
561be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
562be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
563bd481603SBarry Smith 
564bd481603SBarry Smith /*MC
565bd481603SBarry Smith    MatPreallocInitialize - Begins the block of code that will count the number of nonzeros per
566bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
567bd481603SBarry Smith 
568bd481603SBarry Smith    Synopsis:
569c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
570bd481603SBarry Smith 
571bd481603SBarry Smith    Collective on MPI_Comm
572bd481603SBarry Smith 
573bd481603SBarry Smith    Input Parameters:
574bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
575859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
576859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
577bd481603SBarry Smith 
578bd481603SBarry Smith    Output Parameters:
579bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
580bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
581bd481603SBarry Smith 
582bd481603SBarry Smith 
583bd481603SBarry Smith    Level: intermediate
584bd481603SBarry Smith 
585bd481603SBarry Smith    Notes:
586bd481603SBarry Smith    See the chapter in the users manual on performance for more details
587bd481603SBarry Smith 
5881d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
589bd481603SBarry Smith 
590bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
591bd481603SBarry Smith 
592bd481603SBarry Smith   Concepts: preallocation^Matrix
593bd481603SBarry Smith 
594bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
595bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
596bd481603SBarry Smith M*/
597c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
5987c922b88SBarry Smith { \
599c1ac3661SBarry Smith   PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
600c1ac3661SBarry Smith   _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\
601c1ac3661SBarry Smith   _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\
602c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
603c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp;
6047c922b88SBarry Smith 
605bd481603SBarry Smith /*MC
606bd481603SBarry Smith    MatPreallocSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
607bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
608bd481603SBarry Smith 
609bd481603SBarry Smith    Synopsis:
610c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
611bd481603SBarry Smith 
612bd481603SBarry Smith    Collective on MPI_Comm
613bd481603SBarry Smith 
614bd481603SBarry Smith    Input Parameters:
615bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
616859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
617859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
618bd481603SBarry Smith 
619bd481603SBarry Smith    Output Parameters:
620bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
621bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
622bd481603SBarry Smith 
623bd481603SBarry Smith 
624bd481603SBarry Smith    Level: intermediate
625bd481603SBarry Smith 
626bd481603SBarry Smith    Notes:
627bd481603SBarry Smith    See the chapter in the users manual on performance for more details
628bd481603SBarry Smith 
6291d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
630bd481603SBarry Smith 
631bd481603SBarry Smith   Concepts: preallocation^Matrix
632bd481603SBarry Smith 
633bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
634bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
635bd481603SBarry Smith M*/
636222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
637222b16d4SBarry Smith { \
638c1ac3661SBarry Smith   PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__end; \
639c1ac3661SBarry Smith   _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\
640c1ac3661SBarry Smith   _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\
641c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
642c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp;
643222b16d4SBarry Smith 
644bd481603SBarry Smith /*MC
645bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
646bd481603SBarry Smith        inserted using a local number of the rows and columns
647bd481603SBarry Smith 
648bd481603SBarry Smith    Synopsis:
649c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
650bd481603SBarry Smith 
651bd481603SBarry Smith    Not Collective
652bd481603SBarry Smith 
653bd481603SBarry Smith    Input Parameters:
654bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
655bd481603SBarry Smith .  nrows - the number of rows indicated
6561d73ed98SBarry Smith .  rows - the indices of the rows
657bd481603SBarry Smith .  ncols - the number of columns in the matrix
658bd481603SBarry Smith .  cols - the columns indicated
659bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
660bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
661bd481603SBarry Smith 
662bd481603SBarry Smith 
663bd481603SBarry Smith    Level: intermediate
664bd481603SBarry Smith 
665bd481603SBarry Smith    Notes:
666bd481603SBarry Smith    See the chapter in the users manual on performance for more details
667bd481603SBarry Smith 
6681d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
669bd481603SBarry Smith 
670bd481603SBarry Smith   Concepts: preallocation^Matrix
671bd481603SBarry Smith 
672bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
673bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
674bd481603SBarry Smith M*/
675c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
676c4f061fbSSatish Balay {\
677c1ac3661SBarry Smith   PetscInt __l;\
678ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
679ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
680c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
681ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
682c4f061fbSSatish Balay   }\
683c4f061fbSSatish Balay }
684c4f061fbSSatish Balay 
685bd481603SBarry Smith /*MC
686bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
687bd481603SBarry Smith        inserted using a local number of the rows and columns
688bd481603SBarry Smith 
689bd481603SBarry Smith    Synopsis:
690c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
691bd481603SBarry Smith 
692bd481603SBarry Smith    Not Collective
693bd481603SBarry Smith 
694bd481603SBarry Smith    Input Parameters:
695bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
696bd481603SBarry Smith .  nrows - the number of rows indicated
6971d73ed98SBarry Smith .  rows - the indices of the rows
698bd481603SBarry Smith .  ncols - the number of columns in the matrix
699bd481603SBarry Smith .  cols - the columns indicated
700bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
701bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
702bd481603SBarry Smith 
703bd481603SBarry Smith 
704bd481603SBarry Smith    Level: intermediate
705bd481603SBarry Smith 
706bd481603SBarry Smith    Notes:
707bd481603SBarry Smith    See the chapter in the users manual on performance for more details
708bd481603SBarry Smith 
709bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
710bd481603SBarry Smith 
711bd481603SBarry Smith   Concepts: preallocation^Matrix
712bd481603SBarry Smith 
713bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
714bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
715bd481603SBarry Smith M*/
716d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
717d3d32019SBarry Smith {\
718c1ac3661SBarry Smith   PetscInt __l;\
719d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
720d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
721d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
722d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
723d3d32019SBarry Smith   }\
724d3d32019SBarry Smith }
725d3d32019SBarry Smith 
726bd481603SBarry Smith /*MC
727bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
728bd481603SBarry Smith        inserted using a local number of the rows and columns
729bd481603SBarry Smith 
730bd481603SBarry Smith    Synopsis:
731c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
732bd481603SBarry Smith 
733bd481603SBarry Smith    Not Collective
734bd481603SBarry Smith 
735bd481603SBarry Smith    Input Parameters:
736bd481603SBarry Smith +  nrows - the number of rows indicated
7371d73ed98SBarry Smith .  rows - the indices of the rows
738bd481603SBarry Smith .  ncols - the number of columns in the matrix
739a50f8bf6SHong Zhang -  cols - the columns indicated
740a50f8bf6SHong Zhang 
741a50f8bf6SHong Zhang    Output Parameters:
742a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
743bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
744bd481603SBarry Smith 
745bd481603SBarry Smith 
746bd481603SBarry Smith    Level: intermediate
747bd481603SBarry Smith 
748bd481603SBarry Smith    Notes:
749bd481603SBarry Smith    See the chapter in the users manual on performance for more details
750bd481603SBarry Smith 
751bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
752bd481603SBarry Smith 
753bd481603SBarry Smith   Concepts: preallocation^Matrix
754bd481603SBarry Smith 
755bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
756bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
757bd481603SBarry Smith M*/
758c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
759c1ac3661SBarry Smith { PetscInt __i; \
7607c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
7617c922b88SBarry Smith     if (cols[__i] < __start || cols[__i] >= __end) onz[row - __rstart]++; \
7627c922b88SBarry Smith   }\
7637c922b88SBarry Smith   dnz[row - __rstart] = nc - onz[row - __rstart];\
7647c922b88SBarry Smith }
7657c922b88SBarry Smith 
766bd481603SBarry Smith /*MC
767bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
768bd481603SBarry Smith        inserted using a local number of the rows and columns
769bd481603SBarry Smith 
770bd481603SBarry Smith    Synopsis:
771c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
772bd481603SBarry Smith 
773bd481603SBarry Smith    Not Collective
774bd481603SBarry Smith 
775bd481603SBarry Smith    Input Parameters:
776bd481603SBarry Smith +  nrows - the number of rows indicated
7771d73ed98SBarry Smith .  rows - the indices of the rows
778bd481603SBarry Smith .  ncols - the number of columns in the matrix
779bd481603SBarry Smith .  cols - the columns indicated
780bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
781bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
782bd481603SBarry Smith 
783bd481603SBarry Smith 
784bd481603SBarry Smith    Level: intermediate
785bd481603SBarry Smith 
786bd481603SBarry Smith    Notes:
787bd481603SBarry Smith    See the chapter in the users manual on performance for more details
788bd481603SBarry Smith 
789bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
790bd481603SBarry Smith 
791bd481603SBarry Smith   Concepts: preallocation^Matrix
792bd481603SBarry Smith 
793bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
794bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
795bd481603SBarry Smith M*/
796d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
797c1ac3661SBarry Smith { PetscInt __i; \
798d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
799d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
800d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
801d3d32019SBarry Smith   }\
802d3d32019SBarry Smith }
803d3d32019SBarry Smith 
804bd481603SBarry Smith /*MC
805bd481603SBarry Smith    MatPreallocFinalize - Ends the block of code that will count the number of nonzeros per
806bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
807bd481603SBarry Smith 
808bd481603SBarry Smith    Synopsis:
809c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
810bd481603SBarry Smith 
811bd481603SBarry Smith    Collective on MPI_Comm
812bd481603SBarry Smith 
813bd481603SBarry Smith    Input Parameters:
814bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
815bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
816bd481603SBarry Smith 
817bd481603SBarry Smith 
818bd481603SBarry Smith    Level: intermediate
819bd481603SBarry Smith 
820bd481603SBarry Smith    Notes:
821bd481603SBarry Smith    See the chapter in the users manual on performance for more details
822bd481603SBarry Smith 
823bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
824bd481603SBarry Smith 
825bd481603SBarry Smith   Concepts: preallocation^Matrix
826bd481603SBarry Smith 
827bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
828bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
829bd481603SBarry Smith M*/
830ef66eb69SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree(dnz);CHKERRQ(_4_ierr);}
8317c922b88SBarry Smith 
832bd481603SBarry Smith 
833bd481603SBarry Smith 
8347b80b807SBarry Smith /* Routines unique to particular data structures */
835be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
8360c451bc4SBarry Smith PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t);
837be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBDiagGetData(Mat,PetscInt*,PetscInt*,PetscInt*[],PetscInt*[],PetscScalar***);
838698d4c6aSKris Buschelman 
839be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
840be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
841698d4c6aSKris Buschelman 
842be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
843be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
844be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
8457b80b807SBarry Smith 
846a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
847a96a251dSBarry Smith 
848be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
8490c451bc4SBarry Smith PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz));
850be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
8510c451bc4SBarry Smith PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz));
852be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
8530c451bc4SBarry Smith PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz));
854be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDensePreallocation(Mat,PetscScalar[]);
855be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]);
856be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
857273d9f13SBarry Smith 
858be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
8590c451bc4SBarry Smith PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz));
860be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
861be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
862be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
863be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
864be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDensePreallocation(Mat,PetscScalar[]);
865be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]);
866be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
867be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
868be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
869be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
870be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
871be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
872273d9f13SBarry Smith 
873be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
8741b807ce4Svictorle 
875be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
876be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
8772e8a6d31SBarry Smith 
878be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
8793a7fca6bSBarry Smith 
8807b80b807SBarry Smith /*
8817b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
88294b7f48cSBarry Smith   done through the KSP and PC interfaces.
8837b80b807SBarry Smith */
8847b80b807SBarry Smith 
885d9274352SBarry Smith /*E
886d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
887d9274352SBarry Smith        with an optional dynamic library name, for example
888d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
889d9274352SBarry Smith 
890d9274352SBarry Smith    Level: beginner
891d9274352SBarry Smith 
892d9274352SBarry Smith .seealso: MatGetOrdering()
893d9274352SBarry Smith E*/
89449773a63SBarry Smith #define MatOrderingType char*
895b12f92e5SBarry Smith #define MATORDERING_NATURAL   "natural"
896b12f92e5SBarry Smith #define MATORDERING_ND        "nd"
897b12f92e5SBarry Smith #define MATORDERING_1WD       "1wd"
898b12f92e5SBarry Smith #define MATORDERING_RCM       "rcm"
899b12f92e5SBarry Smith #define MATORDERING_QMD       "qmd"
900b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength"
90162152c8bSBarry Smith #define MATORDERING_DSC_ND    "dsc_nd"
90262152c8bSBarry Smith #define MATORDERING_DSC_MMD   "dsc_mmd"
90362152c8bSBarry Smith #define MATORDERING_DSC_MDF   "dsc_mdf"
904c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
905c06d978dSMatthew Knepley #define MATORDERING_IDENTITY  "identity"
906c06d978dSMatthew Knepley #define MATORDERING_REVERSE   "reverse"
907b12f92e5SBarry Smith 
908be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
909be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
910be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
91130de9b25SBarry Smith 
91230de9b25SBarry Smith /*MC
91330de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
91430de9b25SBarry Smith                                matrix package.
91530de9b25SBarry Smith 
91630de9b25SBarry Smith    Synopsis:
917c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
91830de9b25SBarry Smith 
91930de9b25SBarry Smith    Not Collective
92030de9b25SBarry Smith 
92130de9b25SBarry Smith    Input Parameters:
92230de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
92330de9b25SBarry Smith .  path - location of library where creation routine is
92430de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
92530de9b25SBarry Smith -  function - function pointer that creates the ordering
92630de9b25SBarry Smith 
92730de9b25SBarry Smith    Level: developer
92830de9b25SBarry Smith 
92930de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
93030de9b25SBarry Smith    is ignored.
93130de9b25SBarry Smith 
93230de9b25SBarry Smith    Sample usage:
93330de9b25SBarry Smith .vb
93430de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
93530de9b25SBarry Smith                "MyOrder",MyOrder);
93630de9b25SBarry Smith .ve
93730de9b25SBarry Smith 
93830de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
93930de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
94030de9b25SBarry Smith    or at runtime via the option
94130de9b25SBarry Smith $     -pc_ilu_mat_ordering_type my_order
94230de9b25SBarry Smith $     -pc_lu_mat_ordering_type my_order
94330de9b25SBarry Smith 
944ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
94530de9b25SBarry Smith 
94630de9b25SBarry Smith .keywords: matrix, ordering, register
94730de9b25SBarry Smith 
94830de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
94930de9b25SBarry Smith M*/
950aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
951f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
952b12f92e5SBarry Smith #else
953f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
954b12f92e5SBarry Smith #endif
95530de9b25SBarry Smith 
956be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
957be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
9582bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
959b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
960d4fbbf0eSBarry Smith 
961be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
962a2ce50c7SBarry Smith 
963d91e6319SBarry Smith /*S
96415e8a5b3SHong Zhang    MatFactorInfo - Data based into the matrix factorization routines
965b00f7748SHong Zhang 
96615e8a5b3SHong Zhang    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE
967b00f7748SHong Zhang 
96815e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
969b00f7748SHong Zhang 
970b00f7748SHong Zhang    Level: developer
971b00f7748SHong Zhang 
972d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
973d7d82daaSBarry Smith           MatFactorInfoInitialize()
974b00f7748SHong Zhang 
975b00f7748SHong Zhang S*/
976b00f7748SHong Zhang typedef struct {
9770a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
9780a29876aSHong Zhang   PetscTruth    shiftpd;         /* if true, shift until positive pivots */
9792cea7109SBarry Smith   PetscReal     shift_fraction; /* record shift fraction taken */
98015e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
98115e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
982b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
98315e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
984f6275e2eSBarry Smith   PetscReal     fill;           /* expected fill; nonzeros in factored matrix/nonzeros in original matrix*/
985348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
986bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
987bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
98815e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
98915e8a5b3SHong Zhang } MatFactorInfo;
990ffa6d0a5SLois Curfman McInnes 
991be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
992be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,MatFactorInfo*);
993be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
994be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,MatFactorInfo*,Mat*);
995be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,MatFactorInfo*);
996be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,MatFactorInfo*);
997be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
998be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
999be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1000be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,MatFactorInfo*);
1001be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,MatFactorInfo*,Mat*);
1002be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,MatFactorInfo*,Mat *);
1003be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1004be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1005be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1006be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1007be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1008be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1009be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1010be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
10118ed539a5SBarry Smith 
1012be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1013bb5a7306SBarry Smith 
1014d91e6319SBarry Smith /*E
1015d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1016bb1eb677SSatish Balay 
1017d91e6319SBarry Smith     Level: beginner
1018d91e6319SBarry Smith 
1019d9274352SBarry Smith    May be bitwise ORd together
1020d9274352SBarry Smith 
1021d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1022d91e6319SBarry Smith 
10234e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10244e7234bfSBarry Smith 
1025a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax()
1026d91e6319SBarry Smith E*/
1027ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1028ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1029ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
103084cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1031be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1032be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10338ed539a5SBarry Smith 
1034d4fbbf0eSBarry Smith /*
1035639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1036639f9d9dSBarry Smith */
1037b12f92e5SBarry Smith 
1038d9274352SBarry Smith /*E
1039d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1040d9274352SBarry Smith        with an optional dynamic library name, for example
1041d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1042d9274352SBarry Smith 
1043d9274352SBarry Smith    Level: beginner
1044d9274352SBarry Smith 
1045d9274352SBarry Smith .seealso: MatGetColoring()
1046d9274352SBarry Smith E*/
104749773a63SBarry Smith #define MatColoringType char*
1048b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1049b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1050b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1051b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1052b12f92e5SBarry Smith 
1053be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*);
1054be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatColoringType,ISColoring *));
105530de9b25SBarry Smith 
105630de9b25SBarry Smith /*MC
105730de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
105830de9b25SBarry Smith                                matrix package.
105930de9b25SBarry Smith 
106030de9b25SBarry Smith    Synopsis:
1061c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
106230de9b25SBarry Smith 
106330de9b25SBarry Smith    Not Collective
106430de9b25SBarry Smith 
106530de9b25SBarry Smith    Input Parameters:
106630de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
106730de9b25SBarry Smith .  path - location of library where creation routine is
106830de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
106930de9b25SBarry Smith -  function - function pointer that creates the coloring
107030de9b25SBarry Smith 
107130de9b25SBarry Smith    Level: developer
107230de9b25SBarry Smith 
107330de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
107430de9b25SBarry Smith    is ignored.
107530de9b25SBarry Smith 
107630de9b25SBarry Smith    Sample usage:
107730de9b25SBarry Smith .vb
107830de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
107930de9b25SBarry Smith                "MyColor",MyColor);
108030de9b25SBarry Smith .ve
108130de9b25SBarry Smith 
108230de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
108330de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
108430de9b25SBarry Smith    or at runtime via the option
108530de9b25SBarry Smith $     -mat_coloring_type my_color
108630de9b25SBarry Smith 
1087ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
108830de9b25SBarry Smith 
108930de9b25SBarry Smith .keywords: matrix, Coloring, register
109030de9b25SBarry Smith 
109130de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
109230de9b25SBarry Smith M*/
1093aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1094f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1095b12f92e5SBarry Smith #else
1096f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1097b12f92e5SBarry Smith #endif
109830de9b25SBarry Smith 
10992bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1100f1144a30SSatish Balay 
1101be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1102be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1103be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1104639f9d9dSBarry Smith 
1105d9274352SBarry Smith /*S
1106d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1107d9274352SBarry Smith         and coloring
1108639f9d9dSBarry Smith 
1109d9274352SBarry Smith    Level: beginner
1110d9274352SBarry Smith 
1111d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1112d9274352SBarry Smith 
1113d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1114d9274352SBarry Smith S*/
1115e2a1c21fSSatish Balay typedef struct _p_MatFDColoring *MatFDColoring;
1116639f9d9dSBarry Smith 
1117be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1118be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1119be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1120be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1121be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1122be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring,PetscInt);
1123be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring,PetscInt*);
1124be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1125be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1126be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1127be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring);
1128be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1129be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1130639f9d9dSBarry Smith /*
11310752156aSBarry Smith     These routines are for partitioning matrices: currently used only
11323eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
11330752156aSBarry Smith */
1134ca161407SBarry Smith 
1135d9274352SBarry Smith /*S
1136d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1137d9274352SBarry Smith 
1138d9274352SBarry Smith    Level: beginner
1139d9274352SBarry Smith 
1140d9274352SBarry Smith   Concepts: partitioning
1141d9274352SBarry Smith 
1142743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1143d9274352SBarry Smith S*/
114491e9ee9fSBarry Smith typedef struct _p_MatPartitioning *MatPartitioning;
1145d9274352SBarry Smith 
1146d9274352SBarry Smith /*E
11475bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1148d9274352SBarry Smith        with an optional dynamic library name, for example
1149d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1150d9274352SBarry Smith 
1151d9274352SBarry Smith    Level: beginner
1152d9274352SBarry Smith 
1153b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1154d9274352SBarry Smith E*/
115549773a63SBarry Smith #define MatPartitioningType char*
11568ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1157c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
11588ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
115971306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
116071306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
116171306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
116271306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
116371306c60Sjroman 
1164ca161407SBarry Smith 
1165be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1166be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1167be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1168be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1169be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1170be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1171be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1172be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
11732aabb6bbSBarry Smith 
1174be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
117530de9b25SBarry Smith 
117630de9b25SBarry Smith /*MC
117730de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
117830de9b25SBarry Smith    matrix package.
117930de9b25SBarry Smith 
118030de9b25SBarry Smith    Synopsis:
1181c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
118230de9b25SBarry Smith 
118330de9b25SBarry Smith    Not Collective
118430de9b25SBarry Smith 
118530de9b25SBarry Smith    Input Parameters:
118630de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
118730de9b25SBarry Smith .  path - location of library where creation routine is
118830de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
118930de9b25SBarry Smith -  function - function pointer that creates the partitioning type
119030de9b25SBarry Smith 
119130de9b25SBarry Smith    Level: developer
119230de9b25SBarry Smith 
119330de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
119430de9b25SBarry Smith    is ignored.
119530de9b25SBarry Smith 
119630de9b25SBarry Smith    Sample usage:
119730de9b25SBarry Smith .vb
119830de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
119930de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
120030de9b25SBarry Smith .ve
120130de9b25SBarry Smith 
120230de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
120330de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
120430de9b25SBarry Smith    or at runtime via the option
120530de9b25SBarry Smith $     -mat_partitioning_type my_part
120630de9b25SBarry Smith 
1207ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
120830de9b25SBarry Smith 
120930de9b25SBarry Smith .keywords: matrix, partitioning, register
121030de9b25SBarry Smith 
121130de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
121230de9b25SBarry Smith M*/
1213aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1214f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
12152aabb6bbSBarry Smith #else
1216f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
12172aabb6bbSBarry Smith #endif
12182aabb6bbSBarry Smith 
12192bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1220f1144a30SSatish Balay 
1221be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1222be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
12232bad1931SBarry Smith 
1224be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1225be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1226be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1227ca161407SBarry Smith 
1228be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
12290752156aSBarry Smith 
1230be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1231be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
123271306c60Sjroman 
123371306c60Sjroman typedef enum { MP_CHACO_MULTILEVEL_KL, MP_CHACO_SPECTRAL, MP_CHACO_LINEAR,
123471306c60Sjroman     MP_CHACO_RANDOM, MP_CHACO_SCATTERED } MPChacoGlobalType;
1235be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
123671306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1237be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1238be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
123971306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1240be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1241be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1242be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
124371306c60Sjroman 
124471306c60Sjroman #define MP_PARTY_OPT "opt"
124571306c60Sjroman #define MP_PARTY_LIN "lin"
124671306c60Sjroman #define MP_PARTY_SCA "sca"
124771306c60Sjroman #define MP_PARTY_RAN "ran"
124871306c60Sjroman #define MP_PARTY_GBF "gbf"
124971306c60Sjroman #define MP_PARTY_GCF "gcf"
125071306c60Sjroman #define MP_PARTY_BUB "bub"
125171306c60Sjroman #define MP_PARTY_DEF "def"
1252be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
125371306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
125471306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
125571306c60Sjroman #define MP_PARTY_NONE "no"
1256be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1257be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1258be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1259be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
126071306c60Sjroman 
126171306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1262be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1263be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1264be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1265be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1266be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
126771306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1269be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1270be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
127171306c60Sjroman 
12720752156aSBarry Smith /*
12730a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1274d4fbbf0eSBarry Smith */
12751c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
12761c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
12771c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
12781c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
12791c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
12807c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
12817c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
12821c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
12831c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
12847c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
12857c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
12861c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
12871c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
12881c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
12891c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
12901c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
12911c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
12921c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
12931c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
12941c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
12951c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
12961c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
12971c1c02c0SLois Curfman McInnes                MATOP_COMPRESS=22,
12981c1c02c0SLois Curfman McInnes                MATOP_SET_OPTION=23,
12991c1c02c0SLois Curfman McInnes                MATOP_ZERO_ENTRIES=24,
13001c1c02c0SLois Curfman McInnes                MATOP_ZERO_ROWS=25,
13011c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_SYMBOLIC=26,
13021c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_NUMERIC=27,
13031c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_SYMBOLIC=28,
13041c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_NUMERIC=29,
1305d643ce63SMatthew Knepley                MATOP_SETUP_PREALLOCATION=30,
1306d643ce63SMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=31,
1307d643ce63SMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=32,
1308d643ce63SMatthew Knepley                MATOP_GET_ARRAY=33,
1309d643ce63SMatthew Knepley                MATOP_RESTORE_ARRAY=34,
1310d643ce63SMatthew Knepley                MATOP_DUPLCIATE=35,
1311d643ce63SMatthew Knepley                MATOP_FORWARD_SOLVE=36,
1312d643ce63SMatthew Knepley                MATOP_BACKWARD_SOLVE=37,
1313d643ce63SMatthew Knepley                MATOP_ILUFACTOR=38,
1314d643ce63SMatthew Knepley                MATOP_ICCFACTOR=39,
1315d643ce63SMatthew Knepley                MATOP_AXPY=40,
1316d643ce63SMatthew Knepley                MATOP_GET_SUBMATRICES=41,
1317d643ce63SMatthew Knepley                MATOP_INCREASE_OVERLAP=42,
1318d643ce63SMatthew Knepley                MATOP_GET_VALUES=43,
1319d643ce63SMatthew Knepley                MATOP_COPY=44,
1320d643ce63SMatthew Knepley                MATOP_PRINT_HELP=45,
1321d643ce63SMatthew Knepley                MATOP_SCALE=46,
1322d643ce63SMatthew Knepley                MATOP_SHIFT=47,
1323d643ce63SMatthew Knepley                MATOP_DIAGONAL_SHIFT=48,
1324d643ce63SMatthew Knepley                MATOP_ILUDT_FACTOR=49,
1325d643ce63SMatthew Knepley                MATOP_GET_BLOCK_SIZE=50,
1326d643ce63SMatthew Knepley                MATOP_GET_ROW_IJ=51,
1327d643ce63SMatthew Knepley                MATOP_RESTORE_ROW_IJ=52,
1328d643ce63SMatthew Knepley                MATOP_GET_COLUMN_IJ=53,
1329d643ce63SMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=54,
1330d643ce63SMatthew Knepley                MATOP_FDCOLORING_CREATE=55,
1331d643ce63SMatthew Knepley                MATOP_COLORING_PATCH=56,
1332d643ce63SMatthew Knepley                MATOP_SET_UNFACTORED=57,
1333d643ce63SMatthew Knepley                MATOP_PERMUTE=58,
1334d643ce63SMatthew Knepley                MATOP_SET_VALUES_BLOCKED=59,
1335d643ce63SMatthew Knepley                MATOP_GET_SUBMATRIX=60,
1336d643ce63SMatthew Knepley                MATOP_DESTROY=61,
1337d643ce63SMatthew Knepley                MATOP_VIEW=62,
1338d643ce63SMatthew Knepley                MATOP_GET_MAPS=63,
1339d643ce63SMatthew Knepley                MATOP_USE_SCALED_FORM=64,
1340d643ce63SMatthew Knepley                MATOP_SCALE_SYSTEM=65,
1341d643ce63SMatthew Knepley                MATOP_UNSCALE_SYSTEM=66,
1342ba16a8cbSKris Buschelman                MATOP_SET_LOCAL_TO_GLOBAL_MAP=67,
1343d643ce63SMatthew Knepley                MATOP_SET_VALUES_LOCAL=68,
1344d643ce63SMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=69,
1345d643ce63SMatthew Knepley                MATOP_GET_ROW_MAX=70,
1346d643ce63SMatthew Knepley                MATOP_CONVERT=71,
1347d643ce63SMatthew Knepley                MATOP_SET_COLORING=72,
1348d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1349d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1350d643ce63SMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1351d643ce63SMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1352ba16a8cbSKris Buschelman                MATOP_MULT_CON=77,
1353ba16a8cbSKris Buschelman                MATOP_MULT_TRANSPOSE_CON=78,
1354ba16a8cbSKris Buschelman                MATOP_ILU_FACTOR_SYMBOLIC_CON=79,
1355d643ce63SMatthew Knepley                MATOP_PERMUTE_SPARSIFY=80,
1356d643ce63SMatthew Knepley                MATOP_MULT_MULTIPLE=81,
135741acf15aSKris Buschelman                MATOP_SOLVE_MULTIPLE=82,
135841acf15aSKris Buschelman                MATOP_GET_INERTIA=83,
135941acf15aSKris Buschelman                MATOP_LOAD=84,
136041acf15aSKris Buschelman                MATOP_IS_SYMMETRIC=85,
136141acf15aSKris Buschelman                MATOP_IS_HERMITIAN=86,
136241acf15aSKris Buschelman                MATOP_IS_STRUCTURALLY_SYMMETRIC=87,
136341acf15aSKris Buschelman                MATOP_PB_RELAX=88,
136441acf15aSKris Buschelman                MATOP_GET_VECS=89,
136541acf15aSKris Buschelman                MATOP_MAT_MULT=90,
136641acf15aSKris Buschelman                MATOP_MAT_MULT_SYMBOLIC=91,
136741acf15aSKris Buschelman                MATOP_MAT_MULT_NUMERIC=92,
136841acf15aSKris Buschelman                MATOP_PTAP=93,
136941acf15aSKris Buschelman                MATOP_PTAP_SYMBOLIC=94,
1370bc011b1eSHong Zhang                MATOP_PTAP_NUMERIC=95,
1371bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE=96,
1372bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97,
13737ba1a0bfSKris Buschelman                MATOP_MAT_MULTTRANSPOSE_NUMERIC=98,
13747ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_SEQAIJ=99,
13757ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_SEQAIJ=100,
13767ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_MPIAIJ=101,
13777ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_MPIAIJ=102
1378fae171e0SBarry Smith              } MatOperation;
1379be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1382be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1383112a2221SBarry Smith 
138490ace30eSBarry Smith /*
138590ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
138690ace30eSBarry Smith  stored in a universal format. By changing the format with
1387fb9695e5SSatish Balay  PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will
138890ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
138990ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
139090ace30eSBarry Smith  read into matrices of the same time.
139190ace30eSBarry Smith */
139290ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
139390ace30eSBarry Smith 
1394be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1395be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsGetColor(Mat,ISColoring *);
1396860d1616SSatish Balay 
1397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
13981f4e1ec7SBarry Smith 
1399d9274352SBarry Smith /*S
1400d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1401d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1402d9274352SBarry Smith 
1403f7a9e4ceSBarry Smith    Level: advanced
1404d9274352SBarry Smith 
1405d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1406d9274352SBarry Smith 
14076e1639daSBarry Smith   Users manual sections:
14086e1639daSBarry Smith .   sec_singular
14096e1639daSBarry Smith 
1410d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1411d9274352SBarry Smith S*/
141274637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1413d9274352SBarry Smith 
1414be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1415be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1416be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1417be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
1418be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat);
141974637425SBarry Smith 
1420be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1421be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1422be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
142346129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
14243f1d51d7SBarry Smith 
1425be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1426be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1427be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1428c4f061fbSSatish Balay 
1429be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1430b0a32e0cSBarry Smith 
1431be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
143204f1ad80SBarry Smith 
1433be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1434be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
14357dbadf16SMatthew Knepley 
1436e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
14372eac72dbSBarry Smith #endif
1438