xref: /petsc/include/petscmat.h (revision b8b68cfd6dee4e3cf13d1d0263895891c1e42174)
12eac72dbSBarry Smith /*
22eac72dbSBarry Smith      Include file for the matrix component of PETSc
32eac72dbSBarry Smith */
40a835dfdSSatish Balay #ifndef __PETSCMAT_H
50a835dfdSSatish Balay #define __PETSCMAT_H
62c8e378dSBarry Smith #include <petscvec.h>
72eac72dbSBarry Smith 
8d9274352SBarry Smith /*S
98f6c3df8SBarry Smith      Mat - Abstract PETSc matrix object used to manage all linear operators in PETSc, even those without
108f6c3df8SBarry Smith            an explicit sparse representation (such as matrix-free operators)
112eac72dbSBarry Smith 
12d91e6319SBarry Smith    Level: beginner
13d91e6319SBarry Smith 
14d9274352SBarry Smith   Concepts: matrix; linear operator
15d9274352SBarry Smith 
168f6c3df8SBarry Smith .seealso:  MatCreate(), MatType, MatSetType(), MatDestroy()
17d9274352SBarry Smith S*/
18d9274352SBarry Smith typedef struct _p_Mat*           Mat;
19d9274352SBarry Smith 
2076bdecfbSBarry Smith /*J
218f6c3df8SBarry Smith     MatType - String with the name of a PETSc matrix type
22d9274352SBarry Smith 
23d9274352SBarry Smith    Level: beginner
24d9274352SBarry Smith 
258f6c3df8SBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage, MatRegister()
2676bdecfbSBarry Smith J*/
2719fd82e9SBarry Smith typedef const char* MatType;
28273d9f13SBarry Smith #define MATSAME            "same"
295a11e1b2SBarry Smith #define MATMAIJ            "maij"
30273d9f13SBarry Smith #define MATSEQMAIJ         "seqmaij"
31273d9f13SBarry Smith #define MATMPIMAIJ         "mpimaij"
32273d9f13SBarry Smith #define MATIS              "is"
335a11e1b2SBarry Smith #define MATAIJ             "aij"
34273d9f13SBarry Smith #define MATSEQAIJ          "seqaij"
35273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
365a11e1b2SBarry Smith #define MATAIJCRL          "aijcrl"
375a11e1b2SBarry Smith #define MATSEQAIJCRL       "seqaijcrl"
385a11e1b2SBarry Smith #define MATMPIAIJCRL       "mpiaijcrl"
398154be41SBarry Smith #define MATAIJCUSP         "aijcusp"
408154be41SBarry Smith #define MATSEQAIJCUSP      "seqaijcusp"
418154be41SBarry Smith #define MATMPIAIJCUSP      "mpiaijcusp"
429ae82921SPaul Mullowney #define MATAIJCUSPARSE     "aijcusparse"
439ae82921SPaul Mullowney #define MATSEQAIJCUSPARSE  "seqaijcusparse"
449ae82921SPaul Mullowney #define MATMPIAIJCUSPARSE  "mpiaijcusparse"
45d67ff14aSKarl Rupp #define MATAIJVIENNACL     "aijviennacl"
46d67ff14aSKarl Rupp #define MATSEQAIJVIENNACL  "seqaijviennacl"
47d67ff14aSKarl Rupp #define MATMPIAIJVIENNACL  "mpiaijviennacl"
485a11e1b2SBarry Smith #define MATAIJPERM         "aijperm"
495a11e1b2SBarry Smith #define MATSEQAIJPERM      "seqaijperm"
505a11e1b2SBarry Smith #define MATMPIAIJPERM      "mpiaijperm"
51273d9f13SBarry Smith #define MATSHELL           "shell"
525a11e1b2SBarry Smith #define MATDENSE           "dense"
53209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
54273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
55db31f6deSJed Brown #define MATELEMENTAL       "elemental"
565a11e1b2SBarry Smith #define MATBAIJ            "baij"
57273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
58273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
59273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
605a11e1b2SBarry Smith #define MATSBAIJ           "sbaij"
61273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
62273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
63cebc7f6cSBarry Smith #define MATDAAD            "daad"
64cebc7f6cSBarry Smith #define MATMFFD            "mffd"
65c8a8475eSBarry Smith #define MATNORMAL          "normal"
66c5e4d11fSDmitry Karpeev #define MATNORMALHERMITIAN "normalh"
67ab92ecdeSBarry Smith #define MATLRC             "lrc"
682a6744ebSBarry Smith #define MATSCATTER         "scatter"
69421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
70793850ffSBarry Smith #define MATCOMPOSITE       "composite"
711f8c7532SHong Zhang #define MATFFT             "fft"
721f8c7532SHong Zhang #define MATFFTW            "fftw"
73e133240eSMatthew G Knepley #define MATSEQCUFFT        "seqcufft"
74557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
7572ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
761d6018f0SLisandro Dalcin #define MATPYTHON          "python"
77f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
78a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT    "hypresstruct"
79ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
802c0dbf93SJed Brown #define MATLOCALREF        "localref"
81d8588912SDave May #define MATNEST            "nest"
82c094ef40SMatthew G. Knepley #define MATPREALLOCATOR    "preallocator"
83773941d6SBarry Smith 
8476bdecfbSBarry Smith /*J
85c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
869989ab13SBarry Smith 
87c2b89b5dSBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu_dist" the parallel SuperLU_DIST package etc
889989ab13SBarry Smith 
899989ab13SBarry Smith    Level: beginner
909989ab13SBarry Smith 
915c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
9276bdecfbSBarry Smith J*/
93c7393fdbSBarry Smith #define MatSolverPackage char*
942692d6eeSBarry Smith #define MATSOLVERSUPERLU      "superlu"
952692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist"
9608f5efcfSPieter Ghysels #define MATSOLVERSTRUMPACK    "strumpack"
972692d6eeSBarry Smith #define MATSOLVERUMFPACK      "umfpack"
9820db9a53SJed Brown #define MATSOLVERCHOLMOD      "cholmod"
99d89f5e7aSBarry Smith #define MATSOLVERKLU          "klu"
100ecddaf3cSBarry Smith #define MATSOLVERCLIQUE       "clique"
101d89f5e7aSBarry Smith #define MATSOLVERELEMENTAL    "elemental"
1022692d6eeSBarry Smith #define MATSOLVERESSL         "essl"
1032692d6eeSBarry Smith #define MATSOLVERLUSOL        "lusol"
1042692d6eeSBarry Smith #define MATSOLVERMUMPS        "mumps"
105d615f992SSatish Balay #define MATSOLVERMKL_PARDISO  "mkl_pardiso"
106d305a81bSVasiliy Kozyrev #define MATSOLVERMKL_CPARDISO "mkl_cpardiso"
1072692d6eeSBarry Smith #define MATSOLVERPASTIX       "pastix"
1082692d6eeSBarry Smith #define MATSOLVERMATLAB       "matlab"
1092692d6eeSBarry Smith #define MATSOLVERPETSC        "petsc"
1102692d6eeSBarry Smith #define MATSOLVERBAS          "bas"
1119ae82921SPaul Mullowney #define MATSOLVERCUSPARSE     "cusparse"
112c0cdd4a1SDahai Guo 
113b24902e0SBarry Smith /*E
114b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
115b24902e0SBarry Smith 
116b24902e0SBarry Smith     Level: beginner
117b24902e0SBarry Smith 
118af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
119b24902e0SBarry Smith 
120c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
121b24902e0SBarry Smith E*/
122599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
123014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[];
124e92e720dSBarry Smith 
125014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
126014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
127014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
128014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*);
12918713533SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister(const MatSolverPackage,const MatType,MatFactorType,PetscErrorCode(*)(Mat,MatFactorType,Mat*));
13042c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageGet(const MatSolverPackage,const MatType,MatFactorType,PetscBool*,PetscBool*,PetscErrorCode (**)(Mat,MatFactorType,Mat*));
1319989ab13SBarry Smith 
132c06d978dSMatthew Knepley /* Logging support */
1330700a824SBarry Smith #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
134014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID;
135335efc43SPeter Brune PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID;
136014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID;
137014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
138014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID;
139014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID;
140014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID;
141014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID;
142c06d978dSMatthew Knepley 
143ceb03754SKris Buschelman /*E
144ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
145d6eff37eSBarry Smith      or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate
146d6eff37eSBarry Smith      that the input matrix is to be replaced with the converted matrix.
147ceb03754SKris Buschelman 
148ceb03754SKris Buschelman     Level: beginner
149ceb03754SKris Buschelman 
150af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
151ceb03754SKris Buschelman 
1520c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
153ceb03754SKris Buschelman E*/
154511c6705SHong Zhang typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse;
155ceb03754SKris Buschelman 
1565494a064SHong Zhang /*E
1575494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
158829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1595494a064SHong Zhang 
1605494a064SHong Zhang     Level: beginner
1615494a064SHong Zhang 
162829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1635494a064SHong Zhang E*/
1645494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1655494a064SHong Zhang 
166607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void);
167c06d978dSMatthew Knepley 
168014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*);
169014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
17019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType);
171014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat);
172685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
173bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat));
174014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]);
175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]);
176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]);
177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]);
17884d44b13SHong Zhang PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat,PetscBool);
179f69a0ea3SMatthew Knepley 
180140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList;
181140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList;
182140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList;
18328988994SBarry Smith 
1843b224e63SBarry Smith /*E
185aa6c7ce3SBarry Smith     MatStructure - Indicates if two matrices have the same nonzero structure
1863b224e63SBarry Smith 
1873b224e63SBarry Smith     Level: beginner
1883b224e63SBarry Smith 
189af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1903b224e63SBarry Smith 
191aa6c7ce3SBarry Smith .seealso: MatCopy(), MatAXPY()
1923b224e63SBarry Smith E*/
193aa6c7ce3SBarry Smith typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure;
1943b224e63SBarry Smith 
195014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
196014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
197014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
198014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
199014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
200014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2018d7a6e47SBarry Smith 
202014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
205d21a29f3SJed Brown 
206014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
208d21a29f3SJed Brown 
209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
21138f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
2134cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]);
214dfb205c3SBarry Smith 
215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
217c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*);
218014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*);
219c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*);
220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
222c0cdd4a1SDahai Guo 
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2306d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
2326d7c1e57SBarry Smith 
23319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*);
234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
235dedccee8SHong Zhang 
236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
237d0de2241SAndrew Spott PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*);
238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*);
239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS);
240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
2411472f72bSBarry Smith 
242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]);
2431d6018f0SLisandro Dalcin 
244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
246e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*);
24721c89e3eSBarry Smith 
248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
253713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
25499cafbc1SBarry Smith 
2558ed539a5SBarry Smith /* ------------------------------------------------------------*/
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
26173a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);
26284cb2905SBarry Smith 
2632ef4de8bSBarry Smith /*S
2642ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
265be479b30SJed Brown         column of a matrix as indexed on an associated grid.
266be479b30SJed Brown 
267be479b30SJed Brown    Fortran usage is different, see MatSetValuesStencil() for details.
2682ef4de8bSBarry Smith 
2692ef4de8bSBarry Smith    Level: beginner
2702ef4de8bSBarry Smith 
2712ef4de8bSBarry Smith   Concepts: matrix; linear operator
2722ef4de8bSBarry Smith 
273be479b30SJed Brown .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil()
2742ef4de8bSBarry Smith S*/
275435da068SBarry Smith typedef struct {
276c1ac3661SBarry Smith   PetscInt k,j,i,c;
277435da068SBarry Smith } MatStencil;
2782ef4de8bSBarry Smith 
279014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
280014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
281014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
282435da068SBarry Smith 
283014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetColoring(Mat,ISColoring);
284014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*);
2853a7fca6bSBarry Smith 
286d91e6319SBarry Smith /*E
287d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
288d91e6319SBarry Smith      to continue to add values to it
289d91e6319SBarry Smith 
290d91e6319SBarry Smith     Level: beginner
291d91e6319SBarry Smith 
292d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
293d91e6319SBarry Smith E*/
2946d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
295014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
296014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
297014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);
2984f9c727eSBarry Smith 
2991d73ed98SBarry Smith 
30030de9b25SBarry Smith 
301d91e6319SBarry Smith /*E
302d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
303d91e6319SBarry Smith 
304d91e6319SBarry Smith     Level: beginner
305d91e6319SBarry Smith 
306af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
307d91e6319SBarry Smith 
308382164b0SBarry Smith    Developer Notes: Entries that are negative need not be called collectively by all processes.
309382164b0SBarry Smith 
310d91e6319SBarry Smith .seealso: MatSetOption()
311d91e6319SBarry Smith E*/
312c5e4d11fSDmitry Karpeev typedef enum {MAT_OPTION_MIN = -3,
313c5e4d11fSDmitry Karpeev               MAT_UNUSED_NONZERO_LOCATION_ERR = -2,
31492d4d306SBarry Smith               MAT_ROW_ORIENTED = -1,
315382164b0SBarry Smith               MAT_SYMMETRIC = 1,
316382164b0SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC = 2,
317382164b0SBarry Smith               MAT_NEW_DIAGONALS = 3,
318382164b0SBarry Smith               MAT_IGNORE_OFF_PROC_ENTRIES = 4,
319382164b0SBarry Smith               MAT_USE_HASH_TABLE = 5,
320382164b0SBarry Smith               MAT_KEEP_NONZERO_PATTERN = 6,
321382164b0SBarry Smith               MAT_IGNORE_ZERO_ENTRIES = 7,
322382164b0SBarry Smith               MAT_USE_INODES = 8,
323382164b0SBarry Smith               MAT_HERMITIAN = 9,
324382164b0SBarry Smith               MAT_SYMMETRY_ETERNAL = 10,
325c5e4d11fSDmitry Karpeev               MAT_NEW_NONZERO_LOCATION_ERR = 11,
326382164b0SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR = 12,
327382164b0SBarry Smith               MAT_ERROR_LOWER_TRIANGULAR = 13,
328382164b0SBarry Smith               MAT_GETROW_UPPERTRIANGULAR = 14,
329382164b0SBarry Smith               MAT_SPD = 15,
33092d4d306SBarry Smith               MAT_NO_OFF_PROC_ZERO_ROWS = 16,
33192d4d306SBarry Smith               MAT_NO_OFF_PROC_ENTRIES = 17,
33292d4d306SBarry Smith               MAT_NEW_NONZERO_LOCATIONS = 18,
333c5e4d11fSDmitry Karpeev               MAT_NEW_NONZERO_ALLOCATION_ERR = 19,
334c5e4d11fSDmitry Karpeev               MAT_SUBSET_OFF_PROC_ENTRIES = 20,
335c5e4d11fSDmitry Karpeev               MAT_OPTION_MAX = 21} MatOption;
336e057df02SPaul Mullowney 
337014dd563SJed Brown PETSC_EXTERN const char *MatOptions[];
338014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool);
3397d68702bSBarry Smith PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*);
34019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);
34184cb2905SBarry Smith 
342014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
343014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
344014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
345014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
346014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
347014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
348014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
349014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
3508c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
3518c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
35221e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*);
353bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3548c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
3558c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
36033d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat);
3611620fd73SBarry Smith 
362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
363014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
364014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
365014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
366014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
367014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
368014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
369014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
374f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);
375f5edf698SHong Zhang 
376d91e6319SBarry Smith /*E
377d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
378d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
379d91e6319SBarry Smith 
380d91e6319SBarry Smith     Level: beginner
381d91e6319SBarry Smith 
382af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
383d91e6319SBarry Smith 
38470dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
38570dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
38670dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
38770dcbbb9SBarry Smith 
388d91e6319SBarry Smith .seealso: MatDuplicate()
389d91e6319SBarry Smith E*/
39070dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
3912e8a6d31SBarry Smith 
39219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
393014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);
39494a9d846SBarry Smith 
395cb5b572fSBarry Smith 
396014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
397014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
398014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
399014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
400014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
401014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
402014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
403014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
404014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
4057b80b807SBarry Smith 
4061a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4071a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
4081a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4091a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
410d4fbbf0eSBarry Smith 
411d91e6319SBarry Smith /*S
412d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
413d91e6319SBarry Smith 
414d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
415d91e6319SBarry Smith 
416d91e6319SBarry Smith    Level: intermediate
417d91e6319SBarry Smith 
418d91e6319SBarry Smith   Concepts: matrix^nonzero information
419d91e6319SBarry Smith 
420d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
421d91e6319SBarry Smith S*/
4224e220ebcSLois Curfman McInnes typedef struct {
423b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
424b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
425b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
426b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
427b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
428b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
429b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4304e220ebcSLois Curfman McInnes } MatInfo;
4314e220ebcSLois Curfman McInnes 
432d9274352SBarry Smith /*E
433d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
434d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
435d9274352SBarry Smith 
436d9274352SBarry Smith     Level: beginner
437d9274352SBarry Smith 
438af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
439d9274352SBarry Smith 
440d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
441d9274352SBarry Smith E*/
4427b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
443014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
444014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
445014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
446014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
447014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
448014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
449014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
450014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
452014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*);
453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);
455a52676ecSHong Zhang 
456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*);
457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*);
458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*);
459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*);
460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*);
461a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
462a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
4637b80b807SBarry Smith 
464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*);
465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*);
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
4737b80b807SBarry Smith 
474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
480566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);
4815ef9f2a5SBarry Smith 
482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
48353dd7562SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);
4907b80b807SBarry Smith 
491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
4988efafbd8SBarry Smith 
499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
500aa1e27eaSFande Kong PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov);
501d2b2a242SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool);
5027b80b807SBarry Smith 
503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat);
50622440eb1SKris Buschelman 
5077bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
5087bab7c10SHong Zhang 
509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat);
512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat);
51522440eb1SKris Buschelman 
516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat);
519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat);
522bc011b1eSHong Zhang 
523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);
5251c741599SHong Zhang 
526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);
5287b80b807SBarry Smith 
529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
531a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*);
532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
538052efed2SBarry Smith 
539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
54190f02eecSBarry Smith 
542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
5452a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*);
54684cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);}
54753cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
5503a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
5519c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
552bd481603SBarry Smith 
553bd481603SBarry Smith /*MC
5541620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5551620fd73SBarry Smith 
5561620fd73SBarry Smith    Not collective
5571620fd73SBarry Smith 
558a9834a7dSSatish Balay    Synopsis:
559a9834a7dSSatish Balay      #include <petscmat.h>
560a9834a7dSSatish Balay      PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode)
561a9834a7dSSatish Balay 
5621620fd73SBarry Smith    Input Parameters:
5631620fd73SBarry Smith +  m - the matrix
5641620fd73SBarry Smith .  row - the row location of the entry
5651620fd73SBarry Smith .  col - the column location of the entry
5661620fd73SBarry Smith .  value - the value to insert
5671620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5681620fd73SBarry Smith 
5691620fd73SBarry Smith    Notes:
5701620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
5711620fd73SBarry Smith    values simultaneously if possible.
5721620fd73SBarry Smith 
5731620fd73SBarry Smith    Level: beginner
5741620fd73SBarry Smith 
5751620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
5761620fd73SBarry Smith M*/
5771620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSetValue(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValues(v,1,&i,1,&j,&va,mode);}
5781620fd73SBarry Smith 
5791620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
5801620fd73SBarry Smith 
5811620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSetValueLocal(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValuesLocal(v,1,&i,1,&j,&va,mode);}
5821620fd73SBarry Smith 
5831620fd73SBarry Smith /*MC
5840d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
585bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
586bd481603SBarry Smith 
587bd481603SBarry Smith    Synopsis:
588aaa7dc30SBarry Smith    #include <petscmat.h>
589c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
590bd481603SBarry Smith 
591bd481603SBarry Smith    Collective on MPI_Comm
592bd481603SBarry Smith 
593bd481603SBarry Smith    Input Parameters:
594bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
595859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
596859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
597bd481603SBarry Smith 
598bd481603SBarry Smith    Output Parameters:
599bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
600bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
601bd481603SBarry Smith 
602bd481603SBarry Smith    Level: intermediate
603bd481603SBarry Smith 
604bd481603SBarry Smith    Notes:
605a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
606bd481603SBarry Smith 
6071d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
608bd481603SBarry Smith 
609bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
610bd481603SBarry Smith 
6111620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6121620fd73SBarry Smith 
613bd481603SBarry Smith   Concepts: preallocation^Matrix
614bd481603SBarry Smith 
615d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
616d6e23781SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock()
617bd481603SBarry Smith M*/
618c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6197c922b88SBarry Smith { \
620a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
6211795a4d1SJed Brown   _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \
6221795a4d1SJed Brown   __start = 0; __end = __start;                                         \
623c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
624a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6257c922b88SBarry Smith 
626bd481603SBarry Smith /*MC
627bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
628bd481603SBarry Smith        inserted using a local number of the rows and columns
629bd481603SBarry Smith 
630bd481603SBarry Smith    Synopsis:
631aaa7dc30SBarry Smith    #include <petscmat.h>
632c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
633bd481603SBarry Smith 
634bd481603SBarry Smith    Not Collective
635bd481603SBarry Smith 
636bd481603SBarry Smith    Input Parameters:
637784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
638bd481603SBarry Smith .  nrows - the number of rows indicated
6391d73ed98SBarry Smith .  rows - the indices of the rows
640784ac674SJed Brown .  cmap - the column mapping from local to global numbering
641bd481603SBarry Smith .  ncols - the number of columns in the matrix
642bd481603SBarry Smith .  cols - the columns indicated
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    Level: intermediate
647bd481603SBarry Smith 
648bd481603SBarry Smith    Notes:
649a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
650bd481603SBarry Smith 
6511d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
652bd481603SBarry Smith 
653bd481603SBarry Smith   Concepts: preallocation^Matrix
654bd481603SBarry Smith 
655d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
656d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
657bd481603SBarry Smith M*/
658784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
659c4f061fbSSatish Balay {\
660c1ac3661SBarry Smith   PetscInt __l;\
661784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
662784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
663c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
664ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
665c4f061fbSSatish Balay   }\
666c4f061fbSSatish Balay }
667c4f061fbSSatish Balay 
668bd481603SBarry Smith /*MC
669d6e23781SBarry Smith    MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
670bd481603SBarry Smith        inserted using a local number of the rows and columns
671bd481603SBarry Smith 
672bd481603SBarry Smith    Synopsis:
673aaa7dc30SBarry Smith    #include <petscmat.h>
674d6e23781SBarry Smith    PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
675d6e23781SBarry Smith 
676d6e23781SBarry Smith    Not Collective
677d6e23781SBarry Smith 
678d6e23781SBarry Smith    Input Parameters:
679d6e23781SBarry Smith +  map - the row mapping from local numbering to global numbering
680d6e23781SBarry Smith .  nrows - the number of rows indicated
681d6e23781SBarry Smith .  rows - the indices of the rows
682d6e23781SBarry Smith .  cmap - the column mapping from local to global numbering
683d6e23781SBarry Smith .  ncols - the number of columns in the matrix
684d6e23781SBarry Smith .  cols - the columns indicated
685d6e23781SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
686d6e23781SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
687d6e23781SBarry Smith 
688d6e23781SBarry Smith    Level: intermediate
689d6e23781SBarry Smith 
690d6e23781SBarry Smith    Notes:
691d6e23781SBarry Smith     See Users-Manual: ch_performance for more details.
692d6e23781SBarry Smith 
693d6e23781SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
694d6e23781SBarry Smith 
695d6e23781SBarry Smith   Concepts: preallocation^Matrix
696d6e23781SBarry Smith 
697d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
698d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
699d6e23781SBarry Smith M*/
700d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
701d6e23781SBarry Smith {\
702d6e23781SBarry Smith   PetscInt __l;\
703d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
704d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
705d6e23781SBarry Smith   for (__l=0;__l<nrows;__l++) {\
706d6e23781SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
707d6e23781SBarry Smith   }\
708d6e23781SBarry Smith }
709d6e23781SBarry Smith 
710d6e23781SBarry Smith /*MC
711d6e23781SBarry Smith    MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
712d6e23781SBarry Smith        inserted using a local number of the rows and columns
713d6e23781SBarry Smith 
714d6e23781SBarry Smith    Synopsis:
715d6e23781SBarry Smith    #include <petscmat.h>
716d6e23781SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
717bd481603SBarry Smith 
718bd481603SBarry Smith    Not Collective
719bd481603SBarry Smith 
720bd481603SBarry Smith    Input Parameters:
721bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
722bd481603SBarry Smith .  nrows - the number of rows indicated
7231d73ed98SBarry Smith .  rows - the indices of the rows
724bd481603SBarry Smith .  ncols - the number of columns in the matrix
725bd481603SBarry Smith .  cols - the columns indicated
726bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
727bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
728bd481603SBarry Smith 
729bd481603SBarry Smith    Level: intermediate
730bd481603SBarry Smith 
731bd481603SBarry Smith    Notes:
732a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
733bd481603SBarry Smith 
734bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
735bd481603SBarry Smith 
736bd481603SBarry Smith   Concepts: preallocation^Matrix
737bd481603SBarry Smith 
738d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(),
739d6e23781SBarry Smith           MatPreallocateInitialize(),  MatPreallocateSetLocal()
740bd481603SBarry Smith M*/
741d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\
742d3d32019SBarry Smith {\
743c1ac3661SBarry Smith   PetscInt __l;\
744d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
745d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
746d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
747d6e23781SBarry Smith     _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
748d3d32019SBarry Smith   }\
749d3d32019SBarry Smith }
750bd481603SBarry Smith /*MC
751bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
752bd481603SBarry Smith        inserted using a local number of the rows and columns
753bd481603SBarry Smith 
754bd481603SBarry Smith    Synopsis:
755aaa7dc30SBarry Smith    #include <petscmat.h>
756c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
757bd481603SBarry Smith 
758bd481603SBarry Smith    Not Collective
759bd481603SBarry Smith 
760bd481603SBarry Smith    Input Parameters:
76164075487SBarry Smith +  row - the row
762bd481603SBarry Smith .  ncols - the number of columns in the matrix
763a50f8bf6SHong Zhang -  cols - the columns indicated
764a50f8bf6SHong Zhang 
765a50f8bf6SHong Zhang    Output Parameters:
766a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
767bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
768bd481603SBarry Smith 
769bd481603SBarry Smith    Level: intermediate
770bd481603SBarry Smith 
771bd481603SBarry Smith    Notes:
772a7f22e61SSatish Balay     See Users-Manual: ch_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 
7761620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7771620fd73SBarry Smith 
778bd481603SBarry Smith   Concepts: preallocation^Matrix
779bd481603SBarry Smith 
780d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
781d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSetLocal()
782bd481603SBarry Smith M*/
783c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
784c1ac3661SBarry Smith { PetscInt __i; \
785e32f2f54SBarry Smith   if (row < __rstart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
786e32f2f54SBarry Smith   if (row >= __rstart+__nrows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D greater than last local row %D",row,__rstart+__nrows-1);\
7877c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
78864075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
7897cc688d7SBarry Smith     else dnz[row - __rstart]++;\
7907c922b88SBarry Smith   }\
7917c922b88SBarry Smith }
7927c922b88SBarry Smith 
793bd481603SBarry Smith /*MC
794d6e23781SBarry Smith    MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
795bd481603SBarry Smith        inserted using a local number of the rows and columns
796bd481603SBarry Smith 
797bd481603SBarry Smith    Synopsis:
798aaa7dc30SBarry Smith    #include <petscmat.h>
799d6e23781SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
800bd481603SBarry Smith 
801bd481603SBarry Smith    Not Collective
802bd481603SBarry Smith 
803bd481603SBarry Smith    Input Parameters:
804bd481603SBarry Smith +  nrows - the number of rows indicated
8051d73ed98SBarry Smith .  rows - the indices of the rows
806bd481603SBarry Smith .  ncols - the number of columns in the matrix
807bd481603SBarry Smith .  cols - the columns indicated
808bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
809bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
810bd481603SBarry Smith 
811bd481603SBarry Smith    Level: intermediate
812bd481603SBarry Smith 
813bd481603SBarry Smith    Notes:
814a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
815bd481603SBarry Smith 
816bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
817bd481603SBarry Smith 
8181620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8191620fd73SBarry Smith 
820bd481603SBarry Smith   Concepts: preallocation^Matrix
821bd481603SBarry Smith 
822d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(),  MatPreallocateInitialize(),
823d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
824bd481603SBarry Smith M*/
825d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\
826c1ac3661SBarry Smith { PetscInt __i; \
827d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
828d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
829d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
830d3d32019SBarry Smith   }\
831d3d32019SBarry Smith }
832d3d32019SBarry Smith 
833bd481603SBarry Smith /*MC
83416371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
83516371a99SBarry Smith 
83616371a99SBarry Smith    Synopsis:
837aaa7dc30SBarry Smith    #include <petscmat.h>
83816371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
83916371a99SBarry Smith 
84016371a99SBarry Smith    Not Collective
84116371a99SBarry Smith 
84216371a99SBarry Smith    Input Parameters:
84316371a99SBarry Smith .  A - matrix
84416371a99SBarry Smith .  row - row where values exist (must be local to this process)
84516371a99SBarry Smith .  ncols - number of columns
84616371a99SBarry Smith .  cols - columns with nonzeros
84716371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
84816371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
84916371a99SBarry Smith 
85016371a99SBarry Smith    Level: intermediate
85116371a99SBarry Smith 
85216371a99SBarry Smith    Notes:
853a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
85416371a99SBarry Smith 
85516371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
85616371a99SBarry Smith 
85716371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
85816371a99SBarry Smith 
85916371a99SBarry Smith   Concepts: preallocation^Matrix
86016371a99SBarry Smith 
861d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
862d6e23781SBarry Smith           MatPreallocateSymmetricSetLocalBlock()
86316371a99SBarry Smith M*/
8640298fd71SBarry Smith #define MatPreallocateLocation(A,row,ncols,cols,dnz,onz) 0;if (A) {ierr = MatSetValues(A,1,&row,ncols,cols,NULL,INSERT_VALUES);CHKERRQ(ierr);} else {ierr =  MatPreallocateSet(row,ncols,cols,dnz,onz);CHKERRQ(ierr);}
86516371a99SBarry Smith 
86616371a99SBarry Smith 
86716371a99SBarry Smith /*MC
8680d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
869bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
870bd481603SBarry Smith 
871bd481603SBarry Smith    Synopsis:
872aaa7dc30SBarry Smith    #include <petscmat.h>
873c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
874bd481603SBarry Smith 
875bd481603SBarry Smith    Collective on MPI_Comm
876bd481603SBarry Smith 
877bd481603SBarry Smith    Input Parameters:
87816371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
879bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
880bd481603SBarry Smith 
881bd481603SBarry Smith    Level: intermediate
882bd481603SBarry Smith 
883bd481603SBarry Smith    Notes:
884a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
885bd481603SBarry Smith 
886bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
887bd481603SBarry Smith 
8881620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
8891620fd73SBarry Smith 
890bd481603SBarry Smith   Concepts: preallocation^Matrix
891bd481603SBarry Smith 
892d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
893d6e23781SBarry Smith           MatPreallocateSymmetricSetLocalBlock()
894bd481603SBarry Smith M*/
895a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
8967c922b88SBarry Smith 
8977b80b807SBarry Smith /* Routines unique to particular data structures */
898014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
899698d4c6aSKris Buschelman 
900014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
901014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
902698d4c6aSKris Buschelman 
903014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
904014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
905014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
906014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
907014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
908014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
9097b80b807SBarry Smith 
910a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
911a96a251dSBarry Smith 
912014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
913014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
914014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
915273d9f13SBarry Smith 
916014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
917014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
918014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
919014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
920014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
921014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
922014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
923014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
924014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
925014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
9269230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
9279230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
928014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
929273d9f13SBarry Smith 
9302e1947a5SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
931014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt);
932014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
9331b807ce4Svictorle 
934014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
935014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
9362e8a6d31SBarry Smith 
937014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
9383a7fca6bSBarry Smith 
939014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
9407b80b807SBarry Smith /*
9417b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
94294b7f48cSBarry Smith   done through the KSP and PC interfaces.
9437b80b807SBarry Smith */
9447b80b807SBarry Smith 
94576bdecfbSBarry Smith /*J
9468f6c3df8SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering
947d9274352SBarry Smith 
948d9274352SBarry Smith    Level: beginner
949d9274352SBarry Smith 
950d9274352SBarry Smith .seealso: MatGetOrdering()
95176bdecfbSBarry Smith J*/
95219fd82e9SBarry Smith typedef const char* MatOrderingType;
9532692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
9542692d6eeSBarry Smith #define MATORDERINGND          "nd"
9552692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
9562692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
9572692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
9582692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
959510184a7SJed Brown #define MATORDERINGWBM         "wbm"
960fc1bef75SJed Brown #define MATORDERINGSPECTRAL    "spectral"
961312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
962b12f92e5SBarry Smith 
96319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
964140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
965bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
966140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList;
967d4fbbf0eSBarry Smith 
968014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
969fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);
970a2ce50c7SBarry Smith 
971d91e6319SBarry Smith /*S
972d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
973d90ac83dSHong Zhang 
974d90ac83dSHong Zhang    Level: beginner
975d90ac83dSHong Zhang 
976d90ac83dSHong Zhang S*/
977d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
9786a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[];
9795e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];
980d90ac83dSHong Zhang 
981d90ac83dSHong Zhang /*S
9829aa7eafdSHong Zhang     MatFactorError - indicates what type of error in matrix factor
9839aa7eafdSHong Zhang 
9849aa7eafdSHong Zhang     Level: beginner
9859aa7eafdSHong Zhang 
98600e125f8SBarry Smith     Developer Notes: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
98700e125f8SBarry Smith 
98800e125f8SBarry Smith .seealso: MatGetFactor()
9899aa7eafdSHong Zhang S*/
9905cd7cf9dSHong Zhang typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError;
9919aa7eafdSHong Zhang 
99200e125f8SBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*);
993*b8b68cfdSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat);
99400e125f8SBarry Smith 
9959aa7eafdSHong Zhang /*S
996422a814eSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization
997b00f7748SHong Zhang 
99861cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
99961cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1000b00f7748SHong Zhang 
100115e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1002b00f7748SHong Zhang 
100361cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
100461cad860SBarry Smith 
1005b00f7748SHong Zhang    Level: developer
1006b00f7748SHong Zhang 
1007d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1008d7d82daaSBarry Smith           MatFactorInfoInitialize()
1009b00f7748SHong Zhang 
1010b00f7748SHong Zhang S*/
1011b00f7748SHong Zhang typedef struct {
101215e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
101385317021SBarry Smith   PetscReal     usedt;
101415e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1015b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
101615e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
101767e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1018348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1019bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1020bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
102115e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1022f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1023f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
102415e8a5b3SHong Zhang } MatFactorInfo;
1025ffa6d0a5SLois Curfman McInnes 
1026014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
10279a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
10289a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
10299a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
10309a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
10319a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
10329a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
10339a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
10349a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
10359a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
10369a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1037014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1038014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1039014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1040014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1041014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1042014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1043014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1044014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
10458e7ba810SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS);
10465a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*);
10475a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*);
10485a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat);
10495a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*);
10505a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec);
10515a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec);
1052014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
1053bb5a7306SBarry Smith 
1054d91e6319SBarry Smith /*E
1055d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1056bb1eb677SSatish Balay 
1057d91e6319SBarry Smith     Level: beginner
1058d91e6319SBarry Smith 
1059d9274352SBarry Smith    May be bitwise ORd together
1060d9274352SBarry Smith 
1061af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1062d91e6319SBarry Smith 
10634e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10644e7234bfSBarry Smith 
106541f059aeSBarry Smith .seealso: MatSOR()
1066d91e6319SBarry Smith E*/
1067ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1068ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1069ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
107084cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1071014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10728ed539a5SBarry Smith 
1073d4fbbf0eSBarry Smith /*
1074639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1075639f9d9dSBarry Smith */
1076b12f92e5SBarry Smith 
10777086a01eSPeter Brune /*S
10787086a01eSPeter Brune      MatColoring - Object for managing the coloring of matrices.
10797086a01eSPeter Brune 
10807086a01eSPeter Brune    Level: beginner
10817086a01eSPeter Brune 
10827086a01eSPeter Brune   Concepts: matrix, coloring
10837086a01eSPeter Brune 
10847086a01eSPeter Brune .seealso:  MatFDColoringCreate() ISColoring MatFDColoring
10857086a01eSPeter Brune S*/
10867086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring;
108776bdecfbSBarry Smith /*J
10888f6c3df8SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring
1089d9274352SBarry Smith 
1090d9274352SBarry Smith    Level: beginner
1091d9274352SBarry Smith 
10927086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring
109376bdecfbSBarry Smith J*/
10947086a01eSPeter Brune 
109519fd82e9SBarry Smith typedef const  char*           MatColoringType;
10965567d71aSPeter Brune #define MATCOLORINGJP      "jp"
1097b9acb617SPeter Brune #define MATCOLORINGPOWER   "power"
10982692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
10992692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
11002692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
11012692d6eeSBarry Smith #define MATCOLORINGID      "id"
11028121bdceSPeter Brune #define MATCOLORINGGREEDY  "greedy"
1103b12f92e5SBarry Smith 
11048ac6da0aSPeter Brune /*E
11058ac6da0aSPeter Brune    MatColoringWeightType - Type of weight scheme
11068ac6da0aSPeter Brune 
11078ac6da0aSPeter Brune     Not Collective
11088ac6da0aSPeter Brune 
11098ac6da0aSPeter Brune +   MAT_COLORING_RANDOM  - Random weights
11108ac6da0aSPeter Brune .   MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
11118ac6da0aSPeter Brune -   MAT_COLORING_LF      - Last-first weighting.
11128ac6da0aSPeter Brune 
11138ac6da0aSPeter Brune     Level: intermediate
11148ac6da0aSPeter Brune 
1115af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
11168ac6da0aSPeter Brune 
11178ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
11188ac6da0aSPeter Brune E*/
1119301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;
11208ac6da0aSPeter Brune 
1121335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1122d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*);
1123335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1124335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1125335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1126335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1127335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1128335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1129335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1130335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1131335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1132335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1133014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
11348ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType);
1135c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*);
11368ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm);
1137cb394dc5SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringTestValid(MatColoring,ISColoring);
1138639f9d9dSBarry Smith 
1139d9274352SBarry Smith /*S
1140d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1141d9274352SBarry Smith         and coloring
1142639f9d9dSBarry Smith 
1143d9274352SBarry Smith    Level: beginner
1144d9274352SBarry Smith 
1145d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1146d9274352SBarry Smith 
1147d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1148d9274352SBarry Smith S*/
1149e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1150639f9d9dSBarry Smith 
1151014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1152014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1153014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1154014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1155014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1156014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1157014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1158d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *);
1159014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1160014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1161f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1162f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1163f86b9fbaSHong Zhang 
1164b1683b59SHong Zhang 
1165b1683b59SHong Zhang /*S
1166b9af6bddSHong Zhang      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1167b1683b59SHong Zhang 
1168b1683b59SHong Zhang    Level: beginner
1169b1683b59SHong Zhang 
1170b1683b59SHong Zhang   Concepts: coloring, sparse matrix product
1171b1683b59SHong Zhang 
1172b9af6bddSHong Zhang .seealso:  MatTransposeColoringCreate()
1173b1683b59SHong Zhang S*/
1174b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1175b1683b59SHong Zhang 
1176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1178014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1179014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1180b1683b59SHong Zhang 
1181639f9d9dSBarry Smith /*
11820752156aSBarry Smith     These routines are for partitioning matrices: currently used only
11833eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
11840752156aSBarry Smith */
1185ca161407SBarry Smith 
1186d9274352SBarry Smith /*S
1187d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1188d9274352SBarry Smith 
1189d9274352SBarry Smith    Level: beginner
1190d9274352SBarry Smith 
1191d9274352SBarry Smith   Concepts: partitioning
1192d9274352SBarry Smith 
1193743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1194d9274352SBarry Smith S*/
119591e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1196d9274352SBarry Smith 
119776bdecfbSBarry Smith /*J
11988f6c3df8SBarry Smith     MatPartitioningType - String with the name of a PETSc matrix partitioning
1199d9274352SBarry Smith 
1200d9274352SBarry Smith    Level: beginner
12010d04baf8SBarry Smith dm
1202b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
120376bdecfbSBarry Smith J*/
120419fd82e9SBarry Smith typedef const char* MatPartitioningType;
12052692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
1206aa1e27eaSFande Kong #define MATPARTITIONINGAVERAGE   "average"
12072692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
12082692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
12092692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
12102692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
121150d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
121288d2ac2bSFande Kong #define MATPARTITIONINGHIERARCH  "hierarch"
121371306c60Sjroman 
1214ca161407SBarry Smith 
1215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
121619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1218014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1219014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
12232aabb6bbSBarry Smith 
1224bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
122530de9b25SBarry Smith 
1226f1144a30SSatish Balay 
12272bad1931SBarry Smith 
1228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
123019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1231ca161407SBarry Smith 
123227538973SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
12350752156aSBarry Smith 
1236b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
12376a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1238b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
12396a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[];
1240b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
12416a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[];
1242b6956c12SJose Roman 
1243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
125471306c60Sjroman 
125571306c60Sjroman #define MP_PARTY_OPT "opt"
125671306c60Sjroman #define MP_PARTY_LIN "lin"
125771306c60Sjroman #define MP_PARTY_SCA "sca"
125871306c60Sjroman #define MP_PARTY_RAN "ran"
125971306c60Sjroman #define MP_PARTY_GBF "gbf"
126071306c60Sjroman #define MP_PARTY_GCF "gcf"
126171306c60Sjroman #define MP_PARTY_BUB "bub"
126271306c60Sjroman #define MP_PARTY_DEF "def"
1263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
126471306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
126571306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
126671306c60Sjroman #define MP_PARTY_NONE "no"
1267014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1268014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
127171306c60Sjroman 
127250d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
12736a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1274e0f1cffaSJose Roman 
1275014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1276014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1277014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1278014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
127971306c60Sjroman 
1280b43b03e9SMark F. Adams /*
12816294fa2bSFande Kong  * hierarchical partitioning
12826294fa2bSFande Kong  */
12834e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*);
12844e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*);
12854e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt);
12864e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);
12876294fa2bSFande Kong 
1288014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1289014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1290591294e4SBarry Smith 
12910752156aSBarry Smith /*
1292af0996ceSBarry Smith     If you add entries here you must also add them to petsc/finclude/petscmat.h
1293d4fbbf0eSBarry Smith */
12941c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
12951c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
12961c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
12971c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
12981c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
12997c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13007c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13011c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13021c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13037c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13047c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13051c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13061c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1307a714c06dSBarry Smith                MATOP_SOR=13,
13081c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13091c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13101c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13111c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13121c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13131c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13141c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13151c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1316d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1317d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1318d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1319d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1320d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1321d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1322d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1323d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1324d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1325d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
13269d39f69dSJed Brown                /* MATOP_PLACEHOLDER_32=32, */
13279d39f69dSJed Brown                /* MATOP_PLACEHOLDER_33=33, */
1328d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1329d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1330d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1331d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1332d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1333d519adbfSMatthew Knepley                MATOP_AXPY=39,
1334d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1335d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1336d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1337d519adbfSMatthew Knepley                MATOP_COPY=43,
1338d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1339d519adbfSMatthew Knepley                MATOP_SCALE=45,
1340d519adbfSMatthew Knepley                MATOP_SHIFT=46,
134135153367SBarry Smith                MATOP_DIAGONAL_SET=47,
13429d39f69dSJed Brown                MATOP_ZERO_ROWS_COLUMNS=48,
13439d39f69dSJed Brown                MATOP_SET_RANDOM=49,
1344d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1345d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1346d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1347d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1348d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1349d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1350d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1351d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1352d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1353d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1354d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1355d519adbfSMatthew Knepley                MATOP_VIEW=61,
1356d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
13577bab7c10SHong Zhang                MATOP_MATMAT_MULT=63,
13587bab7c10SHong Zhang                MATOP_MATMAT_MULT_SYMBOLIC=64,
13597bab7c10SHong Zhang                MATOP_MATMAT_MULT_NUMERIC=65,
1360d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1361d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1362d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1363d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1364d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1365d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1366d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
13679d39f69dSJed Brown                /* MATOP_PLACEHOLDER_73=73, */
1368d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1369d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1370d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1371bda6c4cbSBarry Smith                MATOP_MULT_CONSTRAINED=77,
1372bda6c4cbSBarry Smith                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
13739d39f69dSJed Brown                MATOP_FIND_ZERO_DIAGONALS=79,
1374d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1375d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1376d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1377d519adbfSMatthew Knepley                MATOP_LOAD=83,
1378d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1379d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1380d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1381820469bcSHong Zhang                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1382d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1383d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1384d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1385d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1386d519adbfSMatthew Knepley                MATOP_PTAP=92,
1387d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1388d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1389bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT=95,
1390bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1391bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
13929d39f69dSJed Brown                /* MATOP_PLACEHOLDER_98=98, */
13939d39f69dSJed Brown                /* MATOP_PLACEHOLDER_99=99, */
13949d39f69dSJed Brown                /* MATOP_PLACEHOLDER_100=100, */
13959d39f69dSJed Brown                /* MATOP_PLACEHOLDER_101=101, */
1396d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
13979d39f69dSJed Brown                /* MATOP_PLACEHOLDER_103=103, */
1398d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1399d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1400bda6c4cbSBarry Smith                MATOP_IMAGINARY_PART=106,
1401bda6c4cbSBarry Smith                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1402bda6c4cbSBarry Smith                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1403bda6c4cbSBarry Smith                MATOP_MAT_SOLVE=109,
14040e8d04f7SBarry Smith                MATOP_GET_REDUNDANT_MATRIX=110,
1405d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
14060e8d04f7SBarry Smith                MATOP_GET_COLUMN_VECTOR=112,
1407d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
14080e8d04f7SBarry Smith                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
140989c6957cSBarry Smith                MATOP_CREATE=115,
141089c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
14110e8d04f7SBarry Smith                MATOP_GET_LOCAL_SUB_MATRIX=117,
14120e8d04f7SBarry Smith                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1413eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
14140e8d04f7SBarry Smith                MATOP_HERMITIAN_TRANSPOSE=120,
14150e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
14160e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
14170e8d04f7SBarry Smith                MATOP_GET_MULTI_PROC_BLOCK=123,
14189d39f69dSJed Brown                MATOP_FIND_NONZERO_ROWS=124,
14190e8d04f7SBarry Smith                MATOP_GET_COLUMN_NORMS=125,
14209d39f69dSJed Brown                MATOP_INVERT_BLOCK_DIAGONAL=126,
14219d39f69dSJed Brown                /* MATOP_PLACEHOLDER_127=127, */
14220e8d04f7SBarry Smith                MATOP_GET_SUB_MATRICES_PARALLE=128,
14232b8ad9a3SHong Zhang                MATOP_SET_VALUES_BATCH=129,
14240e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT=130,
14250e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
1426e9b602ebSSatish Balay                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
14270e8d04f7SBarry Smith                MATOP_TRANSPOSE_COLORING_CREAT=133,
14280e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_SPT=134,
14290e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_DEN=135,
14300e8d04f7SBarry Smith                MATOP_RART=136,
14310e8d04f7SBarry Smith                MATOP_RART_SYMBOLIC=137,
14320e8d04f7SBarry Smith                MATOP_RART_NUMERIC=138,
1433e09a3074SHong Zhang                MATOP_SET_BLOCK_SIZES=139,
1434f9426fe0SMark Adams                MATOP_AYPX=140,
14351919a2e2SJed Brown                MATOP_RESIDUAL=141,
14369c8f2541SHong Zhang                MATOP_FDCOLORING_SETUP=142,
14379c8f2541SHong Zhang                MATOP_MPICONCATENATESEQ=144
1438fae171e0SBarry Smith              } MatOperation;
1439014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *);
1440014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1441014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1442014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1443112a2221SBarry Smith 
144490ace30eSBarry Smith /*
144590ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
144690ace30eSBarry Smith    stored in a universal format. By changing the format with
14476a9046bcSBarry Smith    PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
144890ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
144990ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1450f4403165SShri Abhyankar    read into matrices of the same type.
145190ace30eSBarry Smith */
145290ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
145390ace30eSBarry Smith 
1454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
1457b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*);
14581f4e1ec7SBarry Smith 
1459d9274352SBarry Smith /*S
1460d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1461d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1462d9274352SBarry Smith 
1463f7a9e4ceSBarry Smith    Level: advanced
1464d9274352SBarry Smith 
1465d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1466d9274352SBarry Smith 
14676e1639daSBarry Smith   Users manual sections:
14686e1639daSBarry Smith .   sec_singular
14696e1639daSBarry Smith 
1470d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1471d9274352SBarry Smith S*/
147274637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1473d9274352SBarry Smith 
1474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1477d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
14795fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *);
14805fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace);
1481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
148874637425SBarry Smith 
1489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
14933f1d51d7SBarry Smith 
1494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1497c4f061fbSSatish Balay 
1498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*);
1499b0a32e0cSBarry Smith 
1500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
150104f1ad80SBarry Smith 
1502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1515e884886fSBarry Smith 
15166370053bSBarry Smith /*S
15176370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
15186370053bSBarry Smith               Jacobian vector products
1519e884886fSBarry Smith 
15206370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
15216370053bSBarry Smith 
15226370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
15236370053bSBarry Smith 
15246370053bSBarry Smith     Level: developer
15256370053bSBarry Smith 
15266370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
15276370053bSBarry Smith S*/
1528e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1529e884886fSBarry Smith 
153076bdecfbSBarry Smith /*J
1531e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1532e884886fSBarry Smith 
1533e884886fSBarry Smith    Level: beginner
1534e884886fSBarry Smith 
1535e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
153676bdecfbSBarry Smith J*/
153719fd82e9SBarry Smith typedef const char* MatMFFDType;
1538e884886fSBarry Smith #define MATMFFD_DS  "ds"
1539e884886fSBarry Smith #define MATMFFD_WP  "wp"
1540e884886fSBarry Smith 
154119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1542bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));
1543e884886fSBarry Smith 
1544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool );
1546e884886fSBarry Smith 
154761ab5df0SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType);
154861ab5df0SBarry Smith 
1549014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1550014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
15517dbadf16SMatthew Knepley 
155297969023SHong Zhang /*
155397969023SHong Zhang    PETSc interface to MUMPS
155497969023SHong Zhang */
155597969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1557bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*);
1558b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
1559bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*);
1560bc6112feSHong Zhang 
1561ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1562ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1563ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1564ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
156597969023SHong Zhang #endif
156623a5497aSJed Brown 
1567d954db57SHong Zhang /*
1568d615f992SSatish Balay    PETSc interface to Mkl_Pardiso
1569b500a6b7SJose David Bermeol */
1570b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO
1571d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1572b500a6b7SJose David Bermeol #endif
1573b500a6b7SJose David Bermeol 
1574b500a6b7SJose David Bermeol /*
1575d305a81bSVasiliy Kozyrev    PETSc interface to Mkl_CPardiso
1576d305a81bSVasiliy Kozyrev */
1577d305a81bSVasiliy Kozyrev #ifdef PETSC_HAVE_MKL_CPARDISO
1578d305a81bSVasiliy Kozyrev PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt);
1579d305a81bSVasiliy Kozyrev #endif
1580d305a81bSVasiliy Kozyrev 
1581d305a81bSVasiliy Kozyrev /*
1582d954db57SHong Zhang    PETSc interface to SUPERLU
1583d954db57SHong Zhang */
1584d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1585014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1586d954db57SHong Zhang #endif
1587d954db57SHong Zhang 
1588fb785984SHong Zhang /*
1589fb785984SHong Zhang    PETSc interface to SUPERLU_DIST
1590fb785984SHong Zhang */
1591fb785984SHong Zhang #ifdef PETSC_HAVE_SUPERLU_DIST
1592fb785984SHong Zhang PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*);
1593fb785984SHong Zhang #endif
1594fb785984SHong Zhang 
1595575704cbSPieter Ghysels /*
1596575704cbSPieter Ghysels    PETSc interface to STRUMPACK
1597575704cbSPieter Ghysels */
1598575704cbSPieter Ghysels #ifdef PETSC_HAVE_STRUMPACK
1599575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat,PetscReal);
1600575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat,PetscInt);
1601575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool);
1602575704cbSPieter Ghysels #endif
1603575704cbSPieter Ghysels 
1604575704cbSPieter Ghysels 
1605aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA
16063f9c0db1SPaul Mullowney /*E
1607e057df02SPaul Mullowney     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
16082692e278SPaul Mullowney     matrices.
1609e057df02SPaul Mullowney 
1610e057df02SPaul Mullowney     Not Collective
1611e057df02SPaul Mullowney 
16128468deeeSKarl Rupp +   MAT_CUSPARSE_CSR - Compressed Sparse Row
16132692e278SPaul Mullowney .   MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later).
16142692e278SPaul Mullowney -   MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).
1615e057df02SPaul Mullowney 
1616e057df02SPaul Mullowney     Level: intermediate
1617e057df02SPaul Mullowney 
1618af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1619e057df02SPaul Mullowney 
1620e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1621e057df02SPaul Mullowney E*/
1622e057df02SPaul Mullowney 
16233f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;
1624e057df02SPaul Mullowney 
1625e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1626e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
1627e057df02SPaul Mullowney 
16283f9c0db1SPaul Mullowney /*E
1629e057df02SPaul Mullowney     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
16302692e278SPaul Mullowney     matrices whose operation should use a particular storage format.
1631e057df02SPaul Mullowney 
1632e057df02SPaul Mullowney     Not Collective
1633e057df02SPaul Mullowney 
16348468deeeSKarl Rupp +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
16358468deeeSKarl Rupp .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
16368468deeeSKarl Rupp .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
16378468deeeSKarl Rupp -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices
1638e057df02SPaul Mullowney 
1639e057df02SPaul Mullowney     Level: intermediate
1640e057df02SPaul Mullowney 
1641e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1642e057df02SPaul Mullowney E*/
164336d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;
1644e057df02SPaul Mullowney 
164510e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
164610e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1647e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
16489ae82921SPaul Mullowney #endif
16499ae82921SPaul Mullowney 
165090c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
1651014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1652014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1653e057df02SPaul Mullowney 
16543f9c0db1SPaul Mullowney /*E
1655e057df02SPaul Mullowney     MatCUSPStorageFormat - indicates the storage format for CUSP (GPU)
165636d62e41SPaul Mullowney     matrices.
1657e057df02SPaul Mullowney 
1658e057df02SPaul Mullowney     Not Collective
1659e057df02SPaul Mullowney 
16608468deeeSKarl Rupp +   MAT_CUSP_CSR - Compressed Sparse Row
16618468deeeSKarl Rupp .   MAT_CUSP_DIA - Diagonal
16628468deeeSKarl Rupp -   MAT_CUSP_ELL - Ellpack
1663e057df02SPaul Mullowney 
1664e057df02SPaul Mullowney     Level: intermediate
1665e057df02SPaul Mullowney 
1666af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1667e057df02SPaul Mullowney 
1668e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
1669e057df02SPaul Mullowney E*/
16703f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat;
1671e057df02SPaul Mullowney 
1672e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1673e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[];
1674e057df02SPaul Mullowney 
16753f9c0db1SPaul Mullowney /*E
1676e057df02SPaul Mullowney     MatCUSPFormatOperation - indicates the operation of CUSP (GPU)
167736d62e41SPaul Mullowney     matrices whose operation should use a particular storage format.
1678e057df02SPaul Mullowney 
1679e057df02SPaul Mullowney     Not Collective
1680e057df02SPaul Mullowney 
16818468deeeSKarl Rupp +   MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
16828468deeeSKarl Rupp .   MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
16838468deeeSKarl Rupp .   MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
16848468deeeSKarl Rupp -   MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices
1685e057df02SPaul Mullowney 
1686e057df02SPaul Mullowney     Level: intermediate
1687e057df02SPaul Mullowney 
1688af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1689e057df02SPaul Mullowney 
1690e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat
1691e057df02SPaul Mullowney E*/
169236d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation;
1693e057df02SPaul Mullowney 
1694e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat);
1695e057df02SPaul Mullowney #endif
169690c53902SBarry Smith 
1697d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL)
1698d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1699d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1700d67ff14aSKarl Rupp #endif
1701d67ff14aSKarl Rupp 
170254efbe56SHong Zhang /*
170354efbe56SHong Zhang    PETSc interface to FFTW
170454efbe56SHong Zhang */
170554efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
1706014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1707014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
17082a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*);
170954efbe56SHong Zhang #endif
171054efbe56SHong Zhang 
1711382fd914SHong Zhang /*
1712382fd914SHong Zhang    PETSc interface to ELEMENTAL
1713382fd914SHong Zhang */
1714db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL)
17152ef15aa2SHong Zhang #if defined(PETSC_USE_COMPLEX)
17162ef15aa2SHong Zhang typedef El::Complex<PetscReal> PetscElemScalar;
17172ef15aa2SHong Zhang #else
17182ef15aa2SHong Zhang typedef PetscScalar PetscElemScalar;
1719382fd914SHong Zhang #endif
1720607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
1721db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
1722db31f6deSJed Brown #endif
1723db31f6deSJed Brown 
1724014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1725014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1726014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1727014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1728014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1729014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
173019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1731014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1732014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1733d8588912SDave May 
17344325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
1735e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);
17364325cce7SMatthew G Knepley 
1737b541e6a4SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**);
173853dd7562SDmitry Karpeev 
1739c094ef40SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat);
1740c094ef40SMatthew G. Knepley 
1741539c167fSBarry Smith PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*);
1742539c167fSBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*);
1743539c167fSBarry Smith 
174423a5497aSJed Brown #endif
1745