xref: /petsc/include/petscmat.h (revision abf0e977b6b78a4d15cdf5bdbf4ebae62a1f2f72)
12eac72dbSBarry Smith /*
22eac72dbSBarry Smith      Include file for the matrix component of PETSc
32eac72dbSBarry Smith */
40a835dfdSSatish Balay #ifndef __PETSCMAT_H
50a835dfdSSatish Balay #define __PETSCMAT_H
60a835dfdSSatish Balay #include "petscvec.h"
7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN
82eac72dbSBarry Smith 
9d9274352SBarry Smith /*S
10d9274352SBarry Smith      Mat - Abstract PETSc matrix object
112eac72dbSBarry Smith 
12d91e6319SBarry Smith    Level: beginner
13d91e6319SBarry Smith 
14d9274352SBarry Smith   Concepts: matrix; linear operator
15d9274352SBarry Smith 
16d9274352SBarry Smith .seealso:  MatCreate(), MatType, MatSetType()
17d9274352SBarry Smith S*/
18d9274352SBarry Smith typedef struct _p_Mat*           Mat;
19d9274352SBarry Smith 
20d9274352SBarry Smith /*E
21d9274352SBarry Smith     MatType - String with the name of a PETSc matrix or the creation function
22d9274352SBarry Smith        with an optional dynamic library name, for example
23d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mymatcreate()
24d9274352SBarry Smith 
25d9274352SBarry Smith    Level: beginner
26d9274352SBarry Smith 
27d9274352SBarry Smith .seealso: MatSetType(), Mat
28d91e6319SBarry Smith E*/
29273d9f13SBarry Smith #define MATSAME            "same"
30273d9f13SBarry Smith #define MATSEQMAIJ         "seqmaij"
31273d9f13SBarry Smith #define MATMPIMAIJ         "mpimaij"
32209238afSKris Buschelman #define MATMAIJ            "maij"
33273d9f13SBarry Smith #define MATIS              "is"
34273d9f13SBarry Smith #define MATMPIROWBS        "mpirowbs"
35273d9f13SBarry Smith #define MATSEQAIJ          "seqaij"
36273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
37209238afSKris Buschelman #define MATAIJ             "aij"
38273d9f13SBarry Smith #define MATSHELL           "shell"
39273d9f13SBarry Smith #define MATSEQBDIAG        "seqbdiag"
40273d9f13SBarry Smith #define MATMPIBDIAG        "mpibdiag"
41209238afSKris Buschelman #define MATBDIAG           "bdiag"
42209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
43273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
44209238afSKris Buschelman #define MATDENSE           "dense"
45273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
46273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
47209238afSKris Buschelman #define MATBAIJ            "baij"
48273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
49273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
50273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
51209238afSKris Buschelman #define MATSBAIJ           "sbaij"
52cebc7f6cSBarry Smith #define MATDAAD            "daad"
53cebc7f6cSBarry Smith #define MATMFFD            "mffd"
54c8a8475eSBarry Smith #define MATNORMAL          "normal"
55ab92ecdeSBarry Smith #define MATLRC             "lrc"
56b3a1e11cSKris Buschelman #define MATSEQAIJSPOOLES   "seqaijspooles"
57d10c748bSKris Buschelman #define MATMPIAIJSPOOLES   "mpiaijspooles"
589abb65ffSKris Buschelman #define MATSEQSBAIJSPOOLES "seqsbaijspooles"
5922191285SKris Buschelman #define MATMPISBAIJSPOOLES "mpisbaijspooles"
60bb4d4166SHong Zhang #define MATAIJSPOOLES      "aijspooles"
61bb4d4166SHong Zhang #define MATSBAIJSPOOLES    "sbaijspooles"
6214b396bbSKris Buschelman #define MATSUPERLU         "superlu"
631d96aa28SKris Buschelman #define MATSUPERLU_DIST    "superlu_dist"
641677d0b8SKris Buschelman #define MATUMFPACK         "umfpack"
65e8aa55a4SKris Buschelman #define MATESSL            "essl"
664eb8e494SKris Buschelman #define MATLUSOL           "lusol"
67397b6df1SKris Buschelman #define MATAIJMUMPS        "aijmumps"
68397b6df1SKris Buschelman #define MATSBAIJMUMPS      "sbaijmumps"
698da957c5SKris Buschelman #define MATDSCPACK         "dscpack"
707065e2aaSKris Buschelman #define MATMATLAB          "matlab"
714b8d99adSRichard Tran Mills #define MATSEQCSRPERM      "seqcsrperm"
7218def556SRichard Tran Mills #define MATMPICSRPERM      "mpicsrperm"
73814655b8SBarry Smith #define MATCSRPERM         "csrperm"
74f69a0ea3SMatthew Knepley #define MatType const char*
75d91e6319SBarry Smith 
76c06d978dSMatthew Knepley /* Logging support */
77552e946dSBarry Smith #define    MAT_FILE_COOKIE 1211216    /* used to indicate matrices in binary files */
78be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE;
79be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
80be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE;
81be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE;
826849ba73SBarry Smith extern PetscEvent  MAT_Mult, MAT_MultMatrixFree, MAT_Mults, MAT_MultConstrained, MAT_MultAdd, MAT_MultTranspose;
836849ba73SBarry Smith extern PetscEvent  MAT_MultTransposeConstrained, MAT_MultTransposeAdd, MAT_Solve, MAT_Solves, MAT_SolveAdd, MAT_SolveTranspose;
846849ba73SBarry Smith extern PetscEvent  MAT_SolveTransposeAdd, MAT_Relax, MAT_ForwardSolve, MAT_BackwardSolve, MAT_LUFactor, MAT_LUFactorSymbolic;
856849ba73SBarry Smith extern PetscEvent  MAT_LUFactorNumeric, MAT_CholeskyFactor, MAT_CholeskyFactorSymbolic, MAT_CholeskyFactorNumeric, MAT_ILUFactor;
866849ba73SBarry Smith extern PetscEvent  MAT_ILUFactorSymbolic, MAT_ICCFactorSymbolic, MAT_Copy, MAT_Convert, MAT_Scale, MAT_AssemblyBegin;
876849ba73SBarry Smith extern PetscEvent  MAT_AssemblyEnd, MAT_SetValues, MAT_GetValues, MAT_GetRow, MAT_GetSubMatrices, MAT_GetColoring, MAT_GetOrdering;
886849ba73SBarry Smith extern PetscEvent  MAT_IncreaseOverlap, MAT_Partitioning, MAT_ZeroEntries, MAT_Load, MAT_View, MAT_AXPY, MAT_FDColoringCreate;
896849ba73SBarry Smith extern PetscEvent  MAT_FDColoringApply, MAT_Transpose, MAT_FDColoringFunction;
90cb00f407SKris Buschelman extern PetscEvent  MAT_MatMult, MAT_MatMultSymbolic, MAT_MatMultNumeric;
91534c1384SKris Buschelman extern PetscEvent  MAT_PtAP, MAT_PtAPSymbolic, MAT_PtAPNumeric;
92bc011b1eSHong Zhang extern PetscEvent  MAT_MatMultTranspose, MAT_MatMultTransposeSymbolic, MAT_MatMultTransposeNumeric;
93c06d978dSMatthew Knepley 
94ceb03754SKris Buschelman /*E
95ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
96ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
97ceb03754SKris Buschelman 
98ceb03754SKris Buschelman     Level: beginner
99ceb03754SKris Buschelman 
100ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
101ceb03754SKris Buschelman 
1020c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
103ceb03754SKris Buschelman E*/
104ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse;
105ceb03754SKris Buschelman 
106be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(char *);
107c06d978dSMatthew Knepley 
108f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
109ba966639SSatish Balay PetscPolymorphicFunction(MatCreate,(MPI_Comm comm,void *ctx),(comm,&A),Mat,A)
110ba966639SSatish Balay PetscPolymorphicFunction(MatCreate,(void *ctx),(PETSC_COMM_WORLD,&A),Mat,A)
111f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
112f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,MatType);
113be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
114be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
115be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
116be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
11730de9b25SBarry Smith 
118f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
119f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
120f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
121f69a0ea3SMatthew Knepley 
12230de9b25SBarry Smith /*MC
12330de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
12430de9b25SBarry Smith 
12530de9b25SBarry Smith    Synopsis:
126c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
12730de9b25SBarry Smith 
12830de9b25SBarry Smith    Not Collective
12930de9b25SBarry Smith 
13030de9b25SBarry Smith    Input Parameters:
13130de9b25SBarry Smith +  name - name of a new user-defined matrix type
13230de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
13330de9b25SBarry Smith .  name_create - name of routine to create method context
13430de9b25SBarry Smith -  routine_create - routine to create method context
13530de9b25SBarry Smith 
13630de9b25SBarry Smith    Notes:
13730de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
13830de9b25SBarry Smith 
13930de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
14030de9b25SBarry Smith    is ignored.
14130de9b25SBarry Smith 
14230de9b25SBarry Smith    Sample usage:
14330de9b25SBarry Smith .vb
14430de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
14530de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
14630de9b25SBarry Smith .ve
14730de9b25SBarry Smith 
14830de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
14930de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
15030de9b25SBarry Smith    or at runtime via the option
15130de9b25SBarry Smith $     -mat_type my_mat
15230de9b25SBarry Smith 
15330de9b25SBarry Smith    Level: advanced
15430de9b25SBarry Smith 
155ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
15630de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
15730de9b25SBarry Smith 
15830de9b25SBarry Smith .keywords: Mat, register
15930de9b25SBarry Smith 
16030de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
16130de9b25SBarry Smith 
16230de9b25SBarry Smith M*/
163273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
164273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
165273d9f13SBarry Smith #else
166273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
16730de9b25SBarry Smith #endif
16830de9b25SBarry Smith 
169273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
170b0a32e0cSBarry Smith extern PetscFList MatList;
17128988994SBarry Smith 
172be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
173be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
174be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
175ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
176ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
177ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
178ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
179ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
180ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
181ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
182be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
183ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(comm,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
184ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(comm,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
185ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(comm,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
186ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(comm,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
187ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(comm,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
188ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(comm,m,n,M,N,0,nnz,0,onz,A))
189ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(comm,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
190ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(PETSC_COMM_WORLD,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
191ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(PETSC_COMM_WORLD,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
192ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(PETSC_COMM_WORLD,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
193ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N),(PETSC_COMM_WORLD,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
194ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(PETSC_COMM_WORLD,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
195ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(PETSC_COMM_WORLD,m,n,M,N,0,nnz,0,onz,A))
196ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(PETSC_COMM_WORLD,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
197be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
198be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*);
199be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*);
200be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
201ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,nz,nnz,&A),Mat,A)
202ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
203ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
204ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
205ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
206ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
207ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
208be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
209ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(comm,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
210ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
211ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(comm,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
212ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
213ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
214ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(comm,bs,m,n,M,N,0,nnz,0,onz,A))
215ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
216ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
217ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
218ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
219ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
220ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
221ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,A))
222ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
223be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
224be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
225ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,nz,nnz,&A),Mat,A)
226ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
227ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
228ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
229ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
230ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
231ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
232be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
233ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(comm,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
234ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
235ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(comm,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
236ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
237ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
238ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(comm,bs,m,n,M,N,0,nnz,0,onz,A))
239ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
240ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
241ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
242ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
243ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
244ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
245ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,A))
246ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
247be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
248ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(comm,m,n,M,N,ctx,&A),Mat,A)
249ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
250be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
251be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
252ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
253ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*);
254f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
255be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
25621c89e3eSBarry Smith 
257be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPrintHelp(Mat);
258be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetPetscMaps(Mat,PetscMap*,PetscMap*);
259ec0117caSBarry Smith 
2600c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
2618ed539a5SBarry Smith /* ------------------------------------------------------------*/
262be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
263be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
26484cb2905SBarry Smith 
2652ef4de8bSBarry Smith /*S
2662ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
2672ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
2682ef4de8bSBarry Smith 
2692ef4de8bSBarry Smith    Level: beginner
2702ef4de8bSBarry Smith 
2712ef4de8bSBarry Smith   Concepts: matrix; linear operator
2722ef4de8bSBarry Smith 
273d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
2742ef4de8bSBarry Smith S*/
275435da068SBarry Smith typedef struct {
276c1ac3661SBarry Smith   PetscInt k,j,i,c;
277435da068SBarry Smith } MatStencil;
2782ef4de8bSBarry Smith 
279be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
280be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
281be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
282435da068SBarry Smith 
283be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
284be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
285be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
2863a7fca6bSBarry Smith 
287d91e6319SBarry Smith /*E
288d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
289d91e6319SBarry Smith      to continue to add values to it
290d91e6319SBarry Smith 
291d91e6319SBarry Smith     Level: beginner
292d91e6319SBarry Smith 
293d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
294d91e6319SBarry Smith E*/
2956d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
297be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
298be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
2994f9c727eSBarry Smith 
300be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt    MatSetValue_Row;
301be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt    MatSetValue_Column;
302be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscScalar MatSetValue_Value;
3031d73ed98SBarry Smith 
30430de9b25SBarry Smith /*MC
30530de9b25SBarry Smith    MatSetValue - Set a single entry into a matrix.
30630de9b25SBarry Smith 
30730de9b25SBarry Smith    Synopsis:
308c1ac3661SBarry Smith    PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode);
30930de9b25SBarry Smith 
31030de9b25SBarry Smith    Not collective
31130de9b25SBarry Smith 
31230de9b25SBarry Smith    Input Parameters:
31330de9b25SBarry Smith +  m - the matrix
31430de9b25SBarry Smith .  row - the row location of the entry
31530de9b25SBarry Smith .  col - the column location of the entry
31630de9b25SBarry Smith .  value - the value to insert
31730de9b25SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
31830de9b25SBarry Smith 
31930de9b25SBarry Smith    Notes:
32030de9b25SBarry Smith    For efficiency one should use MatSetValues() and set several or many
32130de9b25SBarry Smith    values simultaneously if possible.
32230de9b25SBarry Smith 
32330de9b25SBarry Smith    Level: beginner
32430de9b25SBarry Smith 
32530de9b25SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
32630de9b25SBarry Smith M*/
327b951964fSBarry Smith #define MatSetValue(v,i,j,va,mode) \
328f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \
3291d73ed98SBarry Smith    MatSetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode))
33030de9b25SBarry Smith 
331ea06a074SBarry Smith #define MatGetValue(v,i,j,va) \
332f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,0) || \
3331d73ed98SBarry Smith    MatGetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&va))
33430de9b25SBarry Smith 
335d91e6319SBarry Smith #define MatSetValueLocal(v,i,j,va,mode) \
336f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \
3371d73ed98SBarry Smith    MatSetValuesLocal(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode))
33830de9b25SBarry Smith 
339d91e6319SBarry Smith /*E
340d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
341d91e6319SBarry Smith 
342d91e6319SBarry Smith     Level: beginner
343d91e6319SBarry Smith 
3440a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
345d91e6319SBarry Smith 
346d91e6319SBarry Smith .seealso: MatSetOption()
347d91e6319SBarry Smith E*/
3486d4a8577SBarry Smith typedef enum {MAT_ROW_ORIENTED=1,MAT_COLUMN_ORIENTED=2,MAT_ROWS_SORTED=4,
3496d4a8577SBarry Smith               MAT_COLUMNS_SORTED=8,MAT_NO_NEW_NONZERO_LOCATIONS=16,
3506d4a8577SBarry Smith               MAT_YES_NEW_NONZERO_LOCATIONS=32,MAT_SYMMETRIC=64,
3516ca9ecd3SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC=65,MAT_NO_NEW_DIAGONALS=66,
3526ca9ecd3SBarry Smith               MAT_YES_NEW_DIAGONALS=67,MAT_INODE_LIMIT_1=68,MAT_INODE_LIMIT_2=69,
3536ca9ecd3SBarry Smith               MAT_INODE_LIMIT_3=70,MAT_INODE_LIMIT_4=71,MAT_INODE_LIMIT_5=72,
3546ca9ecd3SBarry Smith               MAT_IGNORE_OFF_PROC_ENTRIES=73,MAT_ROWS_UNSORTED=74,
3554787f768SSatish Balay               MAT_COLUMNS_UNSORTED=75,MAT_NEW_NONZERO_LOCATION_ERR=76,
3567c922b88SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR=77,MAT_USE_HASH_TABLE=78,
3572bad1931SBarry Smith               MAT_KEEP_ZEROED_ROWS=79,MAT_IGNORE_ZERO_ENTRIES=80,MAT_USE_INODES=81,
358efcf0fc3SBarry Smith               MAT_DO_NOT_USE_INODES=82,MAT_NOT_SYMMETRIC=83,MAT_HERMITIAN=84,
3599a4540c5SBarry Smith               MAT_NOT_STRUCTURALLY_SYMMETRIC=85,MAT_NOT_HERMITIAN=86,
360d487561eSHong Zhang               MAT_SYMMETRY_ETERNAL=87,MAT_NOT_SYMMETRY_ETERNAL=88,
361941593c8SHong Zhang               MAT_USE_COMPRESSEDROW=89,MAT_DO_NOT_USE_COMPRESSEDROW=90,
362941593c8SHong Zhang               MAT_IGNORE_LOWER_TRIANGULAR=91,MAT_ERROR_LOWER_TRIANGULAR=92} MatOption;
363be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption);
364be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,MatType*);
365ba966639SSatish Balay PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),MatType,t)
36684cb2905SBarry Smith 
367be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
368be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
370be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
373be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
374ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
375be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
376be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
377ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
378be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
3797b80b807SBarry Smith 
380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
382ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
384be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
385ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
386ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0.0,&t),PetscTruth,t)
387be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
388ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
389be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
390be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
3912eac72dbSBarry Smith 
392d91e6319SBarry Smith /*E
393d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
394d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
395d91e6319SBarry Smith 
396d91e6319SBarry Smith     Level: beginner
397d91e6319SBarry Smith 
398d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
399d91e6319SBarry Smith 
400d91e6319SBarry Smith .seealso: MatDuplicate()
401d91e6319SBarry Smith E*/
4022e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
4032e8a6d31SBarry Smith 
404be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvertRegister(const char[],const char[],const char[],PetscErrorCode (*)(Mat,MatType,MatReuse,Mat*));
405273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
406273d9f13SBarry Smith #define MatConvertRegisterDynamic(a,b,c,d) MatConvertRegister(a,b,c,0)
407273d9f13SBarry Smith #else
408273d9f13SBarry Smith #define MatConvertRegisterDynamic(a,b,c,d) MatConvertRegister(a,b,c,d)
409273d9f13SBarry Smith #endif
410be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvertRegisterAll(const char[]);
411be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvertRegisterDestroy(void);
412273d9f13SBarry Smith extern PetscTruth MatConvertRegisterAllCalled;
413b0a32e0cSBarry Smith extern PetscFList MatConvertList;
414f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,MatType,MatReuse,Mat*);
415ba966639SSatish Balay PetscPolymorphicFunction(MatConvert,(Mat A,MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
416be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
417ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
418ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
41994a9d846SBarry Smith 
420d91e6319SBarry Smith /*E
421d91e6319SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
422d91e6319SBarry Smith 
423d91e6319SBarry Smith     Level: beginner
424d91e6319SBarry Smith 
425d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
426d91e6319SBarry Smith 
42794b7f48cSBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
428d91e6319SBarry Smith E*/
429c537a176SHong Zhang typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
430cb5b572fSBarry Smith 
431be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
432be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
433be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
434ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
435ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0.0,&t),PetscTruth,t)
436be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
437ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
438be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscTruth*);
439ba966639SSatish Balay PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,&t),PetscTruth,t)
440be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
441be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
442f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,MatType,Mat*);
443ba966639SSatish Balay PetscPolymorphicFunction(MatLoad,(PetscViewer v,MatType t),(v,t,&a),Mat,a)
4447b80b807SBarry Smith 
445be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
446be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
447be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
449d4fbbf0eSBarry Smith 
450d91e6319SBarry Smith /*S
451d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
452d91e6319SBarry Smith 
453d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
454d91e6319SBarry Smith 
455d91e6319SBarry Smith    Level: intermediate
456d91e6319SBarry Smith 
457d91e6319SBarry Smith   Concepts: matrix^nonzero information
458d91e6319SBarry Smith 
459d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
460d91e6319SBarry Smith S*/
4614e220ebcSLois Curfman McInnes typedef struct {
462b0a32e0cSBarry Smith   PetscLogDouble rows_global,columns_global;         /* number of global rows and columns */
463b0a32e0cSBarry Smith   PetscLogDouble rows_local,columns_local;           /* number of local rows and columns */
464b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
465b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
466b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
467b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
468b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
469b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
470b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4714e220ebcSLois Curfman McInnes } MatInfo;
4724e220ebcSLois Curfman McInnes 
473d9274352SBarry Smith /*E
474d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
475d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
476d9274352SBarry Smith 
477d9274352SBarry Smith     Level: beginner
478d9274352SBarry Smith 
479d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
480d9274352SBarry Smith 
481d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
482d9274352SBarry Smith E*/
4837b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
484be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
485be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
486be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
487be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec);
488be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,Mat*);
489ba966639SSatish Balay PetscPolymorphicFunction(MatTranspose,(Mat A),(A,&t),Mat,t)
490be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
491ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
492be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
493be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
494be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
495be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
496ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
497be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
499be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
5017b80b807SBarry Smith 
502be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
503ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
504be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
505f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
506f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
507f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
508f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5097b80b807SBarry Smith 
510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5135ef9f2a5SBarry Smith 
514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
515be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5177b80b807SBarry Smith 
518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
519be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
520be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *);
521*abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *);
522be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
523be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
524be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
525be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
526be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
529be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscScalar**,Mat*);
53143eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
532cd116e26SSatish Balay #include "petscctable.h"
53343eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
53443eb5e2fSMatthew Knepley #else
53543eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
53643eb5e2fSMatthew Knepley #endif
5378efafbd8SBarry Smith 
538be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5397b80b807SBarry Smith 
540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
542be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
54322440eb1SKris Buschelman 
544be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
54722440eb1SKris Buschelman 
548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
551bc011b1eSHong Zhang 
552f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
553f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat);
554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat);
5551c741599SHong Zhang 
556f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
557f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5587b80b807SBarry Smith 
559be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
561f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
562f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
563be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
565052efed2SBarry Smith 
566be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
56890f02eecSBarry Smith 
569be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
5700c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
571ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
572be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
573be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
574bd481603SBarry Smith 
575bd481603SBarry Smith /*MC
576bd481603SBarry Smith    MatPreallocInitialize - Begins the block of code that will count the number of nonzeros per
577bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
578bd481603SBarry Smith 
579bd481603SBarry Smith    Synopsis:
580c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
581bd481603SBarry Smith 
582bd481603SBarry Smith    Collective on MPI_Comm
583bd481603SBarry Smith 
584bd481603SBarry Smith    Input Parameters:
585bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
586859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
587859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
588bd481603SBarry Smith 
589bd481603SBarry Smith    Output Parameters:
590bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
591bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
592bd481603SBarry Smith 
593bd481603SBarry Smith 
594bd481603SBarry Smith    Level: intermediate
595bd481603SBarry Smith 
596bd481603SBarry Smith    Notes:
597bd481603SBarry Smith    See the chapter in the users manual on performance for more details
598bd481603SBarry Smith 
5991d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
600bd481603SBarry Smith 
601bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
602bd481603SBarry Smith 
603bd481603SBarry Smith   Concepts: preallocation^Matrix
604bd481603SBarry Smith 
605bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
606bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
607bd481603SBarry Smith M*/
608c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6097c922b88SBarry Smith { \
610c1ac3661SBarry Smith   PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
611c1ac3661SBarry Smith   _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\
612c1ac3661SBarry Smith   _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\
613c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
614c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp;
6157c922b88SBarry Smith 
616bd481603SBarry Smith /*MC
617bd481603SBarry Smith    MatPreallocSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
618bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
619bd481603SBarry Smith 
620bd481603SBarry Smith    Synopsis:
621c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
622bd481603SBarry Smith 
623bd481603SBarry Smith    Collective on MPI_Comm
624bd481603SBarry Smith 
625bd481603SBarry Smith    Input Parameters:
626bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
627859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
628859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
629bd481603SBarry Smith 
630bd481603SBarry Smith    Output Parameters:
631bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
632bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
633bd481603SBarry Smith 
634bd481603SBarry Smith 
635bd481603SBarry Smith    Level: intermediate
636bd481603SBarry Smith 
637bd481603SBarry Smith    Notes:
638bd481603SBarry Smith    See the chapter in the users manual on performance for more details
639bd481603SBarry Smith 
6401d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
641bd481603SBarry Smith 
642bd481603SBarry Smith   Concepts: preallocation^Matrix
643bd481603SBarry Smith 
644bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
645bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
646bd481603SBarry Smith M*/
647222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
648222b16d4SBarry Smith { \
649c1ac3661SBarry Smith   PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__end; \
650c1ac3661SBarry Smith   _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\
651c1ac3661SBarry Smith   _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\
652c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
653c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp;
654222b16d4SBarry Smith 
655bd481603SBarry Smith /*MC
656bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
657bd481603SBarry Smith        inserted using a local number of the rows and columns
658bd481603SBarry Smith 
659bd481603SBarry Smith    Synopsis:
660c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
661bd481603SBarry Smith 
662bd481603SBarry Smith    Not Collective
663bd481603SBarry Smith 
664bd481603SBarry Smith    Input Parameters:
665bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
666bd481603SBarry Smith .  nrows - the number of rows indicated
6671d73ed98SBarry Smith .  rows - the indices of the rows
668bd481603SBarry Smith .  ncols - the number of columns in the matrix
669bd481603SBarry Smith .  cols - the columns indicated
670bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
671bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
672bd481603SBarry Smith 
673bd481603SBarry Smith 
674bd481603SBarry Smith    Level: intermediate
675bd481603SBarry Smith 
676bd481603SBarry Smith    Notes:
677bd481603SBarry Smith    See the chapter in the users manual on performance for more details
678bd481603SBarry Smith 
6791d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
680bd481603SBarry Smith 
681bd481603SBarry Smith   Concepts: preallocation^Matrix
682bd481603SBarry Smith 
683bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
684bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
685bd481603SBarry Smith M*/
686c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
687c4f061fbSSatish Balay {\
688c1ac3661SBarry Smith   PetscInt __l;\
689ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
690ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
691c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
692ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
693c4f061fbSSatish Balay   }\
694c4f061fbSSatish Balay }
695c4f061fbSSatish Balay 
696bd481603SBarry Smith /*MC
697bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
698bd481603SBarry Smith        inserted using a local number of the rows and columns
699bd481603SBarry Smith 
700bd481603SBarry Smith    Synopsis:
701c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
702bd481603SBarry Smith 
703bd481603SBarry Smith    Not Collective
704bd481603SBarry Smith 
705bd481603SBarry Smith    Input Parameters:
706bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
707bd481603SBarry Smith .  nrows - the number of rows indicated
7081d73ed98SBarry Smith .  rows - the indices of the rows
709bd481603SBarry Smith .  ncols - the number of columns in the matrix
710bd481603SBarry Smith .  cols - the columns indicated
711bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
712bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
713bd481603SBarry Smith 
714bd481603SBarry Smith 
715bd481603SBarry Smith    Level: intermediate
716bd481603SBarry Smith 
717bd481603SBarry Smith    Notes:
718bd481603SBarry Smith    See the chapter in the users manual on performance for more details
719bd481603SBarry Smith 
720bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
721bd481603SBarry Smith 
722bd481603SBarry Smith   Concepts: preallocation^Matrix
723bd481603SBarry Smith 
724bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
725bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
726bd481603SBarry Smith M*/
727d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
728d3d32019SBarry Smith {\
729c1ac3661SBarry Smith   PetscInt __l;\
730d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
731d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
732d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
733d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
734d3d32019SBarry Smith   }\
735d3d32019SBarry Smith }
736d3d32019SBarry Smith 
737bd481603SBarry Smith /*MC
738bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
739bd481603SBarry Smith        inserted using a local number of the rows and columns
740bd481603SBarry Smith 
741bd481603SBarry Smith    Synopsis:
742c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
743bd481603SBarry Smith 
744bd481603SBarry Smith    Not Collective
745bd481603SBarry Smith 
746bd481603SBarry Smith    Input Parameters:
747bd481603SBarry Smith +  nrows - the number of rows indicated
7481d73ed98SBarry Smith .  rows - the indices of the rows
749bd481603SBarry Smith .  ncols - the number of columns in the matrix
750a50f8bf6SHong Zhang -  cols - the columns indicated
751a50f8bf6SHong Zhang 
752a50f8bf6SHong Zhang    Output Parameters:
753a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
754bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
755bd481603SBarry Smith 
756bd481603SBarry Smith 
757bd481603SBarry Smith    Level: intermediate
758bd481603SBarry Smith 
759bd481603SBarry Smith    Notes:
760bd481603SBarry Smith    See the chapter in the users manual on performance for more details
761bd481603SBarry Smith 
762bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
763bd481603SBarry Smith 
764bd481603SBarry Smith   Concepts: preallocation^Matrix
765bd481603SBarry Smith 
766bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
767bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
768bd481603SBarry Smith M*/
769c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
770c1ac3661SBarry Smith { PetscInt __i; \
7717c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
7727c922b88SBarry Smith     if (cols[__i] < __start || cols[__i] >= __end) onz[row - __rstart]++; \
7737c922b88SBarry Smith   }\
7747c922b88SBarry Smith   dnz[row - __rstart] = nc - onz[row - __rstart];\
7757c922b88SBarry Smith }
7767c922b88SBarry Smith 
777bd481603SBarry Smith /*MC
778bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
779bd481603SBarry Smith        inserted using a local number of the rows and columns
780bd481603SBarry Smith 
781bd481603SBarry Smith    Synopsis:
782c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
783bd481603SBarry Smith 
784bd481603SBarry Smith    Not Collective
785bd481603SBarry Smith 
786bd481603SBarry Smith    Input Parameters:
787bd481603SBarry Smith +  nrows - the number of rows indicated
7881d73ed98SBarry Smith .  rows - the indices of the rows
789bd481603SBarry Smith .  ncols - the number of columns in the matrix
790bd481603SBarry Smith .  cols - the columns indicated
791bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
792bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
793bd481603SBarry Smith 
794bd481603SBarry Smith 
795bd481603SBarry Smith    Level: intermediate
796bd481603SBarry Smith 
797bd481603SBarry Smith    Notes:
798bd481603SBarry Smith    See the chapter in the users manual on performance for more details
799bd481603SBarry Smith 
800bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
801bd481603SBarry Smith 
802bd481603SBarry Smith   Concepts: preallocation^Matrix
803bd481603SBarry Smith 
804bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
805bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
806bd481603SBarry Smith M*/
807d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
808c1ac3661SBarry Smith { PetscInt __i; \
809d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
810d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
811d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
812d3d32019SBarry Smith   }\
813d3d32019SBarry Smith }
814d3d32019SBarry Smith 
815bd481603SBarry Smith /*MC
816bd481603SBarry Smith    MatPreallocFinalize - Ends the block of code that will count the number of nonzeros per
817bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
818bd481603SBarry Smith 
819bd481603SBarry Smith    Synopsis:
820c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
821bd481603SBarry Smith 
822bd481603SBarry Smith    Collective on MPI_Comm
823bd481603SBarry Smith 
824bd481603SBarry Smith    Input Parameters:
825bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
826bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
827bd481603SBarry Smith 
828bd481603SBarry Smith 
829bd481603SBarry Smith    Level: intermediate
830bd481603SBarry Smith 
831bd481603SBarry Smith    Notes:
832bd481603SBarry Smith    See the chapter in the users manual on performance for more details
833bd481603SBarry Smith 
834bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
835bd481603SBarry Smith 
836bd481603SBarry Smith   Concepts: preallocation^Matrix
837bd481603SBarry Smith 
838bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
839bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
840bd481603SBarry Smith M*/
841ef66eb69SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree(dnz);CHKERRQ(_4_ierr);}
8427c922b88SBarry Smith 
843bd481603SBarry Smith 
844bd481603SBarry Smith 
8457b80b807SBarry Smith /* Routines unique to particular data structures */
846be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
847ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
848be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBDiagGetData(Mat,PetscInt*,PetscInt*,PetscInt*[],PetscInt*[],PetscScalar***);
849698d4c6aSKris Buschelman 
850be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
851be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
852698d4c6aSKris Buschelman 
853be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
854be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
855be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
8567b80b807SBarry Smith 
857a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
858a96a251dSBarry Smith 
859be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
860ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
861be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
862ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
863be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
864ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
865be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDensePreallocation(Mat,PetscScalar[]);
866be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]);
867be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
868273d9f13SBarry Smith 
869be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
870ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
871be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
872be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
873be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
874be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
875be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDensePreallocation(Mat,PetscScalar[]);
876be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]);
877be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
878be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
879be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
880be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
881be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
882be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
883273d9f13SBarry Smith 
884be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
885ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*);
8861b807ce4Svictorle 
887be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
888be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
8892e8a6d31SBarry Smith 
890be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
8913a7fca6bSBarry Smith 
8927b80b807SBarry Smith /*
8937b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
89494b7f48cSBarry Smith   done through the KSP and PC interfaces.
8957b80b807SBarry Smith */
8967b80b807SBarry Smith 
897d9274352SBarry Smith /*E
898d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
899d9274352SBarry Smith        with an optional dynamic library name, for example
900d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
901d9274352SBarry Smith 
902d9274352SBarry Smith    Level: beginner
903d9274352SBarry Smith 
904d9274352SBarry Smith .seealso: MatGetOrdering()
905d9274352SBarry Smith E*/
90649773a63SBarry Smith #define MatOrderingType char*
907b12f92e5SBarry Smith #define MATORDERING_NATURAL   "natural"
908b12f92e5SBarry Smith #define MATORDERING_ND        "nd"
909b12f92e5SBarry Smith #define MATORDERING_1WD       "1wd"
910b12f92e5SBarry Smith #define MATORDERING_RCM       "rcm"
911b12f92e5SBarry Smith #define MATORDERING_QMD       "qmd"
912b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength"
91362152c8bSBarry Smith #define MATORDERING_DSC_ND    "dsc_nd"
91462152c8bSBarry Smith #define MATORDERING_DSC_MMD   "dsc_mmd"
91562152c8bSBarry Smith #define MATORDERING_DSC_MDF   "dsc_mdf"
916c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
917c06d978dSMatthew Knepley #define MATORDERING_IDENTITY  "identity"
918c06d978dSMatthew Knepley #define MATORDERING_REVERSE   "reverse"
919b12f92e5SBarry Smith 
920be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
921be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
922be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
92330de9b25SBarry Smith 
92430de9b25SBarry Smith /*MC
92530de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
92630de9b25SBarry Smith                                matrix package.
92730de9b25SBarry Smith 
92830de9b25SBarry Smith    Synopsis:
929c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
93030de9b25SBarry Smith 
93130de9b25SBarry Smith    Not Collective
93230de9b25SBarry Smith 
93330de9b25SBarry Smith    Input Parameters:
93430de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
93530de9b25SBarry Smith .  path - location of library where creation routine is
93630de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
93730de9b25SBarry Smith -  function - function pointer that creates the ordering
93830de9b25SBarry Smith 
93930de9b25SBarry Smith    Level: developer
94030de9b25SBarry Smith 
94130de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
94230de9b25SBarry Smith    is ignored.
94330de9b25SBarry Smith 
94430de9b25SBarry Smith    Sample usage:
94530de9b25SBarry Smith .vb
94630de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
94730de9b25SBarry Smith                "MyOrder",MyOrder);
94830de9b25SBarry Smith .ve
94930de9b25SBarry Smith 
95030de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
95130de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
95230de9b25SBarry Smith    or at runtime via the option
95330de9b25SBarry Smith $     -pc_ilu_mat_ordering_type my_order
95430de9b25SBarry Smith $     -pc_lu_mat_ordering_type my_order
95530de9b25SBarry Smith 
956ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
95730de9b25SBarry Smith 
95830de9b25SBarry Smith .keywords: matrix, ordering, register
95930de9b25SBarry Smith 
96030de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
96130de9b25SBarry Smith M*/
962aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
963f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
964b12f92e5SBarry Smith #else
965f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
966b12f92e5SBarry Smith #endif
96730de9b25SBarry Smith 
968be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
969be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
9702bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
971b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
972d4fbbf0eSBarry Smith 
973be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
974a2ce50c7SBarry Smith 
975d91e6319SBarry Smith /*S
97615e8a5b3SHong Zhang    MatFactorInfo - Data based into the matrix factorization routines
977b00f7748SHong Zhang 
97815e8a5b3SHong Zhang    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE
979b00f7748SHong Zhang 
98015e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
981b00f7748SHong Zhang 
982b00f7748SHong Zhang    Level: developer
983b00f7748SHong Zhang 
984d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
985d7d82daaSBarry Smith           MatFactorInfoInitialize()
986b00f7748SHong Zhang 
987b00f7748SHong Zhang S*/
988b00f7748SHong Zhang typedef struct {
9890a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
990fbf22428SSatish Balay   PetscReal     shiftpd;        /* if true, shift until positive pivots */
9912cea7109SBarry Smith   PetscReal     shift_fraction; /* record shift fraction taken */
99215e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
99315e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
994b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
99515e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
996f6275e2eSBarry Smith   PetscReal     fill;           /* expected fill; nonzeros in factored matrix/nonzeros in original matrix*/
997348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
998bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
999bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
100015e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
100115e8a5b3SHong Zhang } MatFactorInfo;
1002ffa6d0a5SLois Curfman McInnes 
1003be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
1004be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,MatFactorInfo*);
1005be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1006be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,MatFactorInfo*,Mat*);
1007be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,MatFactorInfo*);
1008be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,MatFactorInfo*);
1009be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1010be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1011be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1012be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,MatFactorInfo*);
1013be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,MatFactorInfo*,Mat*);
1014be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,MatFactorInfo*,Mat *);
1015be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1016be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1017be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1018be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1019be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1020be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1021be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1022be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
10238ed539a5SBarry Smith 
1024be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1025bb5a7306SBarry Smith 
1026d91e6319SBarry Smith /*E
1027d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1028bb1eb677SSatish Balay 
1029d91e6319SBarry Smith     Level: beginner
1030d91e6319SBarry Smith 
1031d9274352SBarry Smith    May be bitwise ORd together
1032d9274352SBarry Smith 
1033d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1034d91e6319SBarry Smith 
10354e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10364e7234bfSBarry Smith 
1037a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax()
1038d91e6319SBarry Smith E*/
1039ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1040ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1041ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
104284cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1043be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1044be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10458ed539a5SBarry Smith 
1046d4fbbf0eSBarry Smith /*
1047639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1048639f9d9dSBarry Smith */
1049b12f92e5SBarry Smith 
1050d9274352SBarry Smith /*E
1051d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1052d9274352SBarry Smith        with an optional dynamic library name, for example
1053d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1054d9274352SBarry Smith 
1055d9274352SBarry Smith    Level: beginner
1056d9274352SBarry Smith 
1057d9274352SBarry Smith .seealso: MatGetColoring()
1058d9274352SBarry Smith E*/
105949773a63SBarry Smith #define MatColoringType char*
1060b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1061b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1062b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1063b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1064b12f92e5SBarry Smith 
1065be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*);
1066be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatColoringType,ISColoring *));
106730de9b25SBarry Smith 
106830de9b25SBarry Smith /*MC
106930de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
107030de9b25SBarry Smith                                matrix package.
107130de9b25SBarry Smith 
107230de9b25SBarry Smith    Synopsis:
1073c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
107430de9b25SBarry Smith 
107530de9b25SBarry Smith    Not Collective
107630de9b25SBarry Smith 
107730de9b25SBarry Smith    Input Parameters:
107830de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
107930de9b25SBarry Smith .  path - location of library where creation routine is
108030de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
108130de9b25SBarry Smith -  function - function pointer that creates the coloring
108230de9b25SBarry Smith 
108330de9b25SBarry Smith    Level: developer
108430de9b25SBarry Smith 
108530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
108630de9b25SBarry Smith    is ignored.
108730de9b25SBarry Smith 
108830de9b25SBarry Smith    Sample usage:
108930de9b25SBarry Smith .vb
109030de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
109130de9b25SBarry Smith                "MyColor",MyColor);
109230de9b25SBarry Smith .ve
109330de9b25SBarry Smith 
109430de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
109530de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
109630de9b25SBarry Smith    or at runtime via the option
109730de9b25SBarry Smith $     -mat_coloring_type my_color
109830de9b25SBarry Smith 
1099ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
110030de9b25SBarry Smith 
110130de9b25SBarry Smith .keywords: matrix, Coloring, register
110230de9b25SBarry Smith 
110330de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
110430de9b25SBarry Smith M*/
1105aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1106f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1107b12f92e5SBarry Smith #else
1108f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1109b12f92e5SBarry Smith #endif
111030de9b25SBarry Smith 
11112bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1112f1144a30SSatish Balay 
1113be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1114be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1115be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1116639f9d9dSBarry Smith 
1117d9274352SBarry Smith /*S
1118d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1119d9274352SBarry Smith         and coloring
1120639f9d9dSBarry Smith 
1121d9274352SBarry Smith    Level: beginner
1122d9274352SBarry Smith 
1123d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1124d9274352SBarry Smith 
1125d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1126d9274352SBarry Smith S*/
1127e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1128639f9d9dSBarry Smith 
1129be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1130be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1131be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1132be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1133be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1134be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring,PetscInt);
1135be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring,PetscInt*);
1136be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1137be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1138be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1139be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring);
1140be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1141be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1142639f9d9dSBarry Smith /*
11430752156aSBarry Smith     These routines are for partitioning matrices: currently used only
11443eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
11450752156aSBarry Smith */
1146ca161407SBarry Smith 
1147d9274352SBarry Smith /*S
1148d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1149d9274352SBarry Smith 
1150d9274352SBarry Smith    Level: beginner
1151d9274352SBarry Smith 
1152d9274352SBarry Smith   Concepts: partitioning
1153d9274352SBarry Smith 
1154743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1155d9274352SBarry Smith S*/
115691e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1157d9274352SBarry Smith 
1158d9274352SBarry Smith /*E
11595bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1160d9274352SBarry Smith        with an optional dynamic library name, for example
1161d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1162d9274352SBarry Smith 
1163d9274352SBarry Smith    Level: beginner
1164d9274352SBarry Smith 
1165b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1166d9274352SBarry Smith E*/
116749773a63SBarry Smith #define MatPartitioningType char*
11688ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1169c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
11708ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
117171306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
117271306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
117371306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
117471306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
117571306c60Sjroman 
1176ca161407SBarry Smith 
1177be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1178be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1179be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1180be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1181be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1182be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1183be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1184be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
11852aabb6bbSBarry Smith 
1186be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
118730de9b25SBarry Smith 
118830de9b25SBarry Smith /*MC
118930de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
119030de9b25SBarry Smith    matrix package.
119130de9b25SBarry Smith 
119230de9b25SBarry Smith    Synopsis:
1193c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
119430de9b25SBarry Smith 
119530de9b25SBarry Smith    Not Collective
119630de9b25SBarry Smith 
119730de9b25SBarry Smith    Input Parameters:
119830de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
119930de9b25SBarry Smith .  path - location of library where creation routine is
120030de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
120130de9b25SBarry Smith -  function - function pointer that creates the partitioning type
120230de9b25SBarry Smith 
120330de9b25SBarry Smith    Level: developer
120430de9b25SBarry Smith 
120530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
120630de9b25SBarry Smith    is ignored.
120730de9b25SBarry Smith 
120830de9b25SBarry Smith    Sample usage:
120930de9b25SBarry Smith .vb
121030de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
121130de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
121230de9b25SBarry Smith .ve
121330de9b25SBarry Smith 
121430de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
121530de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
121630de9b25SBarry Smith    or at runtime via the option
121730de9b25SBarry Smith $     -mat_partitioning_type my_part
121830de9b25SBarry Smith 
1219ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
122030de9b25SBarry Smith 
122130de9b25SBarry Smith .keywords: matrix, partitioning, register
122230de9b25SBarry Smith 
122330de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
122430de9b25SBarry Smith M*/
1225aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1226f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
12272aabb6bbSBarry Smith #else
1228f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
12292aabb6bbSBarry Smith #endif
12302aabb6bbSBarry Smith 
12312bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1232f1144a30SSatish Balay 
1233be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1234be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
12352bad1931SBarry Smith 
1236be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1237be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1238be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1239ca161407SBarry Smith 
1240be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
12410752156aSBarry Smith 
1242be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1243be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
124471306c60Sjroman 
124571306c60Sjroman typedef enum { MP_CHACO_MULTILEVEL_KL, MP_CHACO_SPECTRAL, MP_CHACO_LINEAR,
124671306c60Sjroman     MP_CHACO_RANDOM, MP_CHACO_SCATTERED } MPChacoGlobalType;
1247be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
124871306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1249be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1250be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
125171306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1252be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1253be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1254be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
125571306c60Sjroman 
125671306c60Sjroman #define MP_PARTY_OPT "opt"
125771306c60Sjroman #define MP_PARTY_LIN "lin"
125871306c60Sjroman #define MP_PARTY_SCA "sca"
125971306c60Sjroman #define MP_PARTY_RAN "ran"
126071306c60Sjroman #define MP_PARTY_GBF "gbf"
126171306c60Sjroman #define MP_PARTY_GCF "gcf"
126271306c60Sjroman #define MP_PARTY_BUB "bub"
126371306c60Sjroman #define MP_PARTY_DEF "def"
1264be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
126571306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
126671306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
126771306c60Sjroman #define MP_PARTY_NONE "no"
1268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1269be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1270be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1271be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
127271306c60Sjroman 
127371306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1274be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1275be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1277be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1278be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
127971306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1280be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1281be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1282be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
128371306c60Sjroman 
12840752156aSBarry Smith /*
12850a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1286d4fbbf0eSBarry Smith */
12871c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
12881c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
12891c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
12901c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
12911c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
12927c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
12937c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
12941c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
12951c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
12967c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
12977c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
12981c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
12991c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
13001c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
13011c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13021c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13031c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13041c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13051c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13061c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13071c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13081c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
13091c1c02c0SLois Curfman McInnes                MATOP_COMPRESS=22,
13101c1c02c0SLois Curfman McInnes                MATOP_SET_OPTION=23,
13111c1c02c0SLois Curfman McInnes                MATOP_ZERO_ENTRIES=24,
13121c1c02c0SLois Curfman McInnes                MATOP_ZERO_ROWS=25,
13131c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_SYMBOLIC=26,
13141c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_NUMERIC=27,
13151c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_SYMBOLIC=28,
13161c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_NUMERIC=29,
1317d643ce63SMatthew Knepley                MATOP_SETUP_PREALLOCATION=30,
1318d643ce63SMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=31,
1319d643ce63SMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=32,
1320d643ce63SMatthew Knepley                MATOP_GET_ARRAY=33,
1321d643ce63SMatthew Knepley                MATOP_RESTORE_ARRAY=34,
1322e26f1c3cSBarry Smith                MATOP_DUPLICATE=35,
1323d643ce63SMatthew Knepley                MATOP_FORWARD_SOLVE=36,
1324d643ce63SMatthew Knepley                MATOP_BACKWARD_SOLVE=37,
1325d643ce63SMatthew Knepley                MATOP_ILUFACTOR=38,
1326d643ce63SMatthew Knepley                MATOP_ICCFACTOR=39,
1327d643ce63SMatthew Knepley                MATOP_AXPY=40,
1328d643ce63SMatthew Knepley                MATOP_GET_SUBMATRICES=41,
1329d643ce63SMatthew Knepley                MATOP_INCREASE_OVERLAP=42,
1330d643ce63SMatthew Knepley                MATOP_GET_VALUES=43,
1331d643ce63SMatthew Knepley                MATOP_COPY=44,
1332d643ce63SMatthew Knepley                MATOP_PRINT_HELP=45,
1333d643ce63SMatthew Knepley                MATOP_SCALE=46,
1334d643ce63SMatthew Knepley                MATOP_SHIFT=47,
1335d643ce63SMatthew Knepley                MATOP_DIAGONAL_SHIFT=48,
1336d643ce63SMatthew Knepley                MATOP_ILUDT_FACTOR=49,
1337d643ce63SMatthew Knepley                MATOP_GET_BLOCK_SIZE=50,
1338d643ce63SMatthew Knepley                MATOP_GET_ROW_IJ=51,
1339d643ce63SMatthew Knepley                MATOP_RESTORE_ROW_IJ=52,
1340d643ce63SMatthew Knepley                MATOP_GET_COLUMN_IJ=53,
1341d643ce63SMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=54,
1342d643ce63SMatthew Knepley                MATOP_FDCOLORING_CREATE=55,
1343d643ce63SMatthew Knepley                MATOP_COLORING_PATCH=56,
1344d643ce63SMatthew Knepley                MATOP_SET_UNFACTORED=57,
1345d643ce63SMatthew Knepley                MATOP_PERMUTE=58,
1346d643ce63SMatthew Knepley                MATOP_SET_VALUES_BLOCKED=59,
1347d643ce63SMatthew Knepley                MATOP_GET_SUBMATRIX=60,
1348d643ce63SMatthew Knepley                MATOP_DESTROY=61,
1349d643ce63SMatthew Knepley                MATOP_VIEW=62,
1350d643ce63SMatthew Knepley                MATOP_GET_MAPS=63,
1351d643ce63SMatthew Knepley                MATOP_USE_SCALED_FORM=64,
1352d643ce63SMatthew Knepley                MATOP_SCALE_SYSTEM=65,
1353d643ce63SMatthew Knepley                MATOP_UNSCALE_SYSTEM=66,
1354ba16a8cbSKris Buschelman                MATOP_SET_LOCAL_TO_GLOBAL_MAP=67,
1355d643ce63SMatthew Knepley                MATOP_SET_VALUES_LOCAL=68,
1356d643ce63SMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=69,
1357d643ce63SMatthew Knepley                MATOP_GET_ROW_MAX=70,
1358d643ce63SMatthew Knepley                MATOP_CONVERT=71,
1359d643ce63SMatthew Knepley                MATOP_SET_COLORING=72,
1360d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1361d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1362d643ce63SMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1363d643ce63SMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1364ba16a8cbSKris Buschelman                MATOP_MULT_CON=77,
1365ba16a8cbSKris Buschelman                MATOP_MULT_TRANSPOSE_CON=78,
1366ba16a8cbSKris Buschelman                MATOP_ILU_FACTOR_SYMBOLIC_CON=79,
1367d643ce63SMatthew Knepley                MATOP_PERMUTE_SPARSIFY=80,
1368d643ce63SMatthew Knepley                MATOP_MULT_MULTIPLE=81,
136941acf15aSKris Buschelman                MATOP_SOLVE_MULTIPLE=82,
137041acf15aSKris Buschelman                MATOP_GET_INERTIA=83,
137141acf15aSKris Buschelman                MATOP_LOAD=84,
137241acf15aSKris Buschelman                MATOP_IS_SYMMETRIC=85,
137341acf15aSKris Buschelman                MATOP_IS_HERMITIAN=86,
137441acf15aSKris Buschelman                MATOP_IS_STRUCTURALLY_SYMMETRIC=87,
137541acf15aSKris Buschelman                MATOP_PB_RELAX=88,
137641acf15aSKris Buschelman                MATOP_GET_VECS=89,
137741acf15aSKris Buschelman                MATOP_MAT_MULT=90,
137841acf15aSKris Buschelman                MATOP_MAT_MULT_SYMBOLIC=91,
137941acf15aSKris Buschelman                MATOP_MAT_MULT_NUMERIC=92,
138041acf15aSKris Buschelman                MATOP_PTAP=93,
138141acf15aSKris Buschelman                MATOP_PTAP_SYMBOLIC=94,
1382bc011b1eSHong Zhang                MATOP_PTAP_NUMERIC=95,
1383bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE=96,
1384bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97,
13857ba1a0bfSKris Buschelman                MATOP_MAT_MULTTRANSPOSE_NUMERIC=98,
13867ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_SEQAIJ=99,
13877ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_SEQAIJ=100,
13887ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_MPIAIJ=101,
13897ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_MPIAIJ=102
1390fae171e0SBarry Smith              } MatOperation;
1391be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1392be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1393be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1394be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1395112a2221SBarry Smith 
139690ace30eSBarry Smith /*
139790ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
139890ace30eSBarry Smith  stored in a universal format. By changing the format with
1399fb9695e5SSatish Balay  PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will
140090ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
140190ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
140290ace30eSBarry Smith  read into matrices of the same time.
140390ace30eSBarry Smith */
140490ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
140590ace30eSBarry Smith 
1406be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1407be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsGetColor(Mat,ISColoring *);
1408860d1616SSatish Balay 
1409be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
14101f4e1ec7SBarry Smith 
1411d9274352SBarry Smith /*S
1412d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1413d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1414d9274352SBarry Smith 
1415f7a9e4ceSBarry Smith    Level: advanced
1416d9274352SBarry Smith 
1417d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1418d9274352SBarry Smith 
14196e1639daSBarry Smith   Users manual sections:
14206e1639daSBarry Smith .   sec_singular
14216e1639daSBarry Smith 
1422d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1423d9274352SBarry Smith S*/
142474637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1425d9274352SBarry Smith 
1426be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1427be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1428be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1429be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
1430be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat);
143174637425SBarry Smith 
1432be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1433be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1434be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
143546129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
14363f1d51d7SBarry Smith 
1437be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1438be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1439be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1440c4f061fbSSatish Balay 
1441be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1442b0a32e0cSBarry Smith 
1443be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
144404f1ad80SBarry Smith 
1445be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1446be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
14477dbadf16SMatthew Knepley 
1448e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
14492eac72dbSBarry Smith #endif
1450