xref: /petsc/include/petscmat.h (revision 81824310db631689e1e81dca1f16488c1a4fb62d)
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"
74*81824310SBarry Smith #define MATSEQCRL          "seqcrl"
75c0aa2d19SHong Zhang #define MATPLAPACK         "plapack"
76f69a0ea3SMatthew Knepley #define MatType const char*
77d91e6319SBarry Smith 
78c06d978dSMatthew Knepley /* Logging support */
79552e946dSBarry Smith #define    MAT_FILE_COOKIE 1211216    /* used to indicate matrices in binary files */
80be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE;
81be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
82be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE;
83be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE;
846849ba73SBarry Smith extern PetscEvent  MAT_Mult, MAT_MultMatrixFree, MAT_Mults, MAT_MultConstrained, MAT_MultAdd, MAT_MultTranspose;
856849ba73SBarry Smith extern PetscEvent  MAT_MultTransposeConstrained, MAT_MultTransposeAdd, MAT_Solve, MAT_Solves, MAT_SolveAdd, MAT_SolveTranspose;
866849ba73SBarry Smith extern PetscEvent  MAT_SolveTransposeAdd, MAT_Relax, MAT_ForwardSolve, MAT_BackwardSolve, MAT_LUFactor, MAT_LUFactorSymbolic;
876849ba73SBarry Smith extern PetscEvent  MAT_LUFactorNumeric, MAT_CholeskyFactor, MAT_CholeskyFactorSymbolic, MAT_CholeskyFactorNumeric, MAT_ILUFactor;
886849ba73SBarry Smith extern PetscEvent  MAT_ILUFactorSymbolic, MAT_ICCFactorSymbolic, MAT_Copy, MAT_Convert, MAT_Scale, MAT_AssemblyBegin;
896849ba73SBarry Smith extern PetscEvent  MAT_AssemblyEnd, MAT_SetValues, MAT_GetValues, MAT_GetRow, MAT_GetSubMatrices, MAT_GetColoring, MAT_GetOrdering;
906849ba73SBarry Smith extern PetscEvent  MAT_IncreaseOverlap, MAT_Partitioning, MAT_ZeroEntries, MAT_Load, MAT_View, MAT_AXPY, MAT_FDColoringCreate;
916849ba73SBarry Smith extern PetscEvent  MAT_FDColoringApply, MAT_Transpose, MAT_FDColoringFunction;
92cb00f407SKris Buschelman extern PetscEvent  MAT_MatMult, MAT_MatMultSymbolic, MAT_MatMultNumeric;
93534c1384SKris Buschelman extern PetscEvent  MAT_PtAP, MAT_PtAPSymbolic, MAT_PtAPNumeric;
94bc011b1eSHong Zhang extern PetscEvent  MAT_MatMultTranspose, MAT_MatMultTransposeSymbolic, MAT_MatMultTransposeNumeric;
95c06d978dSMatthew Knepley 
96ceb03754SKris Buschelman /*E
97ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
98ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
99ceb03754SKris Buschelman 
100ceb03754SKris Buschelman     Level: beginner
101ceb03754SKris Buschelman 
102ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
103ceb03754SKris Buschelman 
1040c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
105ceb03754SKris Buschelman E*/
106ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse;
107ceb03754SKris Buschelman 
108be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(char *);
109c06d978dSMatthew Knepley 
110f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
111a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
112a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
113f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
114f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,MatType);
115be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
116be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
117be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
118be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
11930de9b25SBarry Smith 
120f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
121f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
122f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
123f69a0ea3SMatthew Knepley 
12430de9b25SBarry Smith /*MC
12530de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
12630de9b25SBarry Smith 
12730de9b25SBarry Smith    Synopsis:
128c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
12930de9b25SBarry Smith 
13030de9b25SBarry Smith    Not Collective
13130de9b25SBarry Smith 
13230de9b25SBarry Smith    Input Parameters:
13330de9b25SBarry Smith +  name - name of a new user-defined matrix type
13430de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
13530de9b25SBarry Smith .  name_create - name of routine to create method context
13630de9b25SBarry Smith -  routine_create - routine to create method context
13730de9b25SBarry Smith 
13830de9b25SBarry Smith    Notes:
13930de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
14030de9b25SBarry Smith 
14130de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
14230de9b25SBarry Smith    is ignored.
14330de9b25SBarry Smith 
14430de9b25SBarry Smith    Sample usage:
14530de9b25SBarry Smith .vb
14630de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
14730de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
14830de9b25SBarry Smith .ve
14930de9b25SBarry Smith 
15030de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
15130de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
15230de9b25SBarry Smith    or at runtime via the option
15330de9b25SBarry Smith $     -mat_type my_mat
15430de9b25SBarry Smith 
15530de9b25SBarry Smith    Level: advanced
15630de9b25SBarry Smith 
157ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
15830de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
15930de9b25SBarry Smith 
16030de9b25SBarry Smith .keywords: Mat, register
16130de9b25SBarry Smith 
16230de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
16330de9b25SBarry Smith 
16430de9b25SBarry Smith M*/
165273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
166273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
167273d9f13SBarry Smith #else
168273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
16930de9b25SBarry Smith #endif
17030de9b25SBarry Smith 
171273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
172b0a32e0cSBarry Smith extern PetscFList MatList;
17328988994SBarry Smith 
174be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
175be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
176be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
177ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
178ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
179ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
180ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
181ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
182ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
183ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
184be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
185ba966639SSatish 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)
186ba966639SSatish 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)
187ba966639SSatish 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)
188ba966639SSatish 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)
189ba966639SSatish 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))
190ba966639SSatish 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))
191ba966639SSatish 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))
192ba966639SSatish 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)
193ba966639SSatish 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)
194ba966639SSatish 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)
195ba966639SSatish 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)
196ba966639SSatish 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))
197ba966639SSatish 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))
198ba966639SSatish 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))
199be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
200be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*);
201be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*);
202be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
203ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,nz,nnz,&A),Mat,A)
204ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
205ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
206ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
207ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
208ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
209ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
210be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
211ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(comm,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
212ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
213ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(comm,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
214ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
215ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
216ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(comm,bs,m,n,M,N,0,nnz,0,onz,A))
217ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
218ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
219ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
220ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
221ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
222ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
223ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,A))
224ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
225be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
226be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
227ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,nz,nnz,&A),Mat,A)
228ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
229ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
230ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
231ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
232ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
233ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
234be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
235ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(comm,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
236ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
237ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(comm,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
238ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
239ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
240ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(comm,bs,m,n,M,N,0,nnz,0,onz,A))
241ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
242ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
243ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
244ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
245ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
246ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
247ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,A))
248ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
249be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
250ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(comm,m,n,M,N,ctx,&A),Mat,A)
251ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
252be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
253be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
254ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
255ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*);
256284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
257f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
258be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
25921c89e3eSBarry Smith 
260be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPrintHelp(Mat);
261be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetPetscMaps(Mat,PetscMap*,PetscMap*);
262ec0117caSBarry Smith 
2630c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
26499cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat);
26599cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat);
26699cafbc1SBarry Smith 
2678ed539a5SBarry Smith /* ------------------------------------------------------------*/
268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
269be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
27087d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
271f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
27284cb2905SBarry Smith 
2732ef4de8bSBarry Smith /*S
2742ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
2752ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
2762ef4de8bSBarry Smith 
2772ef4de8bSBarry Smith    Level: beginner
2782ef4de8bSBarry Smith 
2792ef4de8bSBarry Smith   Concepts: matrix; linear operator
2802ef4de8bSBarry Smith 
281d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
2822ef4de8bSBarry Smith S*/
283435da068SBarry Smith typedef struct {
284c1ac3661SBarry Smith   PetscInt k,j,i,c;
285435da068SBarry Smith } MatStencil;
2862ef4de8bSBarry Smith 
287be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
288be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
289be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
290435da068SBarry Smith 
291be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
292be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
2943a7fca6bSBarry Smith 
295d91e6319SBarry Smith /*E
296d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
297d91e6319SBarry Smith      to continue to add values to it
298d91e6319SBarry Smith 
299d91e6319SBarry Smith     Level: beginner
300d91e6319SBarry Smith 
301d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
302d91e6319SBarry Smith E*/
3036d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
304be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
305be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
306be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
3074f9c727eSBarry Smith 
308be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt    MatSetValue_Row;
309be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt    MatSetValue_Column;
310be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscScalar MatSetValue_Value;
3111d73ed98SBarry Smith 
31230de9b25SBarry Smith /*MC
31330de9b25SBarry Smith    MatSetValue - Set a single entry into a matrix.
31430de9b25SBarry Smith 
31530de9b25SBarry Smith    Synopsis:
316c1ac3661SBarry Smith    PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode);
31730de9b25SBarry Smith 
31830de9b25SBarry Smith    Not collective
31930de9b25SBarry Smith 
32030de9b25SBarry Smith    Input Parameters:
32130de9b25SBarry Smith +  m - the matrix
32230de9b25SBarry Smith .  row - the row location of the entry
32330de9b25SBarry Smith .  col - the column location of the entry
32430de9b25SBarry Smith .  value - the value to insert
32530de9b25SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
32630de9b25SBarry Smith 
32730de9b25SBarry Smith    Notes:
32830de9b25SBarry Smith    For efficiency one should use MatSetValues() and set several or many
32930de9b25SBarry Smith    values simultaneously if possible.
33030de9b25SBarry Smith 
33130de9b25SBarry Smith    Level: beginner
33230de9b25SBarry Smith 
33330de9b25SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
33430de9b25SBarry Smith M*/
335b951964fSBarry Smith #define MatSetValue(v,i,j,va,mode) \
336f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \
3371d73ed98SBarry Smith    MatSetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode))
33830de9b25SBarry Smith 
339ea06a074SBarry Smith #define MatGetValue(v,i,j,va) \
340f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,0) || \
3411d73ed98SBarry Smith    MatGetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&va))
34230de9b25SBarry Smith 
343d91e6319SBarry Smith #define MatSetValueLocal(v,i,j,va,mode) \
344f1144a30SSatish Balay   ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \
3451d73ed98SBarry Smith    MatSetValuesLocal(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode))
34630de9b25SBarry Smith 
347d91e6319SBarry Smith /*E
348d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
349d91e6319SBarry Smith 
350d91e6319SBarry Smith     Level: beginner
351d91e6319SBarry Smith 
3520a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
353d91e6319SBarry Smith 
354d91e6319SBarry Smith .seealso: MatSetOption()
355d91e6319SBarry Smith E*/
3566d4a8577SBarry Smith typedef enum {MAT_ROW_ORIENTED=1,MAT_COLUMN_ORIENTED=2,MAT_ROWS_SORTED=4,
3576d4a8577SBarry Smith               MAT_COLUMNS_SORTED=8,MAT_NO_NEW_NONZERO_LOCATIONS=16,
3586d4a8577SBarry Smith               MAT_YES_NEW_NONZERO_LOCATIONS=32,MAT_SYMMETRIC=64,
3596ca9ecd3SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC=65,MAT_NO_NEW_DIAGONALS=66,
3606ca9ecd3SBarry Smith               MAT_YES_NEW_DIAGONALS=67,MAT_INODE_LIMIT_1=68,MAT_INODE_LIMIT_2=69,
3616ca9ecd3SBarry Smith               MAT_INODE_LIMIT_3=70,MAT_INODE_LIMIT_4=71,MAT_INODE_LIMIT_5=72,
3626ca9ecd3SBarry Smith               MAT_IGNORE_OFF_PROC_ENTRIES=73,MAT_ROWS_UNSORTED=74,
3634787f768SSatish Balay               MAT_COLUMNS_UNSORTED=75,MAT_NEW_NONZERO_LOCATION_ERR=76,
3647c922b88SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR=77,MAT_USE_HASH_TABLE=78,
3652bad1931SBarry Smith               MAT_KEEP_ZEROED_ROWS=79,MAT_IGNORE_ZERO_ENTRIES=80,MAT_USE_INODES=81,
366efcf0fc3SBarry Smith               MAT_DO_NOT_USE_INODES=82,MAT_NOT_SYMMETRIC=83,MAT_HERMITIAN=84,
3679a4540c5SBarry Smith               MAT_NOT_STRUCTURALLY_SYMMETRIC=85,MAT_NOT_HERMITIAN=86,
368d487561eSHong Zhang               MAT_SYMMETRY_ETERNAL=87,MAT_NOT_SYMMETRY_ETERNAL=88,
369941593c8SHong Zhang               MAT_USE_COMPRESSEDROW=89,MAT_DO_NOT_USE_COMPRESSEDROW=90,
370f5edf698SHong Zhang               MAT_IGNORE_LOWER_TRIANGULAR=91,MAT_ERROR_LOWER_TRIANGULAR=92,MAT_GETROW_UPPERTRIANGULAR=93} MatOption;
371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption);
372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,MatType*);
373ba966639SSatish Balay PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),MatType,t)
37484cb2905SBarry Smith 
375be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
376be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
378f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat);
379f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat);
380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
382be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
384ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
385be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
386be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
387ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
388be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
3897b80b807SBarry Smith 
390be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
391be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
392ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
393be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
394be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
395ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
396ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0.0,&t),PetscTruth,t)
397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
398ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
400be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
4012eac72dbSBarry Smith 
402f5edf698SHong Zhang 
403f5edf698SHong Zhang 
404d91e6319SBarry Smith /*E
405d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
406d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
407d91e6319SBarry Smith 
408d91e6319SBarry Smith     Level: beginner
409d91e6319SBarry Smith 
410d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
411d91e6319SBarry Smith 
412d91e6319SBarry Smith .seealso: MatDuplicate()
413d91e6319SBarry Smith E*/
4142e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
4152e8a6d31SBarry Smith 
416be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvertRegister(const char[],const char[],const char[],PetscErrorCode (*)(Mat,MatType,MatReuse,Mat*));
417273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
418273d9f13SBarry Smith #define MatConvertRegisterDynamic(a,b,c,d) MatConvertRegister(a,b,c,0)
419273d9f13SBarry Smith #else
420273d9f13SBarry Smith #define MatConvertRegisterDynamic(a,b,c,d) MatConvertRegister(a,b,c,d)
421273d9f13SBarry Smith #endif
422be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvertRegisterAll(const char[]);
423be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvertRegisterDestroy(void);
424273d9f13SBarry Smith extern PetscTruth MatConvertRegisterAllCalled;
425b0a32e0cSBarry Smith extern PetscFList MatConvertList;
426f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,MatType,MatReuse,Mat*);
427ba966639SSatish Balay PetscPolymorphicFunction(MatConvert,(Mat A,MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
428be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
429ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
430ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
43194a9d846SBarry Smith 
432d91e6319SBarry Smith /*E
433d91e6319SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
434d91e6319SBarry Smith 
435d91e6319SBarry Smith     Level: beginner
436d91e6319SBarry Smith 
437d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
438d91e6319SBarry Smith 
43994b7f48cSBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
440d91e6319SBarry Smith E*/
441c537a176SHong Zhang typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
442cb5b572fSBarry Smith 
443be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
444be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
445be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
446ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
447ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0.0,&t),PetscTruth,t)
448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
449ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
450be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscTruth*);
451ba966639SSatish Balay PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,&t),PetscTruth,t)
452be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
453be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
454f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,MatType,Mat*);
455ba966639SSatish Balay PetscPolymorphicFunction(MatLoad,(PetscViewer v,MatType t),(v,t,&a),Mat,a)
4567b80b807SBarry Smith 
457be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
458be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
459be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
460be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
461d4fbbf0eSBarry Smith 
462d91e6319SBarry Smith /*S
463d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
464d91e6319SBarry Smith 
465d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
466d91e6319SBarry Smith 
467d91e6319SBarry Smith    Level: intermediate
468d91e6319SBarry Smith 
469d91e6319SBarry Smith   Concepts: matrix^nonzero information
470d91e6319SBarry Smith 
471d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
472d91e6319SBarry Smith S*/
4734e220ebcSLois Curfman McInnes typedef struct {
474b0a32e0cSBarry Smith   PetscLogDouble rows_global,columns_global;         /* number of global rows and columns */
475b0a32e0cSBarry Smith   PetscLogDouble rows_local,columns_local;           /* number of local rows and columns */
476b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
477b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
478b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
479b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
480b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
481b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
482b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4834e220ebcSLois Curfman McInnes } MatInfo;
4844e220ebcSLois Curfman McInnes 
485d9274352SBarry Smith /*E
486d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
487d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
488d9274352SBarry Smith 
489d9274352SBarry Smith     Level: beginner
490d9274352SBarry Smith 
491d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
492d9274352SBarry Smith 
493d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
494d9274352SBarry Smith E*/
4957b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
496be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
497be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
499be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec);
500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,Mat*);
501ba966639SSatish Balay PetscPolymorphicFunction(MatTranspose,(Mat A),(A,&t),Mat,t)
502be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
503ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
504be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
505be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
506be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
508ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
509be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
5137b80b807SBarry Smith 
514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
515ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
517f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
518f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
519f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
520f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5217b80b807SBarry Smith 
522be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
523be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
524be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5255ef9f2a5SBarry Smith 
526be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5297b80b807SBarry Smith 
530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *);
533abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *);
534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
535be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
536be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
538be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
542be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscScalar**,Mat*);
54343eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
544cd116e26SSatish Balay #include "petscctable.h"
54543eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
54643eb5e2fSMatthew Knepley #else
54743eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
54843eb5e2fSMatthew Knepley #endif
5498efafbd8SBarry Smith 
550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5517b80b807SBarry Smith 
552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
553be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
55522440eb1SKris Buschelman 
556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
557be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
558be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
55922440eb1SKris Buschelman 
560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
561be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
562be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
563bc011b1eSHong Zhang 
564f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
565f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat);
566be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat);
5671c741599SHong Zhang 
568f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
569f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5707b80b807SBarry Smith 
571be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
572be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
573f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
574f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
575be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
576be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
577052efed2SBarry Smith 
578be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
579be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
58090f02eecSBarry Smith 
581be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
5820c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
583ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
584be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
585be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
586bd481603SBarry Smith 
587bd481603SBarry Smith /*MC
588bd481603SBarry Smith    MatPreallocInitialize - Begins the block of code that will count the number of nonzeros per
589bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
590bd481603SBarry Smith 
591bd481603SBarry Smith    Synopsis:
592c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
593bd481603SBarry Smith 
594bd481603SBarry Smith    Collective on MPI_Comm
595bd481603SBarry Smith 
596bd481603SBarry Smith    Input Parameters:
597bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
598859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
599859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
600bd481603SBarry Smith 
601bd481603SBarry Smith    Output Parameters:
602bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
603bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
604bd481603SBarry Smith 
605bd481603SBarry Smith 
606bd481603SBarry Smith    Level: intermediate
607bd481603SBarry Smith 
608bd481603SBarry Smith    Notes:
609bd481603SBarry Smith    See the chapter in the users manual on performance for more details
610bd481603SBarry Smith 
6111d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
612bd481603SBarry Smith 
613bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
614bd481603SBarry Smith 
615bd481603SBarry Smith   Concepts: preallocation^Matrix
616bd481603SBarry Smith 
617bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
618bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
619bd481603SBarry Smith M*/
620c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6217c922b88SBarry Smith { \
622c1ac3661SBarry Smith   PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
623c1ac3661SBarry Smith   _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\
624c1ac3661SBarry Smith   _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\
625c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
626c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp;
6277c922b88SBarry Smith 
628bd481603SBarry Smith /*MC
629bd481603SBarry Smith    MatPreallocSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
630bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
631bd481603SBarry Smith 
632bd481603SBarry Smith    Synopsis:
633c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
634bd481603SBarry Smith 
635bd481603SBarry Smith    Collective on MPI_Comm
636bd481603SBarry Smith 
637bd481603SBarry Smith    Input Parameters:
638bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
639859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
640859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
641bd481603SBarry Smith 
642bd481603SBarry Smith    Output Parameters:
643bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
644bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
645bd481603SBarry Smith 
646bd481603SBarry Smith 
647bd481603SBarry Smith    Level: intermediate
648bd481603SBarry Smith 
649bd481603SBarry Smith    Notes:
650bd481603SBarry Smith    See the chapter in the users manual on performance for more details
651bd481603SBarry Smith 
6521d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
653bd481603SBarry Smith 
654bd481603SBarry Smith   Concepts: preallocation^Matrix
655bd481603SBarry Smith 
656bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
657bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
658bd481603SBarry Smith M*/
659222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
660222b16d4SBarry Smith { \
661c1ac3661SBarry Smith   PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__end; \
662c1ac3661SBarry Smith   _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\
663c1ac3661SBarry Smith   _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\
664c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
665c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp;
666222b16d4SBarry Smith 
667bd481603SBarry Smith /*MC
668bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
669bd481603SBarry Smith        inserted using a local number of the rows and columns
670bd481603SBarry Smith 
671bd481603SBarry Smith    Synopsis:
672c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
673bd481603SBarry Smith 
674bd481603SBarry Smith    Not Collective
675bd481603SBarry Smith 
676bd481603SBarry Smith    Input Parameters:
677bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
678bd481603SBarry Smith .  nrows - the number of rows indicated
6791d73ed98SBarry Smith .  rows - the indices of the rows
680bd481603SBarry Smith .  ncols - the number of columns in the matrix
681bd481603SBarry Smith .  cols - the columns indicated
682bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
683bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
684bd481603SBarry Smith 
685bd481603SBarry Smith 
686bd481603SBarry Smith    Level: intermediate
687bd481603SBarry Smith 
688bd481603SBarry Smith    Notes:
689bd481603SBarry Smith    See the chapter in the users manual on performance for more details
690bd481603SBarry Smith 
6911d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
692bd481603SBarry Smith 
693bd481603SBarry Smith   Concepts: preallocation^Matrix
694bd481603SBarry Smith 
695bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
696bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
697bd481603SBarry Smith M*/
698c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
699c4f061fbSSatish Balay {\
700c1ac3661SBarry Smith   PetscInt __l;\
701ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
702ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
703c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
704ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
705c4f061fbSSatish Balay   }\
706c4f061fbSSatish Balay }
707c4f061fbSSatish Balay 
708bd481603SBarry Smith /*MC
709bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
710bd481603SBarry Smith        inserted using a local number of the rows and columns
711bd481603SBarry Smith 
712bd481603SBarry Smith    Synopsis:
713c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
714bd481603SBarry Smith 
715bd481603SBarry Smith    Not Collective
716bd481603SBarry Smith 
717bd481603SBarry Smith    Input Parameters:
718bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
719bd481603SBarry Smith .  nrows - the number of rows indicated
7201d73ed98SBarry Smith .  rows - the indices of the rows
721bd481603SBarry Smith .  ncols - the number of columns in the matrix
722bd481603SBarry Smith .  cols - the columns indicated
723bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
724bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
725bd481603SBarry Smith 
726bd481603SBarry Smith 
727bd481603SBarry Smith    Level: intermediate
728bd481603SBarry Smith 
729bd481603SBarry Smith    Notes:
730bd481603SBarry Smith    See the chapter in the users manual on performance for more details
731bd481603SBarry Smith 
732bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
733bd481603SBarry Smith 
734bd481603SBarry Smith   Concepts: preallocation^Matrix
735bd481603SBarry Smith 
736bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
737bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
738bd481603SBarry Smith M*/
739d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
740d3d32019SBarry Smith {\
741c1ac3661SBarry Smith   PetscInt __l;\
742d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
743d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
744d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
745d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
746d3d32019SBarry Smith   }\
747d3d32019SBarry Smith }
748d3d32019SBarry Smith 
749bd481603SBarry Smith /*MC
750bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
751bd481603SBarry Smith        inserted using a local number of the rows and columns
752bd481603SBarry Smith 
753bd481603SBarry Smith    Synopsis:
754c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
755bd481603SBarry Smith 
756bd481603SBarry Smith    Not Collective
757bd481603SBarry Smith 
758bd481603SBarry Smith    Input Parameters:
759bd481603SBarry Smith +  nrows - the number of rows indicated
7601d73ed98SBarry Smith .  rows - the indices of the rows
761bd481603SBarry Smith .  ncols - the number of columns in the matrix
762a50f8bf6SHong Zhang -  cols - the columns indicated
763a50f8bf6SHong Zhang 
764a50f8bf6SHong Zhang    Output Parameters:
765a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
766bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
767bd481603SBarry Smith 
768bd481603SBarry Smith 
769bd481603SBarry Smith    Level: intermediate
770bd481603SBarry Smith 
771bd481603SBarry Smith    Notes:
772bd481603SBarry Smith    See the chapter in the users manual on performance for more details
773bd481603SBarry Smith 
774bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
775bd481603SBarry Smith 
776bd481603SBarry Smith   Concepts: preallocation^Matrix
777bd481603SBarry Smith 
778bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
779bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
780bd481603SBarry Smith M*/
781c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
782c1ac3661SBarry Smith { PetscInt __i; \
7837c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
7847c922b88SBarry Smith     if (cols[__i] < __start || cols[__i] >= __end) onz[row - __rstart]++; \
7857c922b88SBarry Smith   }\
7867c922b88SBarry Smith   dnz[row - __rstart] = nc - onz[row - __rstart];\
7877c922b88SBarry Smith }
7887c922b88SBarry Smith 
789bd481603SBarry Smith /*MC
790bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
791bd481603SBarry Smith        inserted using a local number of the rows and columns
792bd481603SBarry Smith 
793bd481603SBarry Smith    Synopsis:
794c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
795bd481603SBarry Smith 
796bd481603SBarry Smith    Not Collective
797bd481603SBarry Smith 
798bd481603SBarry Smith    Input Parameters:
799bd481603SBarry Smith +  nrows - the number of rows indicated
8001d73ed98SBarry Smith .  rows - the indices of the rows
801bd481603SBarry Smith .  ncols - the number of columns in the matrix
802bd481603SBarry Smith .  cols - the columns indicated
803bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
804bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
805bd481603SBarry Smith 
806bd481603SBarry Smith 
807bd481603SBarry Smith    Level: intermediate
808bd481603SBarry Smith 
809bd481603SBarry Smith    Notes:
810bd481603SBarry Smith    See the chapter in the users manual on performance for more details
811bd481603SBarry Smith 
812bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
813bd481603SBarry Smith 
814bd481603SBarry Smith   Concepts: preallocation^Matrix
815bd481603SBarry Smith 
816bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
817bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
818bd481603SBarry Smith M*/
819d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
820c1ac3661SBarry Smith { PetscInt __i; \
821d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
822d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
823d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
824d3d32019SBarry Smith   }\
825d3d32019SBarry Smith }
826d3d32019SBarry Smith 
827bd481603SBarry Smith /*MC
828bd481603SBarry Smith    MatPreallocFinalize - Ends the block of code that will count the number of nonzeros per
829bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
830bd481603SBarry Smith 
831bd481603SBarry Smith    Synopsis:
832c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
833bd481603SBarry Smith 
834bd481603SBarry Smith    Collective on MPI_Comm
835bd481603SBarry Smith 
836bd481603SBarry Smith    Input Parameters:
837bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
838bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
839bd481603SBarry Smith 
840bd481603SBarry Smith 
841bd481603SBarry Smith    Level: intermediate
842bd481603SBarry Smith 
843bd481603SBarry Smith    Notes:
844bd481603SBarry Smith    See the chapter in the users manual on performance for more details
845bd481603SBarry Smith 
846bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
847bd481603SBarry Smith 
848bd481603SBarry Smith   Concepts: preallocation^Matrix
849bd481603SBarry Smith 
850bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
851bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
852bd481603SBarry Smith M*/
853ef66eb69SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree(dnz);CHKERRQ(_4_ierr);}
8547c922b88SBarry Smith 
855bd481603SBarry Smith 
856bd481603SBarry Smith 
8577b80b807SBarry Smith /* Routines unique to particular data structures */
858be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
859ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
860be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBDiagGetData(Mat,PetscInt*,PetscInt*,PetscInt*[],PetscInt*[],PetscScalar***);
861698d4c6aSKris Buschelman 
862be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
863be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
864698d4c6aSKris Buschelman 
865be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
866be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
867be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
8687b80b807SBarry Smith 
869a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
870a96a251dSBarry Smith 
871be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
872ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
873be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
874ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
875be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
876ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
877be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDensePreallocation(Mat,PetscScalar[]);
878be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]);
879273d9f13SBarry Smith 
880be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
881ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
882be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
883be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
884be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
885be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
886be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDensePreallocation(Mat,PetscScalar[]);
887be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]);
888be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
889be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
890284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
891be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
892be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
893be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
894be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
895273d9f13SBarry Smith 
896be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
897ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*);
8981b807ce4Svictorle 
899be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
900be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
9012e8a6d31SBarry Smith 
902be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
9033a7fca6bSBarry Smith 
9047b80b807SBarry Smith /*
9057b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
90694b7f48cSBarry Smith   done through the KSP and PC interfaces.
9077b80b807SBarry Smith */
9087b80b807SBarry Smith 
909d9274352SBarry Smith /*E
910d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
911d9274352SBarry Smith        with an optional dynamic library name, for example
912d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
913d9274352SBarry Smith 
914d9274352SBarry Smith    Level: beginner
915d9274352SBarry Smith 
916d9274352SBarry Smith .seealso: MatGetOrdering()
917d9274352SBarry Smith E*/
91849773a63SBarry Smith #define MatOrderingType char*
919b12f92e5SBarry Smith #define MATORDERING_NATURAL   "natural"
920b12f92e5SBarry Smith #define MATORDERING_ND        "nd"
921b12f92e5SBarry Smith #define MATORDERING_1WD       "1wd"
922b12f92e5SBarry Smith #define MATORDERING_RCM       "rcm"
923b12f92e5SBarry Smith #define MATORDERING_QMD       "qmd"
924b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength"
92562152c8bSBarry Smith #define MATORDERING_DSC_ND    "dsc_nd"
92662152c8bSBarry Smith #define MATORDERING_DSC_MMD   "dsc_mmd"
92762152c8bSBarry Smith #define MATORDERING_DSC_MDF   "dsc_mdf"
928c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
929c06d978dSMatthew Knepley #define MATORDERING_IDENTITY  "identity"
930c06d978dSMatthew Knepley #define MATORDERING_REVERSE   "reverse"
931b12f92e5SBarry Smith 
932be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
933be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
934be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
93530de9b25SBarry Smith 
93630de9b25SBarry Smith /*MC
93730de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
93830de9b25SBarry Smith                                matrix package.
93930de9b25SBarry Smith 
94030de9b25SBarry Smith    Synopsis:
941c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
94230de9b25SBarry Smith 
94330de9b25SBarry Smith    Not Collective
94430de9b25SBarry Smith 
94530de9b25SBarry Smith    Input Parameters:
94630de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
94730de9b25SBarry Smith .  path - location of library where creation routine is
94830de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
94930de9b25SBarry Smith -  function - function pointer that creates the ordering
95030de9b25SBarry Smith 
95130de9b25SBarry Smith    Level: developer
95230de9b25SBarry Smith 
95330de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
95430de9b25SBarry Smith    is ignored.
95530de9b25SBarry Smith 
95630de9b25SBarry Smith    Sample usage:
95730de9b25SBarry Smith .vb
95830de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
95930de9b25SBarry Smith                "MyOrder",MyOrder);
96030de9b25SBarry Smith .ve
96130de9b25SBarry Smith 
96230de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
96330de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
96430de9b25SBarry Smith    or at runtime via the option
96530de9b25SBarry Smith $     -pc_ilu_mat_ordering_type my_order
96630de9b25SBarry Smith $     -pc_lu_mat_ordering_type my_order
96730de9b25SBarry Smith 
968ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
96930de9b25SBarry Smith 
97030de9b25SBarry Smith .keywords: matrix, ordering, register
97130de9b25SBarry Smith 
97230de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
97330de9b25SBarry Smith M*/
974aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
975f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
976b12f92e5SBarry Smith #else
977f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
978b12f92e5SBarry Smith #endif
97930de9b25SBarry Smith 
980be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
981be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
9822bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
983b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
984d4fbbf0eSBarry Smith 
985be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
986a2ce50c7SBarry Smith 
987d91e6319SBarry Smith /*S
98815e8a5b3SHong Zhang    MatFactorInfo - Data based into the matrix factorization routines
989b00f7748SHong Zhang 
99015e8a5b3SHong Zhang    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE
991b00f7748SHong Zhang 
99215e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
993b00f7748SHong Zhang 
994b00f7748SHong Zhang    Level: developer
995b00f7748SHong Zhang 
996d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
997d7d82daaSBarry Smith           MatFactorInfoInitialize()
998b00f7748SHong Zhang 
999b00f7748SHong Zhang S*/
1000b00f7748SHong Zhang typedef struct {
10010a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
1002fbf22428SSatish Balay   PetscReal     shiftpd;        /* if true, shift until positive pivots */
10032cea7109SBarry Smith   PetscReal     shift_fraction; /* record shift fraction taken */
100415e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
100515e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1006b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
100715e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
1008f6275e2eSBarry Smith   PetscReal     fill;           /* expected fill; nonzeros in factored matrix/nonzeros in original matrix*/
1009348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1010bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1011bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
101215e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
101315e8a5b3SHong Zhang } MatFactorInfo;
1014ffa6d0a5SLois Curfman McInnes 
1015be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
1016be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,MatFactorInfo*);
1017be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1018be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,MatFactorInfo*,Mat*);
1019be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,MatFactorInfo*);
1020be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,MatFactorInfo*);
1021be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1022be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1023be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1024be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,MatFactorInfo*);
1025be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,MatFactorInfo*,Mat*);
1026be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,MatFactorInfo*,Mat *);
1027be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1028be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1029be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1030be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1031be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1032be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1033be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1034be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
10358ed539a5SBarry Smith 
1036be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1037bb5a7306SBarry Smith 
1038d91e6319SBarry Smith /*E
1039d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1040bb1eb677SSatish Balay 
1041d91e6319SBarry Smith     Level: beginner
1042d91e6319SBarry Smith 
1043d9274352SBarry Smith    May be bitwise ORd together
1044d9274352SBarry Smith 
1045d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1046d91e6319SBarry Smith 
10474e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10484e7234bfSBarry Smith 
1049a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax()
1050d91e6319SBarry Smith E*/
1051ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1052ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1053ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
105484cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1055be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1056be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10578ed539a5SBarry Smith 
1058d4fbbf0eSBarry Smith /*
1059639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1060639f9d9dSBarry Smith */
1061b12f92e5SBarry Smith 
1062d9274352SBarry Smith /*E
1063d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1064d9274352SBarry Smith        with an optional dynamic library name, for example
1065d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1066d9274352SBarry Smith 
1067d9274352SBarry Smith    Level: beginner
1068d9274352SBarry Smith 
1069d9274352SBarry Smith .seealso: MatGetColoring()
1070d9274352SBarry Smith E*/
107149773a63SBarry Smith #define MatColoringType char*
1072b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1073b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1074b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1075b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1076b12f92e5SBarry Smith 
1077be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*);
1078be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatColoringType,ISColoring *));
107930de9b25SBarry Smith 
108030de9b25SBarry Smith /*MC
108130de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
108230de9b25SBarry Smith                                matrix package.
108330de9b25SBarry Smith 
108430de9b25SBarry Smith    Synopsis:
1085c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
108630de9b25SBarry Smith 
108730de9b25SBarry Smith    Not Collective
108830de9b25SBarry Smith 
108930de9b25SBarry Smith    Input Parameters:
109030de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
109130de9b25SBarry Smith .  path - location of library where creation routine is
109230de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
109330de9b25SBarry Smith -  function - function pointer that creates the coloring
109430de9b25SBarry Smith 
109530de9b25SBarry Smith    Level: developer
109630de9b25SBarry Smith 
109730de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
109830de9b25SBarry Smith    is ignored.
109930de9b25SBarry Smith 
110030de9b25SBarry Smith    Sample usage:
110130de9b25SBarry Smith .vb
110230de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
110330de9b25SBarry Smith                "MyColor",MyColor);
110430de9b25SBarry Smith .ve
110530de9b25SBarry Smith 
110630de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
110730de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
110830de9b25SBarry Smith    or at runtime via the option
110930de9b25SBarry Smith $     -mat_coloring_type my_color
111030de9b25SBarry Smith 
1111ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
111230de9b25SBarry Smith 
111330de9b25SBarry Smith .keywords: matrix, Coloring, register
111430de9b25SBarry Smith 
111530de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
111630de9b25SBarry Smith M*/
1117aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1118f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1119b12f92e5SBarry Smith #else
1120f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1121b12f92e5SBarry Smith #endif
112230de9b25SBarry Smith 
11232bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1124f1144a30SSatish Balay 
1125be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1126be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1127be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1128639f9d9dSBarry Smith 
1129d9274352SBarry Smith /*S
1130d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1131d9274352SBarry Smith         and coloring
1132639f9d9dSBarry Smith 
1133d9274352SBarry Smith    Level: beginner
1134d9274352SBarry Smith 
1135d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1136d9274352SBarry Smith 
1137d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1138d9274352SBarry Smith S*/
1139e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1140639f9d9dSBarry Smith 
1141be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1142be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1143be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1144be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
11454a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1146be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1147be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring,PetscInt);
1148be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring,PetscInt*);
1149be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1150be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1151be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1152be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring);
1153be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1154be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1155639f9d9dSBarry Smith /*
11560752156aSBarry Smith     These routines are for partitioning matrices: currently used only
11573eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
11580752156aSBarry Smith */
1159ca161407SBarry Smith 
1160d9274352SBarry Smith /*S
1161d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1162d9274352SBarry Smith 
1163d9274352SBarry Smith    Level: beginner
1164d9274352SBarry Smith 
1165d9274352SBarry Smith   Concepts: partitioning
1166d9274352SBarry Smith 
1167743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1168d9274352SBarry Smith S*/
116991e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1170d9274352SBarry Smith 
1171d9274352SBarry Smith /*E
11725bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1173d9274352SBarry Smith        with an optional dynamic library name, for example
1174d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1175d9274352SBarry Smith 
1176d9274352SBarry Smith    Level: beginner
1177d9274352SBarry Smith 
1178b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1179d9274352SBarry Smith E*/
118049773a63SBarry Smith #define MatPartitioningType char*
11818ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1182c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
11838ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
118471306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
118571306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
118671306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
118771306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
118871306c60Sjroman 
1189ca161407SBarry Smith 
1190be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1191be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1192be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1193be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1194be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1195be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1196be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1197be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
11982aabb6bbSBarry Smith 
1199be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
120030de9b25SBarry Smith 
120130de9b25SBarry Smith /*MC
120230de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
120330de9b25SBarry Smith    matrix package.
120430de9b25SBarry Smith 
120530de9b25SBarry Smith    Synopsis:
1206c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
120730de9b25SBarry Smith 
120830de9b25SBarry Smith    Not Collective
120930de9b25SBarry Smith 
121030de9b25SBarry Smith    Input Parameters:
121130de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
121230de9b25SBarry Smith .  path - location of library where creation routine is
121330de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
121430de9b25SBarry Smith -  function - function pointer that creates the partitioning type
121530de9b25SBarry Smith 
121630de9b25SBarry Smith    Level: developer
121730de9b25SBarry Smith 
121830de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
121930de9b25SBarry Smith    is ignored.
122030de9b25SBarry Smith 
122130de9b25SBarry Smith    Sample usage:
122230de9b25SBarry Smith .vb
122330de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
122430de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
122530de9b25SBarry Smith .ve
122630de9b25SBarry Smith 
122730de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
122830de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
122930de9b25SBarry Smith    or at runtime via the option
123030de9b25SBarry Smith $     -mat_partitioning_type my_part
123130de9b25SBarry Smith 
1232ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
123330de9b25SBarry Smith 
123430de9b25SBarry Smith .keywords: matrix, partitioning, register
123530de9b25SBarry Smith 
123630de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
123730de9b25SBarry Smith M*/
1238aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1239f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
12402aabb6bbSBarry Smith #else
1241f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
12422aabb6bbSBarry Smith #endif
12432aabb6bbSBarry Smith 
12442bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1245f1144a30SSatish Balay 
1246be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1247be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
12482bad1931SBarry Smith 
1249be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1250be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1251be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1252ca161407SBarry Smith 
1253be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
12540752156aSBarry Smith 
1255be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1256be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
125771306c60Sjroman 
1258954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
1259be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
126071306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1261be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1262be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
126371306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1264be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1265be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1266be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
126771306c60Sjroman 
126871306c60Sjroman #define MP_PARTY_OPT "opt"
126971306c60Sjroman #define MP_PARTY_LIN "lin"
127071306c60Sjroman #define MP_PARTY_SCA "sca"
127171306c60Sjroman #define MP_PARTY_RAN "ran"
127271306c60Sjroman #define MP_PARTY_GBF "gbf"
127371306c60Sjroman #define MP_PARTY_GCF "gcf"
127471306c60Sjroman #define MP_PARTY_BUB "bub"
127571306c60Sjroman #define MP_PARTY_DEF "def"
1276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
127771306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
127871306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
127971306c60Sjroman #define MP_PARTY_NONE "no"
1280be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1281be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1282be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1283be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
128471306c60Sjroman 
128571306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1286be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1287be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1288be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1289be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1290be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
129171306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1292be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1294be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
129571306c60Sjroman 
12960752156aSBarry Smith /*
12970a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1298d4fbbf0eSBarry Smith */
12991c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13001c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13011c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13021c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13031c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13047c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13057c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13061c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13071c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13087c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13097c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13101c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13111c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
13121c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
13131c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13141c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13151c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13161c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13171c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13181c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13191c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13201c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
13211c1c02c0SLois Curfman McInnes                MATOP_COMPRESS=22,
13221c1c02c0SLois Curfman McInnes                MATOP_SET_OPTION=23,
13231c1c02c0SLois Curfman McInnes                MATOP_ZERO_ENTRIES=24,
13241c1c02c0SLois Curfman McInnes                MATOP_ZERO_ROWS=25,
13251c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_SYMBOLIC=26,
13261c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_NUMERIC=27,
13271c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_SYMBOLIC=28,
13281c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_NUMERIC=29,
1329d643ce63SMatthew Knepley                MATOP_SETUP_PREALLOCATION=30,
1330d643ce63SMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=31,
1331d643ce63SMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=32,
1332d643ce63SMatthew Knepley                MATOP_GET_ARRAY=33,
1333d643ce63SMatthew Knepley                MATOP_RESTORE_ARRAY=34,
1334e26f1c3cSBarry Smith                MATOP_DUPLICATE=35,
1335d643ce63SMatthew Knepley                MATOP_FORWARD_SOLVE=36,
1336d643ce63SMatthew Knepley                MATOP_BACKWARD_SOLVE=37,
1337d643ce63SMatthew Knepley                MATOP_ILUFACTOR=38,
1338d643ce63SMatthew Knepley                MATOP_ICCFACTOR=39,
1339d643ce63SMatthew Knepley                MATOP_AXPY=40,
1340d643ce63SMatthew Knepley                MATOP_GET_SUBMATRICES=41,
1341d643ce63SMatthew Knepley                MATOP_INCREASE_OVERLAP=42,
1342d643ce63SMatthew Knepley                MATOP_GET_VALUES=43,
1343d643ce63SMatthew Knepley                MATOP_COPY=44,
1344d643ce63SMatthew Knepley                MATOP_PRINT_HELP=45,
1345d643ce63SMatthew Knepley                MATOP_SCALE=46,
1346d643ce63SMatthew Knepley                MATOP_SHIFT=47,
1347d643ce63SMatthew Knepley                MATOP_DIAGONAL_SHIFT=48,
1348d643ce63SMatthew Knepley                MATOP_ILUDT_FACTOR=49,
1349d643ce63SMatthew Knepley                MATOP_GET_BLOCK_SIZE=50,
1350d643ce63SMatthew Knepley                MATOP_GET_ROW_IJ=51,
1351d643ce63SMatthew Knepley                MATOP_RESTORE_ROW_IJ=52,
1352d643ce63SMatthew Knepley                MATOP_GET_COLUMN_IJ=53,
1353d643ce63SMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=54,
1354d643ce63SMatthew Knepley                MATOP_FDCOLORING_CREATE=55,
1355d643ce63SMatthew Knepley                MATOP_COLORING_PATCH=56,
1356d643ce63SMatthew Knepley                MATOP_SET_UNFACTORED=57,
1357d643ce63SMatthew Knepley                MATOP_PERMUTE=58,
1358d643ce63SMatthew Knepley                MATOP_SET_VALUES_BLOCKED=59,
1359d643ce63SMatthew Knepley                MATOP_GET_SUBMATRIX=60,
1360d643ce63SMatthew Knepley                MATOP_DESTROY=61,
1361d643ce63SMatthew Knepley                MATOP_VIEW=62,
1362d643ce63SMatthew Knepley                MATOP_GET_MAPS=63,
1363d643ce63SMatthew Knepley                MATOP_USE_SCALED_FORM=64,
1364d643ce63SMatthew Knepley                MATOP_SCALE_SYSTEM=65,
1365d643ce63SMatthew Knepley                MATOP_UNSCALE_SYSTEM=66,
1366ba16a8cbSKris Buschelman                MATOP_SET_LOCAL_TO_GLOBAL_MAP=67,
1367d643ce63SMatthew Knepley                MATOP_SET_VALUES_LOCAL=68,
1368d643ce63SMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=69,
1369d643ce63SMatthew Knepley                MATOP_GET_ROW_MAX=70,
1370d643ce63SMatthew Knepley                MATOP_CONVERT=71,
1371d643ce63SMatthew Knepley                MATOP_SET_COLORING=72,
1372d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1373d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1374d643ce63SMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1375d643ce63SMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1376ba16a8cbSKris Buschelman                MATOP_MULT_CON=77,
1377ba16a8cbSKris Buschelman                MATOP_MULT_TRANSPOSE_CON=78,
1378ba16a8cbSKris Buschelman                MATOP_ILU_FACTOR_SYMBOLIC_CON=79,
1379d643ce63SMatthew Knepley                MATOP_PERMUTE_SPARSIFY=80,
1380d643ce63SMatthew Knepley                MATOP_MULT_MULTIPLE=81,
138141acf15aSKris Buschelman                MATOP_SOLVE_MULTIPLE=82,
138241acf15aSKris Buschelman                MATOP_GET_INERTIA=83,
138341acf15aSKris Buschelman                MATOP_LOAD=84,
138441acf15aSKris Buschelman                MATOP_IS_SYMMETRIC=85,
138541acf15aSKris Buschelman                MATOP_IS_HERMITIAN=86,
138641acf15aSKris Buschelman                MATOP_IS_STRUCTURALLY_SYMMETRIC=87,
138741acf15aSKris Buschelman                MATOP_PB_RELAX=88,
138841acf15aSKris Buschelman                MATOP_GET_VECS=89,
138941acf15aSKris Buschelman                MATOP_MAT_MULT=90,
139041acf15aSKris Buschelman                MATOP_MAT_MULT_SYMBOLIC=91,
139141acf15aSKris Buschelman                MATOP_MAT_MULT_NUMERIC=92,
139241acf15aSKris Buschelman                MATOP_PTAP=93,
139341acf15aSKris Buschelman                MATOP_PTAP_SYMBOLIC=94,
1394bc011b1eSHong Zhang                MATOP_PTAP_NUMERIC=95,
1395bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE=96,
1396bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97,
13977ba1a0bfSKris Buschelman                MATOP_MAT_MULTTRANSPOSE_NUMERIC=98,
13987ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_SEQAIJ=99,
13997ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_SEQAIJ=100,
14007ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_MPIAIJ=101,
140187d4246cSBarry Smith                MATOP_PTAP_NUMERIC_MPIAIJ=102,
1402f5edf698SHong Zhang                MATOP_SET_VALUES_ROW = 105,
1403d5c63f83SSatish Balay                MATOP_GET_ROW_UTRIANGULAR=108,
1404d5c63f83SSatish Balay                MATOP_RESTORE_ROW_UTRIANGULAR=109
1405fae171e0SBarry Smith              } MatOperation;
1406be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1407be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1408be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1409be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1410112a2221SBarry Smith 
141190ace30eSBarry Smith /*
141290ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
141390ace30eSBarry Smith  stored in a universal format. By changing the format with
1414fb9695e5SSatish Balay  PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will
141590ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
141690ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
141790ace30eSBarry Smith  read into matrices of the same time.
141890ace30eSBarry Smith */
141990ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
142090ace30eSBarry Smith 
1421be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1422be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsGetColor(Mat,ISColoring *);
1423860d1616SSatish Balay 
1424be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
14251f4e1ec7SBarry Smith 
1426d9274352SBarry Smith /*S
1427d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1428d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1429d9274352SBarry Smith 
1430f7a9e4ceSBarry Smith    Level: advanced
1431d9274352SBarry Smith 
1432d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1433d9274352SBarry Smith 
14346e1639daSBarry Smith   Users manual sections:
14356e1639daSBarry Smith .   sec_singular
14366e1639daSBarry Smith 
1437d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1438d9274352SBarry Smith S*/
143974637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1440d9274352SBarry Smith 
1441be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1442281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*);
1443be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1444be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1445be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
1446be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat);
144774637425SBarry Smith 
1448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1449be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1450be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
145146129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
14523f1d51d7SBarry Smith 
1453be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1454be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1455be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1456c4f061fbSSatish Balay 
1457be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1458b0a32e0cSBarry Smith 
1459be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
146004f1ad80SBarry Smith 
1461be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1462be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
14637dbadf16SMatthew Knepley 
1464e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
14652eac72dbSBarry Smith #endif
1466