xref: /petsc/include/petscmat.h (revision b2bf63709b6a963fd02e2640983e9e8c2a88c2e3)
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);
378f5edf698SHong Zhang 
379d91e6319SBarry Smith /*E
380d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
381d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
382d91e6319SBarry Smith 
383d91e6319SBarry Smith     Level: beginner
384d91e6319SBarry Smith 
385d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
386d91e6319SBarry Smith 
38770dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
38870dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
38970dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
39070dcbbb9SBarry Smith 
391d91e6319SBarry Smith .seealso: MatDuplicate()
392d91e6319SBarry Smith E*/
39370dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
3942e8a6d31SBarry Smith 
39519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
396014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);
39794a9d846SBarry Smith 
398cb5b572fSBarry Smith 
399014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
400014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
401014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
402014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
403014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
404014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
405014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
406014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
407014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
4087b80b807SBarry Smith 
4091a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4101a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
4111a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4121a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
413d4fbbf0eSBarry Smith 
414d91e6319SBarry Smith /*S
415d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
416d91e6319SBarry Smith 
417d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
418d91e6319SBarry Smith 
419d91e6319SBarry Smith    Level: intermediate
420d91e6319SBarry Smith 
421d91e6319SBarry Smith   Concepts: matrix^nonzero information
422d91e6319SBarry Smith 
423d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
424d91e6319SBarry Smith S*/
4254e220ebcSLois Curfman McInnes typedef struct {
426b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
427b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
428b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
429b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
430b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
431b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
432b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4334e220ebcSLois Curfman McInnes } MatInfo;
4344e220ebcSLois Curfman McInnes 
435d9274352SBarry Smith /*E
436d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
437d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
438d9274352SBarry Smith 
439d9274352SBarry Smith     Level: beginner
440d9274352SBarry Smith 
441d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
442d9274352SBarry Smith 
443d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
444d9274352SBarry Smith E*/
4457b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
446014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
447014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
448014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
449014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
450014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
452014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat *);
456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);
458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool *);
459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
4637b80b807SBarry Smith 
464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal *);
465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *);
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
4737b80b807SBarry Smith 
474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
480566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);
4815ef9f2a5SBarry Smith 
482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);
4907b80b807SBarry Smith 
491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,Mat*);
493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJNumeric(MPI_Comm,Mat,PetscInt,Mat);
494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5018efafbd8SBarry Smith 
502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5037b80b807SBarry Smith 
504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat);
50722440eb1SKris Buschelman 
5087bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
5097bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultSymbolic(Mat,Mat,Mat,PetscReal,Mat*);
5107bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultNumeric(Mat,Mat,Mat,Mat);
5117bab7c10SHong Zhang 
512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat);
515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat);
51822440eb1SKris Buschelman 
519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat);
522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat);
525bc011b1eSHong Zhang 
526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);
5281c741599SHong Zhang 
529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);
5317b80b807SBarry Smith 
532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
542052efed2SBarry Smith 
543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
54590f02eecSBarry Smith 
546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecs(Mat,Vec*,Vec*);
550*b2bf6370SHong Zhang PETSC_EXTERN PetscErrorCode MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
553bd481603SBarry Smith 
554bd481603SBarry Smith /*MC
5551620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5561620fd73SBarry Smith 
5571620fd73SBarry Smith    Not collective
5581620fd73SBarry Smith 
5591620fd73SBarry Smith    Input Parameters:
5601620fd73SBarry Smith +  m - the matrix
5611620fd73SBarry Smith .  row - the row location of the entry
5621620fd73SBarry Smith .  col - the column location of the entry
5631620fd73SBarry Smith .  value - the value to insert
5641620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5651620fd73SBarry Smith 
5661620fd73SBarry Smith    Notes:
5671620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
5681620fd73SBarry Smith    values simultaneously if possible.
5691620fd73SBarry Smith 
5701620fd73SBarry Smith    Level: beginner
5711620fd73SBarry Smith 
5721620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
5731620fd73SBarry Smith M*/
5741620fd73SBarry 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);}
5751620fd73SBarry Smith 
5761620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
5771620fd73SBarry Smith 
5781620fd73SBarry 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);}
5791620fd73SBarry Smith 
5801620fd73SBarry Smith /*MC
5810d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
582bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
583bd481603SBarry Smith 
584bd481603SBarry Smith    Synopsis:
585f2ba6396SBarry Smith    #include "petscmat.h"
586c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
587bd481603SBarry Smith 
588bd481603SBarry Smith    Collective on MPI_Comm
589bd481603SBarry Smith 
590bd481603SBarry Smith    Input Parameters:
591bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
592859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
593859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
594bd481603SBarry Smith 
595bd481603SBarry Smith    Output Parameters:
596bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
597bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
598bd481603SBarry Smith 
599bd481603SBarry Smith 
600bd481603SBarry Smith    Level: intermediate
601bd481603SBarry Smith 
602bd481603SBarry Smith    Notes:
6030598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
604bd481603SBarry Smith 
6051d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
606bd481603SBarry Smith 
607bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
608bd481603SBarry Smith 
6091620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6101620fd73SBarry Smith 
611bd481603SBarry Smith   Concepts: preallocation^Matrix
612bd481603SBarry Smith 
613bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
614bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
615bd481603SBarry Smith M*/
616c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6177c922b88SBarry Smith { \
618a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
619a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
620a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
621eabe889fSLisandro Dalcin   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr); __start = 0; __end = __start; \
622c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
623a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6247c922b88SBarry Smith 
625bd481603SBarry Smith /*MC
626bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
627bd481603SBarry Smith        inserted using a local number of the rows and columns
628bd481603SBarry Smith 
629bd481603SBarry Smith    Synopsis:
630f2ba6396SBarry Smith    #include "petscmat.h"
631c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
632bd481603SBarry Smith 
633bd481603SBarry Smith    Not Collective
634bd481603SBarry Smith 
635bd481603SBarry Smith    Input Parameters:
636784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
637bd481603SBarry Smith .  nrows - the number of rows indicated
6381d73ed98SBarry Smith .  rows - the indices of the rows
639784ac674SJed Brown .  cmap - the column mapping from local to global numbering
640bd481603SBarry Smith .  ncols - the number of columns in the matrix
641bd481603SBarry Smith .  cols - the columns indicated
642bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
643bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
644bd481603SBarry Smith 
645bd481603SBarry Smith 
646bd481603SBarry Smith    Level: intermediate
647bd481603SBarry Smith 
648bd481603SBarry Smith    Notes:
6490598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
650bd481603SBarry Smith 
6511d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
652bd481603SBarry Smith 
653bd481603SBarry Smith   Concepts: preallocation^Matrix
654bd481603SBarry Smith 
655bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
656bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
657bd481603SBarry Smith M*/
658784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
659c4f061fbSSatish Balay {\
660c1ac3661SBarry Smith   PetscInt __l;\
661784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
662784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
663c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
664ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
665c4f061fbSSatish Balay   }\
666c4f061fbSSatish Balay }
667c4f061fbSSatish Balay 
668bd481603SBarry Smith /*MC
669bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
670bd481603SBarry Smith        inserted using a local number of the rows and columns
671bd481603SBarry Smith 
672bd481603SBarry Smith    Synopsis:
673f2ba6396SBarry Smith    #include "petscmat.h"
674c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
675bd481603SBarry Smith 
676bd481603SBarry Smith    Not Collective
677bd481603SBarry Smith 
678bd481603SBarry Smith    Input Parameters:
679bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
680bd481603SBarry Smith .  nrows - the number of rows indicated
6811d73ed98SBarry Smith .  rows - the indices of the rows
682bd481603SBarry Smith .  ncols - the number of columns in the matrix
683bd481603SBarry Smith .  cols - the columns indicated
684bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
685bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
686bd481603SBarry Smith 
687bd481603SBarry Smith 
688bd481603SBarry Smith    Level: intermediate
689bd481603SBarry Smith 
690bd481603SBarry Smith    Notes:
6910598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
692bd481603SBarry Smith 
693bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
694bd481603SBarry Smith 
695bd481603SBarry Smith   Concepts: preallocation^Matrix
696bd481603SBarry Smith 
697bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
698bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
699bd481603SBarry Smith M*/
700d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
701d3d32019SBarry Smith {\
702c1ac3661SBarry Smith   PetscInt __l;\
703d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
704d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
705d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
706d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
707d3d32019SBarry Smith   }\
708d3d32019SBarry Smith }
709d3d32019SBarry Smith 
710bd481603SBarry Smith /*MC
711bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
712bd481603SBarry Smith        inserted using a local number of the rows and columns
713bd481603SBarry Smith 
714bd481603SBarry Smith    Synopsis:
715f2ba6396SBarry Smith    #include "petscmat.h"
716c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
717bd481603SBarry Smith 
718bd481603SBarry Smith    Not Collective
719bd481603SBarry Smith 
720bd481603SBarry Smith    Input Parameters:
72164075487SBarry Smith +  row - the row
722bd481603SBarry Smith .  ncols - the number of columns in the matrix
723a50f8bf6SHong Zhang -  cols - the columns indicated
724a50f8bf6SHong Zhang 
725a50f8bf6SHong Zhang    Output Parameters:
726a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
727bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
728bd481603SBarry Smith 
729bd481603SBarry Smith 
730bd481603SBarry Smith    Level: intermediate
731bd481603SBarry Smith 
732bd481603SBarry Smith    Notes:
7330598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
734bd481603SBarry Smith 
735bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
736bd481603SBarry Smith 
7371620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7381620fd73SBarry Smith 
739bd481603SBarry Smith   Concepts: preallocation^Matrix
740bd481603SBarry Smith 
741bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
742bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
743bd481603SBarry Smith M*/
744c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
745c1ac3661SBarry Smith { PetscInt __i; \
746e32f2f54SBarry 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);\
747e32f2f54SBarry 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);\
7487c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
74964075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
7507cc688d7SBarry Smith     else dnz[row - __rstart]++;\
7517c922b88SBarry Smith   }\
7527c922b88SBarry Smith }
7537c922b88SBarry Smith 
754bd481603SBarry Smith /*MC
755bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
756bd481603SBarry Smith        inserted using a local number of the rows and columns
757bd481603SBarry Smith 
758bd481603SBarry Smith    Synopsis:
759f2ba6396SBarry Smith    #include "petscmat.h"
760c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
761bd481603SBarry Smith 
762bd481603SBarry Smith    Not Collective
763bd481603SBarry Smith 
764bd481603SBarry Smith    Input Parameters:
765bd481603SBarry Smith +  nrows - the number of rows indicated
7661d73ed98SBarry Smith .  rows - the indices of the rows
767bd481603SBarry Smith .  ncols - the number of columns in the matrix
768bd481603SBarry Smith .  cols - the columns indicated
769bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
770bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
771bd481603SBarry Smith 
772bd481603SBarry Smith 
773bd481603SBarry Smith    Level: intermediate
774bd481603SBarry Smith 
775bd481603SBarry Smith    Notes:
7760598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
777bd481603SBarry Smith 
778bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
779bd481603SBarry Smith 
7801620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7811620fd73SBarry Smith 
782bd481603SBarry Smith   Concepts: preallocation^Matrix
783bd481603SBarry Smith 
784bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
785bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
786bd481603SBarry Smith M*/
787d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
788c1ac3661SBarry Smith { PetscInt __i; \
789d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
790d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
791d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
792d3d32019SBarry Smith   }\
793d3d32019SBarry Smith }
794d3d32019SBarry Smith 
795bd481603SBarry Smith /*MC
79616371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
79716371a99SBarry Smith 
79816371a99SBarry Smith    Synopsis:
799f2ba6396SBarry Smith    #include "petscmat.h"
80016371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
80116371a99SBarry Smith 
80216371a99SBarry Smith    Not Collective
80316371a99SBarry Smith 
80416371a99SBarry Smith    Input Parameters:
80516371a99SBarry Smith .  A - matrix
80616371a99SBarry Smith .  row - row where values exist (must be local to this process)
80716371a99SBarry Smith .  ncols - number of columns
80816371a99SBarry Smith .  cols - columns with nonzeros
80916371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
81016371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
81116371a99SBarry Smith 
81216371a99SBarry Smith 
81316371a99SBarry Smith    Level: intermediate
81416371a99SBarry Smith 
81516371a99SBarry Smith    Notes:
8160598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
81716371a99SBarry Smith 
81816371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
81916371a99SBarry Smith 
82016371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
82116371a99SBarry Smith 
82216371a99SBarry Smith   Concepts: preallocation^Matrix
82316371a99SBarry Smith 
82416371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
825eabe889fSLisandro Dalcin           MatPreallocateSymmetricSetLocal()
82616371a99SBarry Smith M*/
8270298fd71SBarry 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);}
82816371a99SBarry Smith 
82916371a99SBarry Smith 
83016371a99SBarry Smith /*MC
8310d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
832bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
833bd481603SBarry Smith 
834bd481603SBarry Smith    Synopsis:
835f2ba6396SBarry Smith    #include "petscmat.h"
836c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
837bd481603SBarry Smith 
838bd481603SBarry Smith    Collective on MPI_Comm
839bd481603SBarry Smith 
840bd481603SBarry Smith    Input Parameters:
84116371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
842bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
843bd481603SBarry Smith 
844bd481603SBarry Smith 
845bd481603SBarry Smith    Level: intermediate
846bd481603SBarry Smith 
847bd481603SBarry Smith    Notes:
8480598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
849bd481603SBarry Smith 
850bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
851bd481603SBarry Smith 
8521620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
8531620fd73SBarry Smith 
854bd481603SBarry Smith   Concepts: preallocation^Matrix
855bd481603SBarry Smith 
856bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
857eabe889fSLisandro Dalcin           MatPreallocateSymmetricSetLocal()
858bd481603SBarry Smith M*/
859a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
8607c922b88SBarry Smith 
861bd481603SBarry Smith 
862bd481603SBarry Smith 
8637b80b807SBarry Smith /* Routines unique to particular data structures */
864014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
865698d4c6aSKris Buschelman 
866014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
867014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
868698d4c6aSKris Buschelman 
869014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
870014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
871014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
872014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
873014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
874014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
8757b80b807SBarry Smith 
876a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
877a96a251dSBarry Smith 
878014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
879014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
880014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
881273d9f13SBarry Smith 
882014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
883014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
884014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
885014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
886014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
887014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
888014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
889014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
890014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
891014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
8929230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
8939230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
894014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
895273d9f13SBarry Smith 
896014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt);
897014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
8981b807ce4Svictorle 
899014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
900014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
9012e8a6d31SBarry Smith 
902014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
9033a7fca6bSBarry Smith 
904014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
9057b80b807SBarry Smith /*
9067b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
90794b7f48cSBarry Smith   done through the KSP and PC interfaces.
9087b80b807SBarry Smith */
9097b80b807SBarry Smith 
91076bdecfbSBarry Smith /*J
9118f6c3df8SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering
912d9274352SBarry Smith 
913d9274352SBarry Smith    Level: beginner
914d9274352SBarry Smith 
915d9274352SBarry Smith .seealso: MatGetOrdering()
91676bdecfbSBarry Smith J*/
91719fd82e9SBarry Smith typedef const char* MatOrderingType;
9182692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
9192692d6eeSBarry Smith #define MATORDERINGND          "nd"
9202692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
9212692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
9222692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
9232692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
924312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
925b12f92e5SBarry Smith 
92619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
927140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
928bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
929607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void);
930014dd563SJed Brown PETSC_EXTERN PetscBool         MatOrderingRegisterAllCalled;
931140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList;
932d4fbbf0eSBarry Smith 
933014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
934a2ce50c7SBarry Smith 
935d91e6319SBarry Smith /*S
936d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
937d90ac83dSHong Zhang 
938d90ac83dSHong Zhang    Level: beginner
939d90ac83dSHong Zhang 
940d90ac83dSHong Zhang S*/
941d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
9426a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[];
9435e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];
944d90ac83dSHong Zhang 
945d90ac83dSHong Zhang /*S
9462401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
947b00f7748SHong Zhang 
94861cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
94961cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
950b00f7748SHong Zhang 
95115e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
952b00f7748SHong Zhang 
95361cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
95461cad860SBarry Smith 
955b00f7748SHong Zhang    Level: developer
956b00f7748SHong Zhang 
957d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
958d7d82daaSBarry Smith           MatFactorInfoInitialize()
959b00f7748SHong Zhang 
960b00f7748SHong Zhang S*/
961b00f7748SHong Zhang typedef struct {
96215e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
96385317021SBarry Smith   PetscReal     usedt;
96415e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
965b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
96615e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
96767e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
968348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
969bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
970bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
97115e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
972f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
973f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
97415e8a5b3SHong Zhang } MatFactorInfo;
975ffa6d0a5SLois Curfman McInnes 
976014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
977014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
978014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
979014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
980014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
981014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
982014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
983014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
984014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
985014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
986014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
987014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
988014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
989014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
990014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
991014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
992014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
993014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
994014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
9958ed539a5SBarry Smith 
996014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
997bb5a7306SBarry Smith 
998d91e6319SBarry Smith /*E
999d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1000bb1eb677SSatish Balay 
1001d91e6319SBarry Smith     Level: beginner
1002d91e6319SBarry Smith 
1003d9274352SBarry Smith    May be bitwise ORd together
1004d9274352SBarry Smith 
1005d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1006d91e6319SBarry Smith 
10074e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10084e7234bfSBarry Smith 
100941f059aeSBarry Smith .seealso: MatSOR()
1010d91e6319SBarry Smith E*/
1011ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1012ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1013ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
101484cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1015014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10168ed539a5SBarry Smith 
1017d4fbbf0eSBarry Smith /*
1018639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1019639f9d9dSBarry Smith */
1020b12f92e5SBarry Smith 
102176bdecfbSBarry Smith /*J
10228f6c3df8SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring
1023d9274352SBarry Smith 
1024d9274352SBarry Smith    Level: beginner
1025d9274352SBarry Smith 
10268f6c3df8SBarry Smith .seealso: MatGetColoring(), MatColoring
102776bdecfbSBarry Smith J*/
102819fd82e9SBarry Smith typedef const char* MatColoringType;
10292692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
10302692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
10312692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
10322692d6eeSBarry Smith #define MATCOLORINGID      "id"
1033b12f92e5SBarry Smith 
103419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetColoring(Mat,MatColoringType,ISColoring*);
1035bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
103630de9b25SBarry Smith 
1037014dd563SJed Brown PETSC_EXTERN PetscBool MatColoringRegisterAllCalled;
1038f1144a30SSatish Balay 
1039607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void);
1040014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1041639f9d9dSBarry Smith 
1042d9274352SBarry Smith /*S
1043d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1044d9274352SBarry Smith         and coloring
1045639f9d9dSBarry Smith 
1046d9274352SBarry Smith    Level: beginner
1047d9274352SBarry Smith 
1048d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1049d9274352SBarry Smith 
1050d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1051d9274352SBarry Smith S*/
1052e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1053639f9d9dSBarry Smith 
1054014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1055014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1056014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1057014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1058014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1059014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1060014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1061014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1062014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1063014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1064b1683b59SHong Zhang 
1065b1683b59SHong Zhang /*S
1066b9af6bddSHong Zhang      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1067b1683b59SHong Zhang 
1068b1683b59SHong Zhang    Level: beginner
1069b1683b59SHong Zhang 
1070b1683b59SHong Zhang   Concepts: coloring, sparse matrix product
1071b1683b59SHong Zhang 
1072b9af6bddSHong Zhang .seealso:  MatTransposeColoringCreate()
1073b1683b59SHong Zhang S*/
1074b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1075b1683b59SHong Zhang 
1076014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1077014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1078014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1079014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1080b1683b59SHong Zhang 
1081639f9d9dSBarry Smith /*
10820752156aSBarry Smith     These routines are for partitioning matrices: currently used only
10833eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
10840752156aSBarry Smith */
1085ca161407SBarry Smith 
1086d9274352SBarry Smith /*S
1087d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1088d9274352SBarry Smith 
1089d9274352SBarry Smith    Level: beginner
1090d9274352SBarry Smith 
1091d9274352SBarry Smith   Concepts: partitioning
1092d9274352SBarry Smith 
1093743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1094d9274352SBarry Smith S*/
109591e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1096d9274352SBarry Smith 
109776bdecfbSBarry Smith /*J
10988f6c3df8SBarry Smith     MatPartitioningType - String with the name of a PETSc matrix partitioning
1099d9274352SBarry Smith 
1100d9274352SBarry Smith    Level: beginner
11010d04baf8SBarry Smith dm
1102b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
110376bdecfbSBarry Smith J*/
110419fd82e9SBarry Smith typedef const char* MatPartitioningType;
11052692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
11062692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
11072692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
11082692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
11092692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
111050d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
111171306c60Sjroman 
1112ca161407SBarry Smith 
1113014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
111419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1115014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1116014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1117014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1118014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1119014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1120014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
11212aabb6bbSBarry Smith 
1122bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
112330de9b25SBarry Smith 
1124014dd563SJed Brown PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled;
1125f1144a30SSatish Balay 
1126607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void);
11272bad1931SBarry Smith 
1128014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1129014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
113019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1131ca161407SBarry Smith 
1132014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1133014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
11340752156aSBarry Smith 
1135b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
11366a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1137b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
11386a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[];
1139b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
11406a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[];
1141b6956c12SJose Roman 
1142014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1143014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1144014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1145014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1146014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1147014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1148014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1149014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1150014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1151014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1152014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
115371306c60Sjroman 
115471306c60Sjroman #define MP_PARTY_OPT "opt"
115571306c60Sjroman #define MP_PARTY_LIN "lin"
115671306c60Sjroman #define MP_PARTY_SCA "sca"
115771306c60Sjroman #define MP_PARTY_RAN "ran"
115871306c60Sjroman #define MP_PARTY_GBF "gbf"
115971306c60Sjroman #define MP_PARTY_GCF "gcf"
116071306c60Sjroman #define MP_PARTY_BUB "bub"
116171306c60Sjroman #define MP_PARTY_DEF "def"
1162014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
116371306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
116471306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
116571306c60Sjroman #define MP_PARTY_NONE "no"
1166014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1167014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1168014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1169014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
117071306c60Sjroman 
117150d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
11726a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1173e0f1cffaSJose Roman 
1174014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
117871306c60Sjroman 
1179b43b03e9SMark F. Adams /*
1180b43b03e9SMark F. Adams     These routines are for coarsening matrices:
1181b43b03e9SMark F. Adams */
1182b43b03e9SMark F. Adams 
1183b43b03e9SMark F. Adams /*S
1184b43b03e9SMark F. Adams      MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix)
1185b43b03e9SMark F. Adams 
1186b43b03e9SMark F. Adams    Level: beginner
1187b43b03e9SMark F. Adams 
1188b43b03e9SMark F. Adams   Concepts: coarsen
1189b43b03e9SMark F. Adams 
1190b43b03e9SMark F. Adams .seealso:  MatCoarsenCreate), MatCoarsenType
1191b43b03e9SMark F. Adams S*/
1192b43b03e9SMark F. Adams typedef struct _p_MatCoarsen* MatCoarsen;
1193b43b03e9SMark F. Adams 
1194b43b03e9SMark F. Adams /*J
11958f6c3df8SBarry Smith     MatCoarsenType - String with the name of a PETSc matrix coarsen
1196b43b03e9SMark F. Adams 
1197b43b03e9SMark F. Adams    Level: beginner
11988f6c3df8SBarry Smith 
1199b43b03e9SMark F. Adams .seealso: MatCoarsenCreate(), MatCoarsen
1200b43b03e9SMark F. Adams J*/
120119fd82e9SBarry Smith typedef const char* MatCoarsenType;
1202b43b03e9SMark F. Adams #define MATCOARSENMIS  "mis"
12039057884aSMark F. Adams #define MATCOARSENHEM  "hem"
1204b43b03e9SMark F. Adams 
12050cbbd2e1SMark F. Adams /* linked list for aggregates */
120641b27cdeSMark F. Adams typedef struct _PetscCDIntNd{
120741b27cdeSMark F. Adams   struct _PetscCDIntNd *next;
12080cbbd2e1SMark F. Adams   PetscInt             gid;
120941b27cdeSMark F. Adams }PetscCDIntNd;
12100cbbd2e1SMark F. Adams 
12110cbbd2e1SMark F. Adams /* only used by node pool */
121241b27cdeSMark F. Adams typedef struct _PetscCDArrNd{
121341b27cdeSMark F. Adams   struct _PetscCDArrNd *next;
121441b27cdeSMark F. Adams   struct _PetscCDIntNd *array;
121541b27cdeSMark F. Adams }PetscCDArrNd;
12160cbbd2e1SMark F. Adams 
12170cbbd2e1SMark F. Adams typedef struct _PetscCoarsenData{
121882c86c8fSBarry Smith   PetscCDArrNd pool_list;  /* node pool */
121941b27cdeSMark F. Adams   PetscCDIntNd *new_node;
12200cbbd2e1SMark F. Adams   PetscInt     new_left;
12210cbbd2e1SMark F. Adams   PetscInt     chk_sz;
122241b27cdeSMark F. Adams   PetscCDIntNd *extra_nodes;
122382c86c8fSBarry Smith   PetscCDIntNd **array;  /* Array of lists */
12240cbbd2e1SMark F. Adams   PetscInt     size;
122582c86c8fSBarry Smith   Mat          mat;  /* cache a Mat for communication data */
12260cbbd2e1SMark F. Adams }PetscCoarsenData;
12270cbbd2e1SMark F. Adams 
1228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*);
122919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,MatCoarsenType);
1230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat);
1231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS);
1232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool);
1233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetVerbose(MatCoarsen,PetscInt);
1234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** );
1235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen);
1236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*);
1237b43b03e9SMark F. Adams 
1238bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],PetscErrorCode (*)(MatCoarsen));
1239b43b03e9SMark F. Adams 
1240014dd563SJed Brown PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled;
1241b43b03e9SMark F. Adams 
1242607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void);
1243b43b03e9SMark F. Adams 
1244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer);
1245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen);
124619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,MatCoarsenType*);
1247b43b03e9SMark F. Adams 
1248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1250591294e4SBarry Smith 
12510752156aSBarry Smith /*
12520a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1253d4fbbf0eSBarry Smith */
12541c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
12551c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
12561c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
12571c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
12581c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
12597c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
12607c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
12611c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
12621c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
12637c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
12647c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
12651c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
12661c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1267a714c06dSBarry Smith                MATOP_SOR=13,
12681c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
12691c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
12701c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
12711c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
12721c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
12731c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
12741c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
12751c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1276d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1277d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1278d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1279d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1280d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1281d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1282d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1283d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1284d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1285d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1286d519adbfSMatthew Knepley                MATOP_GET_ARRAY=32,
1287d519adbfSMatthew Knepley                MATOP_RESTORE_ARRAY=33,
1288d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1289d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1290d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1291d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1292d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1293d519adbfSMatthew Knepley                MATOP_AXPY=39,
1294d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1295d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1296d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1297d519adbfSMatthew Knepley                MATOP_COPY=43,
1298d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1299d519adbfSMatthew Knepley                MATOP_SCALE=45,
1300d519adbfSMatthew Knepley                MATOP_SHIFT=46,
130135153367SBarry Smith                MATOP_DIAGONAL_SET=47,
1302d519adbfSMatthew Knepley                MATOP_ILUDT_FACTOR=48,
1303d519adbfSMatthew Knepley                MATOP_SET_BLOCK_SIZE=49,
1304d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1305d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1306d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1307d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1308d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1309d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1310d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1311d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1312d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1313d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1314d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1315d519adbfSMatthew Knepley                MATOP_VIEW=61,
1316d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
13177bab7c10SHong Zhang                MATOP_MATMAT_MULT=63,
13187bab7c10SHong Zhang                MATOP_MATMAT_MULT_SYMBOLIC=64,
13197bab7c10SHong Zhang                MATOP_MATMAT_MULT_NUMERIC=65,
1320d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1321d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1322d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1323d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1324d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1325d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1326d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
1327bcaeba4dSBarry Smith                MATOP_PLACEHOLDER=73,
1328d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1329d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1330d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1331bda6c4cbSBarry Smith                MATOP_MULT_CONSTRAINED=77,
1332bda6c4cbSBarry Smith                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
1333d519adbfSMatthew Knepley                MATOP_PERMUTE_SPARSIFY=79,
1334d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1335d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1336d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1337d519adbfSMatthew Knepley                MATOP_LOAD=83,
1338d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1339d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1340d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1341820469bcSHong Zhang                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1342d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1343d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1344d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1345d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1346d519adbfSMatthew Knepley                MATOP_PTAP=92,
1347d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1348d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1349bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT=95,
1350bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1351bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
1352820469bcSHong Zhang                MATOP_DUMMY98=98,
1353820469bcSHong Zhang                MATOP_DUMMY99=99,
1354820469bcSHong Zhang                MATOP_DUMMY100=100,
1355820469bcSHong Zhang                MATOP_DUMMY101=101,
1356d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
1357d519adbfSMatthew Knepley                MATOP_SET_SIZES=103,
1358d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1359d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1360bda6c4cbSBarry Smith                MATOP_IMAGINARY_PART=106,
1361bda6c4cbSBarry Smith                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1362bda6c4cbSBarry Smith                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1363bda6c4cbSBarry Smith                MATOP_MAT_SOLVE=109,
13640e8d04f7SBarry Smith                MATOP_GET_REDUNDANT_MATRIX=110,
1365d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
13660e8d04f7SBarry Smith                MATOP_GET_COLUMN_VECTOR=112,
1367d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
13680e8d04f7SBarry Smith                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
136989c6957cSBarry Smith                MATOP_CREATE=115,
137089c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
13710e8d04f7SBarry Smith                MATOP_GET_LOCAL_SUB_MATRIX=117,
13720e8d04f7SBarry Smith                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1373eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
13740e8d04f7SBarry Smith                MATOP_HERMITIAN_TRANSPOSE=120,
13750e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
13760e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
13770e8d04f7SBarry Smith                MATOP_GET_MULTI_PROC_BLOCK=123,
13780e8d04f7SBarry Smith                MATOP_GET_COLUMN_NORMS=125,
13790e8d04f7SBarry Smith                MATOP_GET_SUB_MATRICES_PARALLE=128,
13802b8ad9a3SHong Zhang                MATOP_SET_VALUES_BATCH=129,
13810e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT=130,
13820e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
13830e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
13840e8d04f7SBarry Smith                MATOP_TRANSPOSE_COLORING_CREAT=133,
13850e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_SPT=134,
13860e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_DEN=135,
13870e8d04f7SBarry Smith                MATOP_RART=136,
13880e8d04f7SBarry Smith                MATOP_RART_SYMBOLIC=137,
13890e8d04f7SBarry Smith                MATOP_RART_NUMERIC=138,
1390e09a3074SHong Zhang                MATOP_SET_BLOCK_SIZES=139,
1391e09a3074SHong Zhang                MATOP_AYPX=140
1392fae171e0SBarry Smith              } MatOperation;
1393014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *);
1394014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1395014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1396014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1397112a2221SBarry Smith 
139890ace30eSBarry Smith /*
139990ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
140090ace30eSBarry Smith    stored in a universal format. By changing the format with
14017973ad04SBarry Smith    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
140290ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
140390ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1404f4403165SShri Abhyankar    read into matrices of the same type.
140590ace30eSBarry Smith */
140690ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
140790ace30eSBarry Smith 
1408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1409014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1410014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
14111f4e1ec7SBarry Smith 
1412d9274352SBarry Smith /*S
1413d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1414d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1415d9274352SBarry Smith 
1416f7a9e4ceSBarry Smith    Level: advanced
1417d9274352SBarry Smith 
1418d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1419d9274352SBarry Smith 
14206e1639daSBarry Smith   Users manual sections:
14216e1639daSBarry Smith .   sec_singular
14226e1639daSBarry Smith 
1423d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1424d9274352SBarry Smith S*/
142574637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1426d9274352SBarry Smith 
1427014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1428014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1429014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1430d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1431014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1432014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1433014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1434014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1435014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1436014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1437014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1438014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
143974637425SBarry Smith 
1440014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1441014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1442014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1443014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
14443f1d51d7SBarry Smith 
1445014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1446014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1447014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1448c4f061fbSSatish Balay 
1449014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*);
1450b0a32e0cSBarry Smith 
1451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
145204f1ad80SBarry Smith 
1453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace);
1459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1467e884886fSBarry Smith 
14686370053bSBarry Smith /*S
14696370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
14706370053bSBarry Smith               Jacobian vector products
1471e884886fSBarry Smith 
14726370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
14736370053bSBarry Smith 
14746370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
14756370053bSBarry Smith 
14766370053bSBarry Smith     Level: developer
14776370053bSBarry Smith 
14786370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
14796370053bSBarry Smith S*/
1480e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1481e884886fSBarry Smith 
148276bdecfbSBarry Smith /*J
1483e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1484e884886fSBarry Smith 
1485e884886fSBarry Smith    Level: beginner
1486e884886fSBarry Smith 
1487e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
148876bdecfbSBarry Smith J*/
148919fd82e9SBarry Smith typedef const char* MatMFFDType;
1490e884886fSBarry Smith #define MATMFFD_DS  "ds"
1491e884886fSBarry Smith #define MATMFFD_WP  "wp"
1492e884886fSBarry Smith 
149319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1494bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));
1495e884886fSBarry Smith 
1496607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegisterAll(void);
1497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool );
1499e884886fSBarry Smith 
1500014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1501014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
15027dbadf16SMatthew Knepley 
150397969023SHong Zhang /*
150497969023SHong Zhang    PETSc interface to MUMPS
150597969023SHong Zhang */
150697969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1508b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
150997969023SHong Zhang #endif
151023a5497aSJed Brown 
1511d954db57SHong Zhang /*
1512d954db57SHong Zhang    PETSc interface to SUPERLU
1513d954db57SHong Zhang */
1514d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1516d954db57SHong Zhang #endif
1517d954db57SHong Zhang 
1518aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA
15193f9c0db1SPaul Mullowney /*E
1520e057df02SPaul Mullowney     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
15212692e278SPaul Mullowney     matrices.
1522e057df02SPaul Mullowney 
1523e057df02SPaul Mullowney     Not Collective
1524e057df02SPaul Mullowney 
15258468deeeSKarl Rupp +   MAT_CUSPARSE_CSR - Compressed Sparse Row
15262692e278SPaul Mullowney .   MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later).
15272692e278SPaul Mullowney -   MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).
1528e057df02SPaul Mullowney 
1529e057df02SPaul Mullowney     Level: intermediate
1530e057df02SPaul Mullowney 
1531e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1532e057df02SPaul Mullowney 
1533e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1534e057df02SPaul Mullowney E*/
1535e057df02SPaul Mullowney 
15363f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;
1537e057df02SPaul Mullowney 
1538e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1539e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
1540e057df02SPaul Mullowney 
15413f9c0db1SPaul Mullowney /*E
1542e057df02SPaul Mullowney     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
15432692e278SPaul Mullowney     matrices whose operation should use a particular storage format.
1544e057df02SPaul Mullowney 
1545e057df02SPaul Mullowney     Not Collective
1546e057df02SPaul Mullowney 
15478468deeeSKarl Rupp +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
15488468deeeSKarl Rupp .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
15498468deeeSKarl Rupp .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
15508468deeeSKarl Rupp .   MAT_CUSPARSE_SOLVE - sets the storage format for the triangular factors in the serial (single GPU) MatSolve
15518468deeeSKarl Rupp -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices
1552e057df02SPaul Mullowney 
1553e057df02SPaul Mullowney     Level: intermediate
1554e057df02SPaul Mullowney 
1555e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1556e057df02SPaul Mullowney E*/
15573f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_SOLVE, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;
1558e057df02SPaul Mullowney 
155910e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
156010e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1561e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
15629ae82921SPaul Mullowney #endif
15639ae82921SPaul Mullowney 
156490c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
1565014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1567e057df02SPaul Mullowney 
15683f9c0db1SPaul Mullowney /*E
1569e057df02SPaul Mullowney     MatCUSPStorageFormat - indicates the storage format for CUSP (GPU)
15708468deeeSKarl Rupp     matrices. Requires the txpetscgpu package to use. Configure with
15718468deeeSKarl Rupp     --download-txpetscgpu to build/install petsc with the txpetscgpu library.
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)
15928468deeeSKarl Rupp     matrices whose operation should use a particular storage format. Requires
15938468deeeSKarl Rupp     the txpetscgpu package to use. Configure with --download-txpetscgpu to
15948468deeeSKarl Rupp     build/install petsc with the txpetscgpu library.
1595e057df02SPaul Mullowney 
1596e057df02SPaul Mullowney     Not Collective
1597e057df02SPaul Mullowney 
15988468deeeSKarl Rupp +   MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
15998468deeeSKarl Rupp .   MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
16008468deeeSKarl Rupp .   MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
16018468deeeSKarl Rupp .   MAT_CUSP_SOLVE - sets the storage format for the triangular factors in the serial (single GPU) MatSolve
16028468deeeSKarl Rupp -   MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices
1603e057df02SPaul Mullowney 
1604e057df02SPaul Mullowney     Level: intermediate
1605e057df02SPaul Mullowney 
1606e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1607e057df02SPaul Mullowney 
1608e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat
1609e057df02SPaul Mullowney E*/
16103f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_SOLVE, MAT_CUSP_ALL} MatCUSPFormatOperation;
1611e057df02SPaul Mullowney 
1612e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat);
1613e057df02SPaul Mullowney #endif
161490c53902SBarry Smith 
1615d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL)
1616d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1617d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1618d67ff14aSKarl Rupp #endif
1619d67ff14aSKarl Rupp 
162054efbe56SHong Zhang /*
162154efbe56SHong Zhang    PETSc interface to FFTW
162254efbe56SHong Zhang */
162354efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
1624014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1625014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
1626014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*);
162754efbe56SHong Zhang #endif
162854efbe56SHong Zhang 
1629db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL)
1630607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
1631db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
1632db31f6deSJed Brown #endif
1633db31f6deSJed Brown 
1634014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1635014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1636014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1637014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1638014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1639014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
164019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1641014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1642014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1643d8588912SDave May 
16444325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
16454325cce7SMatthew G Knepley 
164623a5497aSJed Brown #endif
1647