xref: /petsc/include/petscmat.h (revision f9426fe092dba0ba2fdf65dfec8d938c4b10a31c)
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"
357d6a0e61SBarry Smith #define MATSEQAIJPTHREAD   "seqaijpthread"
36bf2c1783SBarry Smith #define MATAIJPTHREAD      "aijpthread"
37273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
385a11e1b2SBarry Smith #define MATAIJCRL          "aijcrl"
395a11e1b2SBarry Smith #define MATSEQAIJCRL       "seqaijcrl"
405a11e1b2SBarry Smith #define MATMPIAIJCRL       "mpiaijcrl"
418154be41SBarry Smith #define MATAIJCUSP         "aijcusp"
428154be41SBarry Smith #define MATSEQAIJCUSP      "seqaijcusp"
438154be41SBarry Smith #define MATMPIAIJCUSP      "mpiaijcusp"
449ae82921SPaul Mullowney #define MATAIJCUSPARSE     "aijcusparse"
459ae82921SPaul Mullowney #define MATSEQAIJCUSPARSE  "seqaijcusparse"
469ae82921SPaul Mullowney #define MATMPIAIJCUSPARSE  "mpiaijcusparse"
47d67ff14aSKarl Rupp #define MATAIJVIENNACL     "aijviennacl"
48d67ff14aSKarl Rupp #define MATSEQAIJVIENNACL  "seqaijviennacl"
49d67ff14aSKarl Rupp #define MATMPIAIJVIENNACL  "mpiaijviennacl"
505a11e1b2SBarry Smith #define MATAIJPERM         "aijperm"
515a11e1b2SBarry Smith #define MATSEQAIJPERM      "seqaijperm"
525a11e1b2SBarry Smith #define MATMPIAIJPERM      "mpiaijperm"
53273d9f13SBarry Smith #define MATSHELL           "shell"
545a11e1b2SBarry Smith #define MATDENSE           "dense"
55209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
56273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
57db31f6deSJed Brown #define MATELEMENTAL       "elemental"
585a11e1b2SBarry Smith #define MATBAIJ            "baij"
59273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
60273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
61273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
625a11e1b2SBarry Smith #define MATSBAIJ           "sbaij"
63273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
64273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
65c0cdd4a1SDahai Guo #define MATSEQBSTRM        "seqbstrm"
66c0cdd4a1SDahai Guo #define MATMPIBSTRM        "mpibstrm"
67c0cdd4a1SDahai Guo #define MATBSTRM           "bstrm"
68aa5a9175SDahai Guo #define MATSEQSBSTRM       "seqsbstrm"
69aa5a9175SDahai Guo #define MATMPISBSTRM       "mpisbstrm"
70aa5a9175SDahai Guo #define MATSBSTRM          "sbstrm"
71cebc7f6cSBarry Smith #define MATDAAD            "daad"
72cebc7f6cSBarry Smith #define MATMFFD            "mffd"
73c8a8475eSBarry Smith #define MATNORMAL          "normal"
74ab92ecdeSBarry Smith #define MATLRC             "lrc"
752a6744ebSBarry Smith #define MATSCATTER         "scatter"
76421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
77793850ffSBarry Smith #define MATCOMPOSITE       "composite"
781f8c7532SHong Zhang #define MATFFT             "fft"
791f8c7532SHong Zhang #define MATFFTW            "fftw"
80e133240eSMatthew G Knepley #define MATSEQCUFFT        "seqcufft"
81557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
8272ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
831d6018f0SLisandro Dalcin #define MATPYTHON          "python"
84f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
85a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT    "hypresstruct"
86ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
872c0dbf93SJed Brown #define MATLOCALREF        "localref"
88d8588912SDave May #define MATNEST            "nest"
89773941d6SBarry Smith 
9076bdecfbSBarry Smith /*J
91c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
929989ab13SBarry Smith 
939989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
949989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
959989ab13SBarry Smith 
969989ab13SBarry Smith 
979989ab13SBarry Smith    Level: beginner
989989ab13SBarry Smith 
995c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
10076bdecfbSBarry Smith J*/
101c7393fdbSBarry Smith #define MatSolverPackage char*
1022692d6eeSBarry Smith #define MATSOLVERSUPERLU      "superlu"
1032692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist"
1042692d6eeSBarry Smith #define MATSOLVERUMFPACK      "umfpack"
10520db9a53SJed Brown #define MATSOLVERCHOLMOD      "cholmod"
1062692d6eeSBarry Smith #define MATSOLVERESSL         "essl"
1072692d6eeSBarry Smith #define MATSOLVERLUSOL        "lusol"
1082692d6eeSBarry Smith #define MATSOLVERMUMPS        "mumps"
1092692d6eeSBarry Smith #define MATSOLVERPASTIX       "pastix"
1102692d6eeSBarry Smith #define MATSOLVERMATLAB       "matlab"
1112692d6eeSBarry Smith #define MATSOLVERPETSC        "petsc"
1122692d6eeSBarry Smith #define MATSOLVERBAS          "bas"
1139ae82921SPaul Mullowney #define MATSOLVERCUSPARSE     "cusparse"
114aa5a9175SDahai Guo #define MATSOLVERBSTRM        "bstrm"
115aa5a9175SDahai Guo #define MATSOLVERSBSTRM       "sbstrm"
11615767789SHong Zhang #define MATSOLVERELEMENTAL    "elemental"
11717f1a0eaSHong Zhang #define MATSOLVERCLIQUE       "clique"
118c0cdd4a1SDahai Guo 
119b24902e0SBarry Smith /*E
120b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
121b24902e0SBarry Smith 
122b24902e0SBarry Smith     Level: beginner
123b24902e0SBarry Smith 
124b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
125b24902e0SBarry Smith 
126c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
127b24902e0SBarry Smith E*/
128599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
129014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[];
130e92e720dSBarry Smith 
131014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*);
1359989ab13SBarry Smith 
136c06d978dSMatthew Knepley /* Logging support */
1370700a824SBarry Smith #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
138014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID;
139014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID;
140014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
141014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID;
142014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID;
143014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID;
144014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID;
145c06d978dSMatthew Knepley 
146ceb03754SKris Buschelman /*E
147ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
148d6eff37eSBarry Smith      or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate
149d6eff37eSBarry Smith      that the input matrix is to be replaced with the converted matrix.
150ceb03754SKris Buschelman 
151ceb03754SKris Buschelman     Level: beginner
152ceb03754SKris Buschelman 
153ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
154ceb03754SKris Buschelman 
1550c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
156ceb03754SKris Buschelman E*/
157dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse;
158ceb03754SKris Buschelman 
1595494a064SHong Zhang /*E
1605494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
161829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1625494a064SHong Zhang 
1635494a064SHong Zhang     Level: beginner
1645494a064SHong Zhang 
165829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1665494a064SHong Zhang E*/
1675494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1685494a064SHong Zhang 
169607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void);
170c06d978dSMatthew Knepley 
171014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*);
172014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
17319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType);
174014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat);
175146574abSBarry Smith PETSC_EXTERN PetscErrorCode MatViewFromOptions(Mat,const char[],const char[]);
176607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatRegisterAll(void);
177bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat));
178014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]);
179014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]);
180014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]);
181014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]);
182f69a0ea3SMatthew Knepley 
183014dd563SJed Brown PETSC_EXTERN PetscBool         MatRegisterAllCalled;
184140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList;
185140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList;
186140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList;
187140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatCoarsenList;
18828988994SBarry Smith 
1893b224e63SBarry Smith /*E
1903b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
1913b224e63SBarry Smith 
1923b224e63SBarry Smith     Level: beginner
1933b224e63SBarry Smith 
1943b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1953b224e63SBarry Smith 
1963b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
1973b224e63SBarry Smith E*/
198214bc6a2SJed Brown typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN,SAME_PRECONDITIONER} MatStructure;
1993b224e63SBarry Smith 
200014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
201014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
202014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2068d7a6e47SBarry Smith 
207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
208014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
210d21a29f3SJed Brown 
211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
213d21a29f3SJed Brown 
214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*);
218dfb205c3SBarry Smith 
219014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*);
222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
225c0cdd4a1SDahai Guo 
226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
230c0cdd4a1SDahai Guo 
231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2386d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
2406d7c1e57SBarry Smith 
24119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*);
242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
243dedccee8SHong Zhang 
244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*);
246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS);
247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
2481472f72bSBarry Smith 
249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]);
2501d6018f0SLisandro Dalcin 
251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
25321c89e3eSBarry Smith 
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
259713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
26099cafbc1SBarry Smith 
2618ed539a5SBarry Smith /* ------------------------------------------------------------*/
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
26773a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);
26884cb2905SBarry Smith 
2692ef4de8bSBarry Smith /*S
2702ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
271be479b30SJed Brown         column of a matrix as indexed on an associated grid.
272be479b30SJed Brown 
273be479b30SJed Brown    Fortran usage is different, see MatSetValuesStencil() for details.
2742ef4de8bSBarry Smith 
2752ef4de8bSBarry Smith    Level: beginner
2762ef4de8bSBarry Smith 
2772ef4de8bSBarry Smith   Concepts: matrix; linear operator
2782ef4de8bSBarry Smith 
279be479b30SJed Brown .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil()
2802ef4de8bSBarry Smith S*/
281435da068SBarry Smith typedef struct {
282c1ac3661SBarry Smith   PetscInt k,j,i,c;
283435da068SBarry Smith } MatStencil;
2842ef4de8bSBarry Smith 
285014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
286014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
287014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
288435da068SBarry Smith 
289014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetColoring(Mat,ISColoring);
290014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*);
2913a7fca6bSBarry Smith 
292d91e6319SBarry Smith /*E
293d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
294d91e6319SBarry Smith      to continue to add values to it
295d91e6319SBarry Smith 
296d91e6319SBarry Smith     Level: beginner
297d91e6319SBarry Smith 
298d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
299d91e6319SBarry Smith E*/
3006d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
301014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
302014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
303014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);
3044f9c727eSBarry Smith 
3051d73ed98SBarry Smith 
30630de9b25SBarry Smith 
307d91e6319SBarry Smith /*E
308d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
309d91e6319SBarry Smith 
310d91e6319SBarry Smith     Level: beginner
311d91e6319SBarry Smith 
3120a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
313d91e6319SBarry Smith 
314382164b0SBarry Smith    Developer Notes: Entries that are negative need not be called collectively by all processes.
315382164b0SBarry Smith 
316d91e6319SBarry Smith .seealso: MatSetOption()
317d91e6319SBarry Smith E*/
318382164b0SBarry Smith typedef enum {MAT_OPTION_MIN = -8,
319382164b0SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR = -7,
320382164b0SBarry Smith               MAT_NO_OFF_PROC_ZERO_ROWS = -6,
321382164b0SBarry Smith               MAT_NO_OFF_PROC_ENTRIES = -5,
322382164b0SBarry Smith               MAT_UNUSED_NONZERO_LOCATION_ERR = -4,
323382164b0SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR = -3,
324382164b0SBarry Smith               MAT_ROW_ORIENTED = -2,
325382164b0SBarry Smith               MAT_NEW_NONZERO_LOCATIONS = -1,
326382164b0SBarry Smith               MAT_SYMMETRIC = 1,
327382164b0SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC = 2,
328382164b0SBarry Smith               MAT_NEW_DIAGONALS = 3,
329382164b0SBarry Smith               MAT_IGNORE_OFF_PROC_ENTRIES = 4,
330382164b0SBarry Smith               MAT_USE_HASH_TABLE = 5,
331382164b0SBarry Smith               MAT_KEEP_NONZERO_PATTERN = 6,
332382164b0SBarry Smith               MAT_IGNORE_ZERO_ENTRIES = 7,
333382164b0SBarry Smith               MAT_USE_INODES = 8,
334382164b0SBarry Smith               MAT_HERMITIAN = 9,
335382164b0SBarry Smith               MAT_SYMMETRY_ETERNAL = 10,
336382164b0SBarry Smith               MAT_CHECK_COMPRESSED_ROW = 11,
337382164b0SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR = 12,
338382164b0SBarry Smith               MAT_ERROR_LOWER_TRIANGULAR = 13,
339382164b0SBarry Smith               MAT_GETROW_UPPERTRIANGULAR = 14,
340382164b0SBarry Smith               MAT_SPD = 15,
341382164b0SBarry Smith               MAT_OPTION_MAX = 16} MatOption;
342e057df02SPaul Mullowney 
343014dd563SJed Brown PETSC_EXTERN const char *MatOptions[];
344014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool );
34519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);
34684cb2905SBarry Smith 
347014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
348014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
349014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
350014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
351014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
3558c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
3568c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
3578c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
3588c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
363014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNThreads(Mat,PetscInt);
364014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNThreads(Mat,PetscInt*);
3651620fd73SBarry Smith 
366014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
367014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
368014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
369014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
378*f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);
379f5edf698SHong Zhang 
380d91e6319SBarry Smith /*E
381d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
382d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
383d91e6319SBarry Smith 
384d91e6319SBarry Smith     Level: beginner
385d91e6319SBarry Smith 
386d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
387d91e6319SBarry Smith 
38870dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
38970dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
39070dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
39170dcbbb9SBarry Smith 
392d91e6319SBarry Smith .seealso: MatDuplicate()
393d91e6319SBarry Smith E*/
39470dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
3952e8a6d31SBarry Smith 
39619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
397014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);
39894a9d846SBarry Smith 
399cb5b572fSBarry Smith 
400014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
401014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
402014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
403014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
404014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
405014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
406014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
407014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
4097b80b807SBarry Smith 
4101a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4111a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
4121a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4131a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
414d4fbbf0eSBarry Smith 
415d91e6319SBarry Smith /*S
416d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
417d91e6319SBarry Smith 
418d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
419d91e6319SBarry Smith 
420d91e6319SBarry Smith    Level: intermediate
421d91e6319SBarry Smith 
422d91e6319SBarry Smith   Concepts: matrix^nonzero information
423d91e6319SBarry Smith 
424d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
425d91e6319SBarry Smith S*/
4264e220ebcSLois Curfman McInnes typedef struct {
427b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
428b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
429b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
430b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
431b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
432b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
433b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4344e220ebcSLois Curfman McInnes } MatInfo;
4354e220ebcSLois Curfman McInnes 
436d9274352SBarry Smith /*E
437d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
438d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
439d9274352SBarry Smith 
440d9274352SBarry Smith     Level: beginner
441d9274352SBarry Smith 
442d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
443d9274352SBarry Smith 
444d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
445d9274352SBarry Smith E*/
4467b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
447014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
448014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
449014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
450014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
452014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat *);
457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);
459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool *);
460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
4647b80b807SBarry Smith 
465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal *);
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *);
467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
4747b80b807SBarry Smith 
475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
481566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);
4825ef9f2a5SBarry Smith 
483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);
4917b80b807SBarry Smith 
492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,Mat*);
494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJNumeric(MPI_Comm,Mat,PetscInt,Mat);
495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5028efafbd8SBarry Smith 
503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5047b80b807SBarry Smith 
505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat);
50822440eb1SKris Buschelman 
5097bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
5107bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultSymbolic(Mat,Mat,Mat,PetscReal,Mat*);
5117bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultNumeric(Mat,Mat,Mat,Mat);
5127bab7c10SHong Zhang 
513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat);
516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat);
51922440eb1SKris Buschelman 
520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat);
523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat);
526bc011b1eSHong Zhang 
527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);
5291c741599SHong Zhang 
530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);
5327b80b807SBarry Smith 
533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
543052efed2SBarry Smith 
544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
54690f02eecSBarry Smith 
547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecs(Mat,Vec*,Vec*);
551b2bf6370SHong Zhang PETSC_EXTERN PetscErrorCode MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
554bd481603SBarry Smith 
555bd481603SBarry Smith /*MC
5561620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5571620fd73SBarry Smith 
5581620fd73SBarry Smith    Not collective
5591620fd73SBarry Smith 
5601620fd73SBarry Smith    Input Parameters:
5611620fd73SBarry Smith +  m - the matrix
5621620fd73SBarry Smith .  row - the row location of the entry
5631620fd73SBarry Smith .  col - the column location of the entry
5641620fd73SBarry Smith .  value - the value to insert
5651620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5661620fd73SBarry Smith 
5671620fd73SBarry Smith    Notes:
5681620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
5691620fd73SBarry Smith    values simultaneously if possible.
5701620fd73SBarry Smith 
5711620fd73SBarry Smith    Level: beginner
5721620fd73SBarry Smith 
5731620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
5741620fd73SBarry Smith M*/
5751620fd73SBarry 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);}
5761620fd73SBarry Smith 
5771620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
5781620fd73SBarry Smith 
5791620fd73SBarry 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);}
5801620fd73SBarry Smith 
5811620fd73SBarry Smith /*MC
5820d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
583bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
584bd481603SBarry Smith 
585bd481603SBarry Smith    Synopsis:
586f2ba6396SBarry Smith    #include "petscmat.h"
587c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
588bd481603SBarry Smith 
589bd481603SBarry Smith    Collective on MPI_Comm
590bd481603SBarry Smith 
591bd481603SBarry Smith    Input Parameters:
592bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
593859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
594859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
595bd481603SBarry Smith 
596bd481603SBarry Smith    Output Parameters:
597bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
598bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
599bd481603SBarry Smith 
600bd481603SBarry Smith 
601bd481603SBarry Smith    Level: intermediate
602bd481603SBarry Smith 
603bd481603SBarry Smith    Notes:
6040598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
605bd481603SBarry Smith 
6061d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
607bd481603SBarry Smith 
608bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
609bd481603SBarry Smith 
6101620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6111620fd73SBarry Smith 
612bd481603SBarry Smith   Concepts: preallocation^Matrix
613bd481603SBarry Smith 
614bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
615bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
616bd481603SBarry Smith M*/
617c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6187c922b88SBarry Smith { \
619a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
620a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
621a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
622eabe889fSLisandro Dalcin   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr); __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:
631f2ba6396SBarry 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 
647bd481603SBarry Smith    Level: intermediate
648bd481603SBarry Smith 
649bd481603SBarry Smith    Notes:
6500598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
651bd481603SBarry Smith 
6521d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
653bd481603SBarry Smith 
654bd481603SBarry Smith   Concepts: preallocation^Matrix
655bd481603SBarry Smith 
656bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
657bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
658bd481603SBarry Smith M*/
659784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
660c4f061fbSSatish Balay {\
661c1ac3661SBarry Smith   PetscInt __l;\
662784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
663784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
664c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
665ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
666c4f061fbSSatish Balay   }\
667c4f061fbSSatish Balay }
668c4f061fbSSatish Balay 
669bd481603SBarry Smith /*MC
670bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
671bd481603SBarry Smith        inserted using a local number of the rows and columns
672bd481603SBarry Smith 
673bd481603SBarry Smith    Synopsis:
674f2ba6396SBarry Smith    #include "petscmat.h"
675c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
676bd481603SBarry Smith 
677bd481603SBarry Smith    Not Collective
678bd481603SBarry Smith 
679bd481603SBarry Smith    Input Parameters:
680bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
681bd481603SBarry Smith .  nrows - the number of rows indicated
6821d73ed98SBarry Smith .  rows - the indices of the rows
683bd481603SBarry Smith .  ncols - the number of columns in the matrix
684bd481603SBarry Smith .  cols - the columns indicated
685bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
686bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
687bd481603SBarry Smith 
688bd481603SBarry Smith 
689bd481603SBarry Smith    Level: intermediate
690bd481603SBarry Smith 
691bd481603SBarry Smith    Notes:
6920598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
693bd481603SBarry Smith 
694bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
695bd481603SBarry Smith 
696bd481603SBarry Smith   Concepts: preallocation^Matrix
697bd481603SBarry Smith 
698bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
699bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
700bd481603SBarry Smith M*/
701d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
702d3d32019SBarry Smith {\
703c1ac3661SBarry Smith   PetscInt __l;\
704d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
705d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
706d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
707d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
708d3d32019SBarry Smith   }\
709d3d32019SBarry Smith }
710d3d32019SBarry Smith 
711bd481603SBarry Smith /*MC
712bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
713bd481603SBarry Smith        inserted using a local number of the rows and columns
714bd481603SBarry Smith 
715bd481603SBarry Smith    Synopsis:
716f2ba6396SBarry Smith    #include "petscmat.h"
717c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
718bd481603SBarry Smith 
719bd481603SBarry Smith    Not Collective
720bd481603SBarry Smith 
721bd481603SBarry Smith    Input Parameters:
72264075487SBarry Smith +  row - the row
723bd481603SBarry Smith .  ncols - the number of columns in the matrix
724a50f8bf6SHong Zhang -  cols - the columns indicated
725a50f8bf6SHong Zhang 
726a50f8bf6SHong Zhang    Output Parameters:
727a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
728bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
729bd481603SBarry Smith 
730bd481603SBarry Smith 
731bd481603SBarry Smith    Level: intermediate
732bd481603SBarry Smith 
733bd481603SBarry Smith    Notes:
7340598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
735bd481603SBarry Smith 
736bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
737bd481603SBarry Smith 
7381620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7391620fd73SBarry Smith 
740bd481603SBarry Smith   Concepts: preallocation^Matrix
741bd481603SBarry Smith 
742bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
743bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
744bd481603SBarry Smith M*/
745c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
746c1ac3661SBarry Smith { PetscInt __i; \
747e32f2f54SBarry 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);\
748e32f2f54SBarry 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);\
7497c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
75064075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
7517cc688d7SBarry Smith     else dnz[row - __rstart]++;\
7527c922b88SBarry Smith   }\
7537c922b88SBarry Smith }
7547c922b88SBarry Smith 
755bd481603SBarry Smith /*MC
756bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
757bd481603SBarry Smith        inserted using a local number of the rows and columns
758bd481603SBarry Smith 
759bd481603SBarry Smith    Synopsis:
760f2ba6396SBarry Smith    #include "petscmat.h"
761c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
762bd481603SBarry Smith 
763bd481603SBarry Smith    Not Collective
764bd481603SBarry Smith 
765bd481603SBarry Smith    Input Parameters:
766bd481603SBarry Smith +  nrows - the number of rows indicated
7671d73ed98SBarry Smith .  rows - the indices of the rows
768bd481603SBarry Smith .  ncols - the number of columns in the matrix
769bd481603SBarry Smith .  cols - the columns indicated
770bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
771bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
772bd481603SBarry Smith 
773bd481603SBarry Smith 
774bd481603SBarry Smith    Level: intermediate
775bd481603SBarry Smith 
776bd481603SBarry Smith    Notes:
7770598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
778bd481603SBarry Smith 
779bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
780bd481603SBarry Smith 
7811620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7821620fd73SBarry Smith 
783bd481603SBarry Smith   Concepts: preallocation^Matrix
784bd481603SBarry Smith 
785bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
786bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
787bd481603SBarry Smith M*/
788d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
789c1ac3661SBarry Smith { PetscInt __i; \
790d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
791d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
792d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
793d3d32019SBarry Smith   }\
794d3d32019SBarry Smith }
795d3d32019SBarry Smith 
796bd481603SBarry Smith /*MC
79716371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
79816371a99SBarry Smith 
79916371a99SBarry Smith    Synopsis:
800f2ba6396SBarry Smith    #include "petscmat.h"
80116371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
80216371a99SBarry Smith 
80316371a99SBarry Smith    Not Collective
80416371a99SBarry Smith 
80516371a99SBarry Smith    Input Parameters:
80616371a99SBarry Smith .  A - matrix
80716371a99SBarry Smith .  row - row where values exist (must be local to this process)
80816371a99SBarry Smith .  ncols - number of columns
80916371a99SBarry Smith .  cols - columns with nonzeros
81016371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
81116371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
81216371a99SBarry Smith 
81316371a99SBarry Smith 
81416371a99SBarry Smith    Level: intermediate
81516371a99SBarry Smith 
81616371a99SBarry Smith    Notes:
8170598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
81816371a99SBarry Smith 
81916371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
82016371a99SBarry Smith 
82116371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
82216371a99SBarry Smith 
82316371a99SBarry Smith   Concepts: preallocation^Matrix
82416371a99SBarry Smith 
82516371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
826eabe889fSLisandro Dalcin           MatPreallocateSymmetricSetLocal()
82716371a99SBarry Smith M*/
8280298fd71SBarry 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);}
82916371a99SBarry Smith 
83016371a99SBarry Smith 
83116371a99SBarry Smith /*MC
8320d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
833bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
834bd481603SBarry Smith 
835bd481603SBarry Smith    Synopsis:
836f2ba6396SBarry Smith    #include "petscmat.h"
837c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
838bd481603SBarry Smith 
839bd481603SBarry Smith    Collective on MPI_Comm
840bd481603SBarry Smith 
841bd481603SBarry Smith    Input Parameters:
84216371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
843bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
844bd481603SBarry Smith 
845bd481603SBarry Smith 
846bd481603SBarry Smith    Level: intermediate
847bd481603SBarry Smith 
848bd481603SBarry Smith    Notes:
8490598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
850bd481603SBarry Smith 
851bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
852bd481603SBarry Smith 
8531620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
8541620fd73SBarry Smith 
855bd481603SBarry Smith   Concepts: preallocation^Matrix
856bd481603SBarry Smith 
857bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
858eabe889fSLisandro Dalcin           MatPreallocateSymmetricSetLocal()
859bd481603SBarry Smith M*/
860a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
8617c922b88SBarry Smith 
862bd481603SBarry Smith 
863bd481603SBarry Smith 
8647b80b807SBarry Smith /* Routines unique to particular data structures */
865014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
866698d4c6aSKris Buschelman 
867014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
868014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
869698d4c6aSKris Buschelman 
870014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
871014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
872014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
873014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
874014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
875014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
8767b80b807SBarry Smith 
877a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
878a96a251dSBarry Smith 
879014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
880014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
881014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
882273d9f13SBarry Smith 
883014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
884014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
885014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
886014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
887014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
888014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
889014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
890014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
891014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
892014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
8939230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
8949230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
895014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
896273d9f13SBarry Smith 
897014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt);
898014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
8991b807ce4Svictorle 
900014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
901014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
9022e8a6d31SBarry Smith 
903014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
9043a7fca6bSBarry Smith 
905014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
9067b80b807SBarry Smith /*
9077b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
90894b7f48cSBarry Smith   done through the KSP and PC interfaces.
9097b80b807SBarry Smith */
9107b80b807SBarry Smith 
91176bdecfbSBarry Smith /*J
9128f6c3df8SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering
913d9274352SBarry Smith 
914d9274352SBarry Smith    Level: beginner
915d9274352SBarry Smith 
916d9274352SBarry Smith .seealso: MatGetOrdering()
91776bdecfbSBarry Smith J*/
91819fd82e9SBarry Smith typedef const char* MatOrderingType;
9192692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
9202692d6eeSBarry Smith #define MATORDERINGND          "nd"
9212692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
9222692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
9232692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
9242692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
925312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
926b12f92e5SBarry Smith 
92719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
928140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
929bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
930607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void);
931014dd563SJed Brown PETSC_EXTERN PetscBool         MatOrderingRegisterAllCalled;
932140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList;
933d4fbbf0eSBarry Smith 
934014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
935a2ce50c7SBarry Smith 
936d91e6319SBarry Smith /*S
937d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
938d90ac83dSHong Zhang 
939d90ac83dSHong Zhang    Level: beginner
940d90ac83dSHong Zhang 
941d90ac83dSHong Zhang S*/
942d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
9436a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[];
9445e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];
945d90ac83dSHong Zhang 
946d90ac83dSHong Zhang /*S
9472401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
948b00f7748SHong Zhang 
94961cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
95061cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
951b00f7748SHong Zhang 
95215e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
953b00f7748SHong Zhang 
95461cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
95561cad860SBarry Smith 
956b00f7748SHong Zhang    Level: developer
957b00f7748SHong Zhang 
958d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
959d7d82daaSBarry Smith           MatFactorInfoInitialize()
960b00f7748SHong Zhang 
961b00f7748SHong Zhang S*/
962b00f7748SHong Zhang typedef struct {
96315e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
96485317021SBarry Smith   PetscReal     usedt;
96515e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
966b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
96715e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
96867e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
969348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
970bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
971bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
97215e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
973f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
974f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
97515e8a5b3SHong Zhang } MatFactorInfo;
976ffa6d0a5SLois Curfman McInnes 
977014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
978014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
979014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
980014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
981014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
982014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
983014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
984014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
985014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
986014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
987014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
988014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
989014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
990014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
991014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
992014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
993014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
994014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
995014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
9968ed539a5SBarry Smith 
997014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
998bb5a7306SBarry Smith 
999d91e6319SBarry Smith /*E
1000d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1001bb1eb677SSatish Balay 
1002d91e6319SBarry Smith     Level: beginner
1003d91e6319SBarry Smith 
1004d9274352SBarry Smith    May be bitwise ORd together
1005d9274352SBarry Smith 
1006d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1007d91e6319SBarry Smith 
10084e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10094e7234bfSBarry Smith 
101041f059aeSBarry Smith .seealso: MatSOR()
1011d91e6319SBarry Smith E*/
1012ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1013ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1014ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
101584cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1016014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10178ed539a5SBarry Smith 
1018d4fbbf0eSBarry Smith /*
1019639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1020639f9d9dSBarry Smith */
1021b12f92e5SBarry Smith 
102276bdecfbSBarry Smith /*J
10238f6c3df8SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring
1024d9274352SBarry Smith 
1025d9274352SBarry Smith    Level: beginner
1026d9274352SBarry Smith 
10278f6c3df8SBarry Smith .seealso: MatGetColoring(), MatColoring
102876bdecfbSBarry Smith J*/
102919fd82e9SBarry Smith typedef const char* MatColoringType;
10302692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
10312692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
10322692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
10332692d6eeSBarry Smith #define MATCOLORINGID      "id"
1034b12f92e5SBarry Smith 
103519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetColoring(Mat,MatColoringType,ISColoring*);
1036bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
103730de9b25SBarry Smith 
1038014dd563SJed Brown PETSC_EXTERN PetscBool MatColoringRegisterAllCalled;
1039f1144a30SSatish Balay 
1040607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void);
1041014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1042639f9d9dSBarry Smith 
1043d9274352SBarry Smith /*S
1044d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1045d9274352SBarry Smith         and coloring
1046639f9d9dSBarry Smith 
1047d9274352SBarry Smith    Level: beginner
1048d9274352SBarry Smith 
1049d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1050d9274352SBarry Smith 
1051d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1052d9274352SBarry Smith S*/
1053e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1054639f9d9dSBarry Smith 
1055014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1056014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1057014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1058014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1059014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1060014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1061014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1062014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1063014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1064014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1065b1683b59SHong Zhang 
1066b1683b59SHong Zhang /*S
1067b9af6bddSHong Zhang      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1068b1683b59SHong Zhang 
1069b1683b59SHong Zhang    Level: beginner
1070b1683b59SHong Zhang 
1071b1683b59SHong Zhang   Concepts: coloring, sparse matrix product
1072b1683b59SHong Zhang 
1073b9af6bddSHong Zhang .seealso:  MatTransposeColoringCreate()
1074b1683b59SHong Zhang S*/
1075b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1076b1683b59SHong Zhang 
1077014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1078014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1079014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1080014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1081b1683b59SHong Zhang 
1082639f9d9dSBarry Smith /*
10830752156aSBarry Smith     These routines are for partitioning matrices: currently used only
10843eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
10850752156aSBarry Smith */
1086ca161407SBarry Smith 
1087d9274352SBarry Smith /*S
1088d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1089d9274352SBarry Smith 
1090d9274352SBarry Smith    Level: beginner
1091d9274352SBarry Smith 
1092d9274352SBarry Smith   Concepts: partitioning
1093d9274352SBarry Smith 
1094743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1095d9274352SBarry Smith S*/
109691e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1097d9274352SBarry Smith 
109876bdecfbSBarry Smith /*J
10998f6c3df8SBarry Smith     MatPartitioningType - String with the name of a PETSc matrix partitioning
1100d9274352SBarry Smith 
1101d9274352SBarry Smith    Level: beginner
11020d04baf8SBarry Smith dm
1103b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
110476bdecfbSBarry Smith J*/
110519fd82e9SBarry Smith typedef const char* MatPartitioningType;
11062692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
11072692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
11082692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
11092692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
11102692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
111150d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
111271306c60Sjroman 
1113ca161407SBarry Smith 
1114014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
111519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1116014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1117014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1118014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1119014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1120014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1121014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
11222aabb6bbSBarry Smith 
1123bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
112430de9b25SBarry Smith 
1125014dd563SJed Brown PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled;
1126f1144a30SSatish Balay 
1127607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void);
11282bad1931SBarry Smith 
1129014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1130014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
113119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1132ca161407SBarry Smith 
1133014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1134014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
11350752156aSBarry Smith 
1136b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
11376a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1138b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
11396a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[];
1140b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
11416a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[];
1142b6956c12SJose Roman 
1143014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1144014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1145014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1146014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1147014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1148014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1149014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1150014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1151014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1152014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1153014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
115471306c60Sjroman 
115571306c60Sjroman #define MP_PARTY_OPT "opt"
115671306c60Sjroman #define MP_PARTY_LIN "lin"
115771306c60Sjroman #define MP_PARTY_SCA "sca"
115871306c60Sjroman #define MP_PARTY_RAN "ran"
115971306c60Sjroman #define MP_PARTY_GBF "gbf"
116071306c60Sjroman #define MP_PARTY_GCF "gcf"
116171306c60Sjroman #define MP_PARTY_BUB "bub"
116271306c60Sjroman #define MP_PARTY_DEF "def"
1163014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
116471306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
116571306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
116671306c60Sjroman #define MP_PARTY_NONE "no"
1167014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1168014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1169014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1170014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
117171306c60Sjroman 
117250d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
11736a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1174e0f1cffaSJose Roman 
1175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1178014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
117971306c60Sjroman 
1180b43b03e9SMark F. Adams /*
1181b43b03e9SMark F. Adams     These routines are for coarsening matrices:
1182b43b03e9SMark F. Adams */
1183b43b03e9SMark F. Adams 
1184b43b03e9SMark F. Adams /*S
1185b43b03e9SMark F. Adams      MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix)
1186b43b03e9SMark F. Adams 
1187b43b03e9SMark F. Adams    Level: beginner
1188b43b03e9SMark F. Adams 
1189b43b03e9SMark F. Adams   Concepts: coarsen
1190b43b03e9SMark F. Adams 
1191b43b03e9SMark F. Adams .seealso:  MatCoarsenCreate), MatCoarsenType
1192b43b03e9SMark F. Adams S*/
1193b43b03e9SMark F. Adams typedef struct _p_MatCoarsen* MatCoarsen;
1194b43b03e9SMark F. Adams 
1195b43b03e9SMark F. Adams /*J
11968f6c3df8SBarry Smith     MatCoarsenType - String with the name of a PETSc matrix coarsen
1197b43b03e9SMark F. Adams 
1198b43b03e9SMark F. Adams    Level: beginner
11998f6c3df8SBarry Smith 
1200b43b03e9SMark F. Adams .seealso: MatCoarsenCreate(), MatCoarsen
1201b43b03e9SMark F. Adams J*/
120219fd82e9SBarry Smith typedef const char* MatCoarsenType;
1203b43b03e9SMark F. Adams #define MATCOARSENMIS  "mis"
12049057884aSMark F. Adams #define MATCOARSENHEM  "hem"
1205b43b03e9SMark F. Adams 
12060cbbd2e1SMark F. Adams /* linked list for aggregates */
120741b27cdeSMark F. Adams typedef struct _PetscCDIntNd{
120841b27cdeSMark F. Adams   struct _PetscCDIntNd *next;
12090cbbd2e1SMark F. Adams   PetscInt             gid;
121041b27cdeSMark F. Adams }PetscCDIntNd;
12110cbbd2e1SMark F. Adams 
12120cbbd2e1SMark F. Adams /* only used by node pool */
121341b27cdeSMark F. Adams typedef struct _PetscCDArrNd{
121441b27cdeSMark F. Adams   struct _PetscCDArrNd *next;
121541b27cdeSMark F. Adams   struct _PetscCDIntNd *array;
121641b27cdeSMark F. Adams }PetscCDArrNd;
12170cbbd2e1SMark F. Adams 
12180cbbd2e1SMark F. Adams typedef struct _PetscCoarsenData{
121982c86c8fSBarry Smith   PetscCDArrNd pool_list;  /* node pool */
122041b27cdeSMark F. Adams   PetscCDIntNd *new_node;
12210cbbd2e1SMark F. Adams   PetscInt     new_left;
12220cbbd2e1SMark F. Adams   PetscInt     chk_sz;
122341b27cdeSMark F. Adams   PetscCDIntNd *extra_nodes;
122482c86c8fSBarry Smith   PetscCDIntNd **array;  /* Array of lists */
12250cbbd2e1SMark F. Adams   PetscInt     size;
122682c86c8fSBarry Smith   Mat          mat;  /* cache a Mat for communication data */
12270cbbd2e1SMark F. Adams }PetscCoarsenData;
12280cbbd2e1SMark F. Adams 
1229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*);
123019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,MatCoarsenType);
1231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat);
1232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS);
1233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool);
1234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetVerbose(MatCoarsen,PetscInt);
1235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** );
1236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen);
1237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*);
1238b43b03e9SMark F. Adams 
1239bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],PetscErrorCode (*)(MatCoarsen));
1240b43b03e9SMark F. Adams 
1241014dd563SJed Brown PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled;
1242b43b03e9SMark F. Adams 
1243607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void);
1244b43b03e9SMark F. Adams 
1245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer);
1246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen);
124719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,MatCoarsenType*);
1248b43b03e9SMark F. Adams 
1249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1251591294e4SBarry Smith 
12520752156aSBarry Smith /*
12530a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1254d4fbbf0eSBarry Smith */
12551c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
12561c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
12571c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
12581c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
12591c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
12607c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
12617c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
12621c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
12631c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
12647c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
12657c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
12661c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
12671c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1268a714c06dSBarry Smith                MATOP_SOR=13,
12691c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
12701c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
12711c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
12721c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
12731c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
12741c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
12751c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
12761c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1277d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1278d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1279d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1280d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1281d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1282d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1283d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1284d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1285d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1286d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1287d519adbfSMatthew Knepley                MATOP_GET_ARRAY=32,
1288d519adbfSMatthew Knepley                MATOP_RESTORE_ARRAY=33,
1289d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1290d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1291d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1292d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1293d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1294d519adbfSMatthew Knepley                MATOP_AXPY=39,
1295d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1296d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1297d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1298d519adbfSMatthew Knepley                MATOP_COPY=43,
1299d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1300d519adbfSMatthew Knepley                MATOP_SCALE=45,
1301d519adbfSMatthew Knepley                MATOP_SHIFT=46,
130235153367SBarry Smith                MATOP_DIAGONAL_SET=47,
1303d519adbfSMatthew Knepley                MATOP_ILUDT_FACTOR=48,
1304d519adbfSMatthew Knepley                MATOP_SET_BLOCK_SIZE=49,
1305d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1306d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1307d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1308d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1309d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1310d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1311d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1312d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1313d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1314d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1315d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1316d519adbfSMatthew Knepley                MATOP_VIEW=61,
1317d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
13187bab7c10SHong Zhang                MATOP_MATMAT_MULT=63,
13197bab7c10SHong Zhang                MATOP_MATMAT_MULT_SYMBOLIC=64,
13207bab7c10SHong Zhang                MATOP_MATMAT_MULT_NUMERIC=65,
1321d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1322d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1323d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1324d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1325d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1326d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1327d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
1328bcaeba4dSBarry Smith                MATOP_PLACEHOLDER=73,
1329d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1330d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1331d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1332bda6c4cbSBarry Smith                MATOP_MULT_CONSTRAINED=77,
1333bda6c4cbSBarry Smith                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
1334d519adbfSMatthew Knepley                MATOP_PERMUTE_SPARSIFY=79,
1335d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1336d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1337d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1338d519adbfSMatthew Knepley                MATOP_LOAD=83,
1339d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1340d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1341d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1342820469bcSHong Zhang                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1343d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1344d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1345d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1346d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1347d519adbfSMatthew Knepley                MATOP_PTAP=92,
1348d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1349d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1350bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT=95,
1351bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1352bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
1353820469bcSHong Zhang                MATOP_DUMMY98=98,
1354820469bcSHong Zhang                MATOP_DUMMY99=99,
1355820469bcSHong Zhang                MATOP_DUMMY100=100,
1356820469bcSHong Zhang                MATOP_DUMMY101=101,
1357d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
1358d519adbfSMatthew Knepley                MATOP_SET_SIZES=103,
1359d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1360d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1361bda6c4cbSBarry Smith                MATOP_IMAGINARY_PART=106,
1362bda6c4cbSBarry Smith                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1363bda6c4cbSBarry Smith                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1364bda6c4cbSBarry Smith                MATOP_MAT_SOLVE=109,
13650e8d04f7SBarry Smith                MATOP_GET_REDUNDANT_MATRIX=110,
1366d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
13670e8d04f7SBarry Smith                MATOP_GET_COLUMN_VECTOR=112,
1368d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
13690e8d04f7SBarry Smith                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
137089c6957cSBarry Smith                MATOP_CREATE=115,
137189c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
13720e8d04f7SBarry Smith                MATOP_GET_LOCAL_SUB_MATRIX=117,
13730e8d04f7SBarry Smith                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1374eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
13750e8d04f7SBarry Smith                MATOP_HERMITIAN_TRANSPOSE=120,
13760e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
13770e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
13780e8d04f7SBarry Smith                MATOP_GET_MULTI_PROC_BLOCK=123,
13790e8d04f7SBarry Smith                MATOP_GET_COLUMN_NORMS=125,
13800e8d04f7SBarry Smith                MATOP_GET_SUB_MATRICES_PARALLE=128,
13812b8ad9a3SHong Zhang                MATOP_SET_VALUES_BATCH=129,
13820e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT=130,
13830e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
13840e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
13850e8d04f7SBarry Smith                MATOP_TRANSPOSE_COLORING_CREAT=133,
13860e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_SPT=134,
13870e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_DEN=135,
13880e8d04f7SBarry Smith                MATOP_RART=136,
13890e8d04f7SBarry Smith                MATOP_RART_SYMBOLIC=137,
13900e8d04f7SBarry Smith                MATOP_RART_NUMERIC=138,
1391e09a3074SHong Zhang                MATOP_SET_BLOCK_SIZES=139,
1392*f9426fe0SMark Adams                MATOP_AYPX=140,
1393*f9426fe0SMark Adams                MATOP_RESIDUAL=141
1394fae171e0SBarry Smith              } MatOperation;
1395014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *);
1396014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1397014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1398014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1399112a2221SBarry Smith 
140090ace30eSBarry Smith /*
140190ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
140290ace30eSBarry Smith    stored in a universal format. By changing the format with
14037973ad04SBarry Smith    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
140490ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
140590ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1406f4403165SShri Abhyankar    read into matrices of the same type.
140790ace30eSBarry Smith */
140890ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
140990ace30eSBarry Smith 
1410014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1411014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
14131f4e1ec7SBarry Smith 
1414d9274352SBarry Smith /*S
1415d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1416d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1417d9274352SBarry Smith 
1418f7a9e4ceSBarry Smith    Level: advanced
1419d9274352SBarry Smith 
1420d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1421d9274352SBarry Smith 
14226e1639daSBarry Smith   Users manual sections:
14236e1639daSBarry Smith .   sec_singular
14246e1639daSBarry Smith 
1425d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1426d9274352SBarry Smith S*/
142774637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1428d9274352SBarry Smith 
1429014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1430014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1431014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1432d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1433014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1434014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1435014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1436014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1437014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1438014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1439014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1440014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
144174637425SBarry Smith 
1442014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1443014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1444014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1445014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
14463f1d51d7SBarry Smith 
1447014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1448014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1449014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1450c4f061fbSSatish Balay 
1451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*);
1452b0a32e0cSBarry Smith 
1453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
145404f1ad80SBarry Smith 
1455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace);
1461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1469e884886fSBarry Smith 
14706370053bSBarry Smith /*S
14716370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
14726370053bSBarry Smith               Jacobian vector products
1473e884886fSBarry Smith 
14746370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
14756370053bSBarry Smith 
14766370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
14776370053bSBarry Smith 
14786370053bSBarry Smith     Level: developer
14796370053bSBarry Smith 
14806370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
14816370053bSBarry Smith S*/
1482e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1483e884886fSBarry Smith 
148476bdecfbSBarry Smith /*J
1485e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1486e884886fSBarry Smith 
1487e884886fSBarry Smith    Level: beginner
1488e884886fSBarry Smith 
1489e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
149076bdecfbSBarry Smith J*/
149119fd82e9SBarry Smith typedef const char* MatMFFDType;
1492e884886fSBarry Smith #define MATMFFD_DS  "ds"
1493e884886fSBarry Smith #define MATMFFD_WP  "wp"
1494e884886fSBarry Smith 
149519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1496bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));
1497e884886fSBarry Smith 
1498607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegisterAll(void);
1499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool );
1501e884886fSBarry Smith 
1502014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1503014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
15047dbadf16SMatthew Knepley 
150597969023SHong Zhang /*
150697969023SHong Zhang    PETSc interface to MUMPS
150797969023SHong Zhang */
150897969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1510b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
151197969023SHong Zhang #endif
151223a5497aSJed Brown 
1513d954db57SHong Zhang /*
1514d954db57SHong Zhang    PETSc interface to SUPERLU
1515d954db57SHong Zhang */
1516d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1518d954db57SHong Zhang #endif
1519d954db57SHong Zhang 
1520aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA
15213f9c0db1SPaul Mullowney /*E
1522e057df02SPaul Mullowney     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
15232692e278SPaul Mullowney     matrices.
1524e057df02SPaul Mullowney 
1525e057df02SPaul Mullowney     Not Collective
1526e057df02SPaul Mullowney 
15278468deeeSKarl Rupp +   MAT_CUSPARSE_CSR - Compressed Sparse Row
15282692e278SPaul Mullowney .   MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later).
15292692e278SPaul Mullowney -   MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).
1530e057df02SPaul Mullowney 
1531e057df02SPaul Mullowney     Level: intermediate
1532e057df02SPaul Mullowney 
1533e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1534e057df02SPaul Mullowney 
1535e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1536e057df02SPaul Mullowney E*/
1537e057df02SPaul Mullowney 
15383f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;
1539e057df02SPaul Mullowney 
1540e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1541e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
1542e057df02SPaul Mullowney 
15433f9c0db1SPaul Mullowney /*E
1544e057df02SPaul Mullowney     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
15452692e278SPaul Mullowney     matrices whose operation should use a particular storage format.
1546e057df02SPaul Mullowney 
1547e057df02SPaul Mullowney     Not Collective
1548e057df02SPaul Mullowney 
15498468deeeSKarl Rupp +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
15508468deeeSKarl Rupp .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
15518468deeeSKarl Rupp .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
15528468deeeSKarl Rupp -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices
1553e057df02SPaul Mullowney 
1554e057df02SPaul Mullowney     Level: intermediate
1555e057df02SPaul Mullowney 
1556e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1557e057df02SPaul Mullowney E*/
155836d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;
1559e057df02SPaul Mullowney 
156010e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
156110e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1562e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
15639ae82921SPaul Mullowney #endif
15649ae82921SPaul Mullowney 
156590c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
1566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1567014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1568e057df02SPaul Mullowney 
15693f9c0db1SPaul Mullowney /*E
1570e057df02SPaul Mullowney     MatCUSPStorageFormat - indicates the storage format for CUSP (GPU)
157136d62e41SPaul Mullowney     matrices.
1572e057df02SPaul Mullowney 
1573e057df02SPaul Mullowney     Not Collective
1574e057df02SPaul Mullowney 
15758468deeeSKarl Rupp +   MAT_CUSP_CSR - Compressed Sparse Row
15768468deeeSKarl Rupp .   MAT_CUSP_DIA - Diagonal
15778468deeeSKarl Rupp -   MAT_CUSP_ELL - Ellpack
1578e057df02SPaul Mullowney 
1579e057df02SPaul Mullowney     Level: intermediate
1580e057df02SPaul Mullowney 
1581e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1582e057df02SPaul Mullowney 
1583e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
1584e057df02SPaul Mullowney E*/
15853f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat;
1586e057df02SPaul Mullowney 
1587e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1588e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[];
1589e057df02SPaul Mullowney 
15903f9c0db1SPaul Mullowney /*E
1591e057df02SPaul Mullowney     MatCUSPFormatOperation - indicates the operation of CUSP (GPU)
159236d62e41SPaul Mullowney     matrices whose operation should use a particular storage format.
1593e057df02SPaul Mullowney 
1594e057df02SPaul Mullowney     Not Collective
1595e057df02SPaul Mullowney 
15968468deeeSKarl Rupp +   MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
15978468deeeSKarl Rupp .   MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
15988468deeeSKarl Rupp .   MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
15998468deeeSKarl Rupp -   MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices
1600e057df02SPaul Mullowney 
1601e057df02SPaul Mullowney     Level: intermediate
1602e057df02SPaul Mullowney 
1603e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1604e057df02SPaul Mullowney 
1605e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat
1606e057df02SPaul Mullowney E*/
160736d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation;
1608e057df02SPaul Mullowney 
1609e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat);
1610e057df02SPaul Mullowney #endif
161190c53902SBarry Smith 
1612d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL)
1613d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1614d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1615d67ff14aSKarl Rupp #endif
1616d67ff14aSKarl Rupp 
161754efbe56SHong Zhang /*
161854efbe56SHong Zhang    PETSc interface to FFTW
161954efbe56SHong Zhang */
162054efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
1621014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1622014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
1623014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*);
162454efbe56SHong Zhang #endif
162554efbe56SHong Zhang 
1626db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL)
1627607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
1628db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
1629db31f6deSJed Brown #endif
1630db31f6deSJed Brown 
1631014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1632014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1633014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1634014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1635014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1636014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
163719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1638014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1639014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1640d8588912SDave May 
16414325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
16424325cce7SMatthew G Knepley 
164323a5497aSJed Brown #endif
1644