xref: /petsc/include/petscmat.h (revision 38f409eb36822ff3419b5ed9e4be1699a2d6fe1d)
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 *);
216*38f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
2184cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]);
219dfb205c3SBarry Smith 
220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*);
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
226c0cdd4a1SDahai Guo 
227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
231c0cdd4a1SDahai Guo 
232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2396d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
2416d7c1e57SBarry Smith 
24219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*);
243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
244dedccee8SHong Zhang 
245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*);
247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS);
248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
2491472f72bSBarry Smith 
250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]);
2511d6018f0SLisandro Dalcin 
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
25421c89e3eSBarry Smith 
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
260713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
26199cafbc1SBarry Smith 
2628ed539a5SBarry Smith /* ------------------------------------------------------------*/
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
267014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
26873a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);
26984cb2905SBarry Smith 
2702ef4de8bSBarry Smith /*S
2712ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
272be479b30SJed Brown         column of a matrix as indexed on an associated grid.
273be479b30SJed Brown 
274be479b30SJed Brown    Fortran usage is different, see MatSetValuesStencil() for details.
2752ef4de8bSBarry Smith 
2762ef4de8bSBarry Smith    Level: beginner
2772ef4de8bSBarry Smith 
2782ef4de8bSBarry Smith   Concepts: matrix; linear operator
2792ef4de8bSBarry Smith 
280be479b30SJed Brown .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil()
2812ef4de8bSBarry Smith S*/
282435da068SBarry Smith typedef struct {
283c1ac3661SBarry Smith   PetscInt k,j,i,c;
284435da068SBarry Smith } MatStencil;
2852ef4de8bSBarry Smith 
286014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
287014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
288014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
289435da068SBarry Smith 
290014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetColoring(Mat,ISColoring);
291014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*);
2923a7fca6bSBarry Smith 
293d91e6319SBarry Smith /*E
294d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
295d91e6319SBarry Smith      to continue to add values to it
296d91e6319SBarry Smith 
297d91e6319SBarry Smith     Level: beginner
298d91e6319SBarry Smith 
299d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
300d91e6319SBarry Smith E*/
3016d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
302014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
303014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
304014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);
3054f9c727eSBarry Smith 
3061d73ed98SBarry Smith 
30730de9b25SBarry Smith 
308d91e6319SBarry Smith /*E
309d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
310d91e6319SBarry Smith 
311d91e6319SBarry Smith     Level: beginner
312d91e6319SBarry Smith 
3130a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
314d91e6319SBarry Smith 
315382164b0SBarry Smith    Developer Notes: Entries that are negative need not be called collectively by all processes.
316382164b0SBarry Smith 
317d91e6319SBarry Smith .seealso: MatSetOption()
318d91e6319SBarry Smith E*/
319382164b0SBarry Smith typedef enum {MAT_OPTION_MIN = -8,
320382164b0SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR = -7,
321382164b0SBarry Smith               MAT_NO_OFF_PROC_ZERO_ROWS = -6,
322382164b0SBarry Smith               MAT_NO_OFF_PROC_ENTRIES = -5,
323382164b0SBarry Smith               MAT_UNUSED_NONZERO_LOCATION_ERR = -4,
324382164b0SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR = -3,
325382164b0SBarry Smith               MAT_ROW_ORIENTED = -2,
326382164b0SBarry Smith               MAT_NEW_NONZERO_LOCATIONS = -1,
327382164b0SBarry Smith               MAT_SYMMETRIC = 1,
328382164b0SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC = 2,
329382164b0SBarry Smith               MAT_NEW_DIAGONALS = 3,
330382164b0SBarry Smith               MAT_IGNORE_OFF_PROC_ENTRIES = 4,
331382164b0SBarry Smith               MAT_USE_HASH_TABLE = 5,
332382164b0SBarry Smith               MAT_KEEP_NONZERO_PATTERN = 6,
333382164b0SBarry Smith               MAT_IGNORE_ZERO_ENTRIES = 7,
334382164b0SBarry Smith               MAT_USE_INODES = 8,
335382164b0SBarry Smith               MAT_HERMITIAN = 9,
336382164b0SBarry Smith               MAT_SYMMETRY_ETERNAL = 10,
33711e456e1SBarry Smith               MAT_DUMMY = 11,
338382164b0SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR = 12,
339382164b0SBarry Smith               MAT_ERROR_LOWER_TRIANGULAR = 13,
340382164b0SBarry Smith               MAT_GETROW_UPPERTRIANGULAR = 14,
341382164b0SBarry Smith               MAT_SPD = 15,
342382164b0SBarry Smith               MAT_OPTION_MAX = 16} MatOption;
343e057df02SPaul Mullowney 
344014dd563SJed Brown PETSC_EXTERN const char *MatOptions[];
345014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool );
34619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);
34784cb2905SBarry Smith 
348014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
349014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
350014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
351014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
3568c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
3578c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
3588c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
3598c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
363014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
364014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNThreads(Mat,PetscInt);
365014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNThreads(Mat,PetscInt*);
3661620fd73SBarry Smith 
367014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
368014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
369014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
379f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);
380f5edf698SHong Zhang 
381d91e6319SBarry Smith /*E
382d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
383d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
384d91e6319SBarry Smith 
385d91e6319SBarry Smith     Level: beginner
386d91e6319SBarry Smith 
387d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
388d91e6319SBarry Smith 
38970dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
39070dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
39170dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
39270dcbbb9SBarry Smith 
393d91e6319SBarry Smith .seealso: MatDuplicate()
394d91e6319SBarry Smith E*/
39570dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
3962e8a6d31SBarry Smith 
39719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
398014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);
39994a9d846SBarry Smith 
400cb5b572fSBarry Smith 
401014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
402014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
403014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
404014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
405014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
406014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
407014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
409014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
4107b80b807SBarry Smith 
4111a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4121a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
4131a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4141a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
415d4fbbf0eSBarry Smith 
416d91e6319SBarry Smith /*S
417d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
418d91e6319SBarry Smith 
419d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
420d91e6319SBarry Smith 
421d91e6319SBarry Smith    Level: intermediate
422d91e6319SBarry Smith 
423d91e6319SBarry Smith   Concepts: matrix^nonzero information
424d91e6319SBarry Smith 
425d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
426d91e6319SBarry Smith S*/
4274e220ebcSLois Curfman McInnes typedef struct {
428b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
429b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
430b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
431b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
432b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
433b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
434b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4354e220ebcSLois Curfman McInnes } MatInfo;
4364e220ebcSLois Curfman McInnes 
437d9274352SBarry Smith /*E
438d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
439d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
440d9274352SBarry Smith 
441d9274352SBarry Smith     Level: beginner
442d9274352SBarry Smith 
443d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
444d9274352SBarry Smith 
445d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
446d9274352SBarry Smith E*/
4477b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
448014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
449014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
450014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
452014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat *);
458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);
460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool *);
461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
4657b80b807SBarry Smith 
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal *);
467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *);
468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
4757b80b807SBarry Smith 
476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
482566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);
4835ef9f2a5SBarry Smith 
484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);
4927b80b807SBarry Smith 
493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,Mat*);
495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJNumeric(MPI_Comm,Mat,PetscInt,Mat);
496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5038efafbd8SBarry Smith 
504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5057b80b807SBarry Smith 
506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat);
50922440eb1SKris Buschelman 
5107bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
5117bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultSymbolic(Mat,Mat,Mat,PetscReal,Mat*);
5127bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultNumeric(Mat,Mat,Mat,Mat);
5137bab7c10SHong Zhang 
514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat);
517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat);
52022440eb1SKris Buschelman 
521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat);
524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat);
527bc011b1eSHong Zhang 
528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);
5301c741599SHong Zhang 
531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);
5337b80b807SBarry Smith 
534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
544052efed2SBarry Smith 
545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
54790f02eecSBarry Smith 
548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecs(Mat,Vec*,Vec*);
552b2bf6370SHong Zhang PETSC_EXTERN PetscErrorCode MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
555bd481603SBarry Smith 
556bd481603SBarry Smith /*MC
5571620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5581620fd73SBarry Smith 
5591620fd73SBarry Smith    Not collective
5601620fd73SBarry Smith 
5611620fd73SBarry Smith    Input Parameters:
5621620fd73SBarry Smith +  m - the matrix
5631620fd73SBarry Smith .  row - the row location of the entry
5641620fd73SBarry Smith .  col - the column location of the entry
5651620fd73SBarry Smith .  value - the value to insert
5661620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5671620fd73SBarry Smith 
5681620fd73SBarry Smith    Notes:
5691620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
5701620fd73SBarry Smith    values simultaneously if possible.
5711620fd73SBarry Smith 
5721620fd73SBarry Smith    Level: beginner
5731620fd73SBarry Smith 
5741620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
5751620fd73SBarry Smith M*/
5761620fd73SBarry 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);}
5771620fd73SBarry Smith 
5781620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
5791620fd73SBarry Smith 
5801620fd73SBarry 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);}
5811620fd73SBarry Smith 
5821620fd73SBarry Smith /*MC
5830d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
584bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
585bd481603SBarry Smith 
586bd481603SBarry Smith    Synopsis:
587f2ba6396SBarry Smith    #include "petscmat.h"
588c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
589bd481603SBarry Smith 
590bd481603SBarry Smith    Collective on MPI_Comm
591bd481603SBarry Smith 
592bd481603SBarry Smith    Input Parameters:
593bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
594859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
595859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
596bd481603SBarry Smith 
597bd481603SBarry Smith    Output Parameters:
598bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
599bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
600bd481603SBarry Smith 
601bd481603SBarry Smith 
602bd481603SBarry Smith    Level: intermediate
603bd481603SBarry Smith 
604bd481603SBarry Smith    Notes:
6050598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
606bd481603SBarry Smith 
6071d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
608bd481603SBarry Smith 
609bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
610bd481603SBarry Smith 
6111620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6121620fd73SBarry Smith 
613bd481603SBarry Smith   Concepts: preallocation^Matrix
614bd481603SBarry Smith 
615bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
616bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
617bd481603SBarry Smith M*/
618c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6197c922b88SBarry Smith { \
620a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
621a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
622a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
623eabe889fSLisandro Dalcin   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr); __start = 0; __end = __start; \
624c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
625a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6267c922b88SBarry Smith 
627bd481603SBarry Smith /*MC
628bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
629bd481603SBarry Smith        inserted using a local number of the rows and columns
630bd481603SBarry Smith 
631bd481603SBarry Smith    Synopsis:
632f2ba6396SBarry Smith    #include "petscmat.h"
633c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
634bd481603SBarry Smith 
635bd481603SBarry Smith    Not Collective
636bd481603SBarry Smith 
637bd481603SBarry Smith    Input Parameters:
638784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
639bd481603SBarry Smith .  nrows - the number of rows indicated
6401d73ed98SBarry Smith .  rows - the indices of the rows
641784ac674SJed Brown .  cmap - the column mapping from local to global numbering
642bd481603SBarry Smith .  ncols - the number of columns in the matrix
643bd481603SBarry Smith .  cols - the columns indicated
644bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
645bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
646bd481603SBarry Smith 
647bd481603SBarry Smith 
648bd481603SBarry Smith    Level: intermediate
649bd481603SBarry Smith 
650bd481603SBarry Smith    Notes:
6510598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
652bd481603SBarry Smith 
6531d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
654bd481603SBarry Smith 
655bd481603SBarry Smith   Concepts: preallocation^Matrix
656bd481603SBarry Smith 
657bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
658bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
659bd481603SBarry Smith M*/
660784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
661c4f061fbSSatish Balay {\
662c1ac3661SBarry Smith   PetscInt __l;\
663784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
664784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
665c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
666ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
667c4f061fbSSatish Balay   }\
668c4f061fbSSatish Balay }
669c4f061fbSSatish Balay 
670bd481603SBarry Smith /*MC
671bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
672bd481603SBarry Smith        inserted using a local number of the rows and columns
673bd481603SBarry Smith 
674bd481603SBarry Smith    Synopsis:
675f2ba6396SBarry Smith    #include "petscmat.h"
676c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
677bd481603SBarry Smith 
678bd481603SBarry Smith    Not Collective
679bd481603SBarry Smith 
680bd481603SBarry Smith    Input Parameters:
681bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
682bd481603SBarry Smith .  nrows - the number of rows indicated
6831d73ed98SBarry Smith .  rows - the indices of the rows
684bd481603SBarry Smith .  ncols - the number of columns in the matrix
685bd481603SBarry Smith .  cols - the columns indicated
686bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
687bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
688bd481603SBarry Smith 
689bd481603SBarry Smith 
690bd481603SBarry Smith    Level: intermediate
691bd481603SBarry Smith 
692bd481603SBarry Smith    Notes:
6930598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
694bd481603SBarry Smith 
695bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
696bd481603SBarry Smith 
697bd481603SBarry Smith   Concepts: preallocation^Matrix
698bd481603SBarry Smith 
699bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
700bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
701bd481603SBarry Smith M*/
702d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
703d3d32019SBarry Smith {\
704c1ac3661SBarry Smith   PetscInt __l;\
705d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
706d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
707d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
708d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
709d3d32019SBarry Smith   }\
710d3d32019SBarry Smith }
711d3d32019SBarry Smith 
712bd481603SBarry Smith /*MC
713bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
714bd481603SBarry Smith        inserted using a local number of the rows and columns
715bd481603SBarry Smith 
716bd481603SBarry Smith    Synopsis:
717f2ba6396SBarry Smith    #include "petscmat.h"
718c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
719bd481603SBarry Smith 
720bd481603SBarry Smith    Not Collective
721bd481603SBarry Smith 
722bd481603SBarry Smith    Input Parameters:
72364075487SBarry Smith +  row - the row
724bd481603SBarry Smith .  ncols - the number of columns in the matrix
725a50f8bf6SHong Zhang -  cols - the columns indicated
726a50f8bf6SHong Zhang 
727a50f8bf6SHong Zhang    Output Parameters:
728a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
729bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
730bd481603SBarry Smith 
731bd481603SBarry Smith 
732bd481603SBarry Smith    Level: intermediate
733bd481603SBarry Smith 
734bd481603SBarry Smith    Notes:
7350598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
736bd481603SBarry Smith 
737bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
738bd481603SBarry Smith 
7391620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7401620fd73SBarry Smith 
741bd481603SBarry Smith   Concepts: preallocation^Matrix
742bd481603SBarry Smith 
743bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
744bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
745bd481603SBarry Smith M*/
746c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
747c1ac3661SBarry Smith { PetscInt __i; \
748e32f2f54SBarry 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);\
749e32f2f54SBarry 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);\
7507c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
75164075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
7527cc688d7SBarry Smith     else dnz[row - __rstart]++;\
7537c922b88SBarry Smith   }\
7547c922b88SBarry Smith }
7557c922b88SBarry Smith 
756bd481603SBarry Smith /*MC
757bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
758bd481603SBarry Smith        inserted using a local number of the rows and columns
759bd481603SBarry Smith 
760bd481603SBarry Smith    Synopsis:
761f2ba6396SBarry Smith    #include "petscmat.h"
762c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
763bd481603SBarry Smith 
764bd481603SBarry Smith    Not Collective
765bd481603SBarry Smith 
766bd481603SBarry Smith    Input Parameters:
767bd481603SBarry Smith +  nrows - the number of rows indicated
7681d73ed98SBarry Smith .  rows - the indices of the rows
769bd481603SBarry Smith .  ncols - the number of columns in the matrix
770bd481603SBarry Smith .  cols - the columns indicated
771bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
772bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
773bd481603SBarry Smith 
774bd481603SBarry Smith 
775bd481603SBarry Smith    Level: intermediate
776bd481603SBarry Smith 
777bd481603SBarry Smith    Notes:
7780598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
779bd481603SBarry Smith 
780bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
781bd481603SBarry Smith 
7821620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7831620fd73SBarry Smith 
784bd481603SBarry Smith   Concepts: preallocation^Matrix
785bd481603SBarry Smith 
786bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
787bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
788bd481603SBarry Smith M*/
789d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
790c1ac3661SBarry Smith { PetscInt __i; \
791d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
792d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
793d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
794d3d32019SBarry Smith   }\
795d3d32019SBarry Smith }
796d3d32019SBarry Smith 
797bd481603SBarry Smith /*MC
79816371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
79916371a99SBarry Smith 
80016371a99SBarry Smith    Synopsis:
801f2ba6396SBarry Smith    #include "petscmat.h"
80216371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
80316371a99SBarry Smith 
80416371a99SBarry Smith    Not Collective
80516371a99SBarry Smith 
80616371a99SBarry Smith    Input Parameters:
80716371a99SBarry Smith .  A - matrix
80816371a99SBarry Smith .  row - row where values exist (must be local to this process)
80916371a99SBarry Smith .  ncols - number of columns
81016371a99SBarry Smith .  cols - columns with nonzeros
81116371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
81216371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
81316371a99SBarry Smith 
81416371a99SBarry Smith 
81516371a99SBarry Smith    Level: intermediate
81616371a99SBarry Smith 
81716371a99SBarry Smith    Notes:
8180598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
81916371a99SBarry Smith 
82016371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
82116371a99SBarry Smith 
82216371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
82316371a99SBarry Smith 
82416371a99SBarry Smith   Concepts: preallocation^Matrix
82516371a99SBarry Smith 
82616371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
827eabe889fSLisandro Dalcin           MatPreallocateSymmetricSetLocal()
82816371a99SBarry Smith M*/
8290298fd71SBarry 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);}
83016371a99SBarry Smith 
83116371a99SBarry Smith 
83216371a99SBarry Smith /*MC
8330d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
834bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
835bd481603SBarry Smith 
836bd481603SBarry Smith    Synopsis:
837f2ba6396SBarry Smith    #include "petscmat.h"
838c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
839bd481603SBarry Smith 
840bd481603SBarry Smith    Collective on MPI_Comm
841bd481603SBarry Smith 
842bd481603SBarry Smith    Input Parameters:
84316371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
844bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
845bd481603SBarry Smith 
846bd481603SBarry Smith 
847bd481603SBarry Smith    Level: intermediate
848bd481603SBarry Smith 
849bd481603SBarry Smith    Notes:
8500598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
851bd481603SBarry Smith 
852bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
853bd481603SBarry Smith 
8541620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
8551620fd73SBarry Smith 
856bd481603SBarry Smith   Concepts: preallocation^Matrix
857bd481603SBarry Smith 
858bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
859eabe889fSLisandro Dalcin           MatPreallocateSymmetricSetLocal()
860bd481603SBarry Smith M*/
861a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
8627c922b88SBarry Smith 
863bd481603SBarry Smith 
864bd481603SBarry Smith 
8657b80b807SBarry Smith /* Routines unique to particular data structures */
866014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
867698d4c6aSKris Buschelman 
868014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
869014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
870698d4c6aSKris Buschelman 
871014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
872014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
873014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
874014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
875014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
876014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
8777b80b807SBarry Smith 
878a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
879a96a251dSBarry Smith 
880014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
881014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
882014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
883273d9f13SBarry Smith 
884014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
885014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
886014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
887014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
888014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
889014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
890014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
891014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
892014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
893014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
8949230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
8959230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
896014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
897273d9f13SBarry Smith 
898014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt);
899014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
9001b807ce4Svictorle 
901014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
902014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
9032e8a6d31SBarry Smith 
904014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
9053a7fca6bSBarry Smith 
906014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
9077b80b807SBarry Smith /*
9087b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
90994b7f48cSBarry Smith   done through the KSP and PC interfaces.
9107b80b807SBarry Smith */
9117b80b807SBarry Smith 
91276bdecfbSBarry Smith /*J
9138f6c3df8SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering
914d9274352SBarry Smith 
915d9274352SBarry Smith    Level: beginner
916d9274352SBarry Smith 
917d9274352SBarry Smith .seealso: MatGetOrdering()
91876bdecfbSBarry Smith J*/
91919fd82e9SBarry Smith typedef const char* MatOrderingType;
9202692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
9212692d6eeSBarry Smith #define MATORDERINGND          "nd"
9222692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
9232692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
9242692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
9252692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
926312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
927b12f92e5SBarry Smith 
92819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
929140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
930bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
931607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void);
932014dd563SJed Brown PETSC_EXTERN PetscBool         MatOrderingRegisterAllCalled;
933140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList;
934d4fbbf0eSBarry Smith 
935014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
936a2ce50c7SBarry Smith 
937d91e6319SBarry Smith /*S
938d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
939d90ac83dSHong Zhang 
940d90ac83dSHong Zhang    Level: beginner
941d90ac83dSHong Zhang 
942d90ac83dSHong Zhang S*/
943d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
9446a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[];
9455e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];
946d90ac83dSHong Zhang 
947d90ac83dSHong Zhang /*S
9482401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
949b00f7748SHong Zhang 
95061cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
95161cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
952b00f7748SHong Zhang 
95315e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
954b00f7748SHong Zhang 
95561cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
95661cad860SBarry Smith 
957b00f7748SHong Zhang    Level: developer
958b00f7748SHong Zhang 
959d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
960d7d82daaSBarry Smith           MatFactorInfoInitialize()
961b00f7748SHong Zhang 
962b00f7748SHong Zhang S*/
963b00f7748SHong Zhang typedef struct {
96415e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
96585317021SBarry Smith   PetscReal     usedt;
96615e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
967b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
96815e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
96967e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
970348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
971bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
972bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
97315e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
974f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
975f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
97615e8a5b3SHong Zhang } MatFactorInfo;
977ffa6d0a5SLois Curfman McInnes 
978014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
979014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
980014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
981014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
982014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
983014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
984014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
985014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
986014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
987014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
988014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
989014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
990014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
991014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
992014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
993014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
994014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
995014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
996014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
9978ed539a5SBarry Smith 
998014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
999bb5a7306SBarry Smith 
1000d91e6319SBarry Smith /*E
1001d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1002bb1eb677SSatish Balay 
1003d91e6319SBarry Smith     Level: beginner
1004d91e6319SBarry Smith 
1005d9274352SBarry Smith    May be bitwise ORd together
1006d9274352SBarry Smith 
1007d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1008d91e6319SBarry Smith 
10094e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10104e7234bfSBarry Smith 
101141f059aeSBarry Smith .seealso: MatSOR()
1012d91e6319SBarry Smith E*/
1013ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1014ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1015ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
101684cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1017014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10188ed539a5SBarry Smith 
1019d4fbbf0eSBarry Smith /*
1020639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1021639f9d9dSBarry Smith */
1022b12f92e5SBarry Smith 
102376bdecfbSBarry Smith /*J
10248f6c3df8SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring
1025d9274352SBarry Smith 
1026d9274352SBarry Smith    Level: beginner
1027d9274352SBarry Smith 
10288f6c3df8SBarry Smith .seealso: MatGetColoring(), MatColoring
102976bdecfbSBarry Smith J*/
103019fd82e9SBarry Smith typedef const char* MatColoringType;
10312692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
10322692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
10332692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
10342692d6eeSBarry Smith #define MATCOLORINGID      "id"
1035b12f92e5SBarry Smith 
103619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetColoring(Mat,MatColoringType,ISColoring*);
1037bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
103830de9b25SBarry Smith 
1039014dd563SJed Brown PETSC_EXTERN PetscBool MatColoringRegisterAllCalled;
1040f1144a30SSatish Balay 
1041607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void);
1042014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1043639f9d9dSBarry Smith 
1044d9274352SBarry Smith /*S
1045d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1046d9274352SBarry Smith         and coloring
1047639f9d9dSBarry Smith 
1048d9274352SBarry Smith    Level: beginner
1049d9274352SBarry Smith 
1050d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1051d9274352SBarry Smith 
1052d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1053d9274352SBarry Smith S*/
1054e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1055639f9d9dSBarry Smith 
1056014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1057014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1058014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1059014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1060014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1061014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1062014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1063014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1064014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1065014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1066f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1067f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1068f86b9fbaSHong Zhang 
1069b1683b59SHong Zhang 
1070b1683b59SHong Zhang /*S
1071b9af6bddSHong Zhang      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1072b1683b59SHong Zhang 
1073b1683b59SHong Zhang    Level: beginner
1074b1683b59SHong Zhang 
1075b1683b59SHong Zhang   Concepts: coloring, sparse matrix product
1076b1683b59SHong Zhang 
1077b9af6bddSHong Zhang .seealso:  MatTransposeColoringCreate()
1078b1683b59SHong Zhang S*/
1079b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1080b1683b59SHong Zhang 
1081014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1082014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1083014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1084014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1085b1683b59SHong Zhang 
1086639f9d9dSBarry Smith /*
10870752156aSBarry Smith     These routines are for partitioning matrices: currently used only
10883eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
10890752156aSBarry Smith */
1090ca161407SBarry Smith 
1091d9274352SBarry Smith /*S
1092d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1093d9274352SBarry Smith 
1094d9274352SBarry Smith    Level: beginner
1095d9274352SBarry Smith 
1096d9274352SBarry Smith   Concepts: partitioning
1097d9274352SBarry Smith 
1098743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1099d9274352SBarry Smith S*/
110091e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1101d9274352SBarry Smith 
110276bdecfbSBarry Smith /*J
11038f6c3df8SBarry Smith     MatPartitioningType - String with the name of a PETSc matrix partitioning
1104d9274352SBarry Smith 
1105d9274352SBarry Smith    Level: beginner
11060d04baf8SBarry Smith dm
1107b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
110876bdecfbSBarry Smith J*/
110919fd82e9SBarry Smith typedef const char* MatPartitioningType;
11102692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
11112692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
11122692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
11132692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
11142692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
111550d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
111671306c60Sjroman 
1117ca161407SBarry Smith 
1118014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
111919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1120014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1121014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1122014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1123014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1124014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1125014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
11262aabb6bbSBarry Smith 
1127bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
112830de9b25SBarry Smith 
1129014dd563SJed Brown PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled;
1130f1144a30SSatish Balay 
1131607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void);
11322bad1931SBarry Smith 
1133014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1134014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
113519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1136ca161407SBarry Smith 
1137014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1138014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
11390752156aSBarry Smith 
1140b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
11416a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1142b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
11436a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[];
1144b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
11456a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[];
1146b6956c12SJose Roman 
1147014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1148014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1149014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1150014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1151014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1152014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1153014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1154014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1155014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1156014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1157014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
115871306c60Sjroman 
115971306c60Sjroman #define MP_PARTY_OPT "opt"
116071306c60Sjroman #define MP_PARTY_LIN "lin"
116171306c60Sjroman #define MP_PARTY_SCA "sca"
116271306c60Sjroman #define MP_PARTY_RAN "ran"
116371306c60Sjroman #define MP_PARTY_GBF "gbf"
116471306c60Sjroman #define MP_PARTY_GCF "gcf"
116571306c60Sjroman #define MP_PARTY_BUB "bub"
116671306c60Sjroman #define MP_PARTY_DEF "def"
1167014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
116871306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
116971306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
117071306c60Sjroman #define MP_PARTY_NONE "no"
1171014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1172014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1173014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1174014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
117571306c60Sjroman 
117650d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
11776a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1178e0f1cffaSJose Roman 
1179014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1180014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1181014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1182014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
118371306c60Sjroman 
1184b43b03e9SMark F. Adams /*
1185b43b03e9SMark F. Adams     These routines are for coarsening matrices:
1186b43b03e9SMark F. Adams */
1187b43b03e9SMark F. Adams 
1188b43b03e9SMark F. Adams /*S
1189b43b03e9SMark F. Adams      MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix)
1190b43b03e9SMark F. Adams 
1191b43b03e9SMark F. Adams    Level: beginner
1192b43b03e9SMark F. Adams 
1193b43b03e9SMark F. Adams   Concepts: coarsen
1194b43b03e9SMark F. Adams 
1195b43b03e9SMark F. Adams .seealso:  MatCoarsenCreate), MatCoarsenType
1196b43b03e9SMark F. Adams S*/
1197b43b03e9SMark F. Adams typedef struct _p_MatCoarsen* MatCoarsen;
1198b43b03e9SMark F. Adams 
1199b43b03e9SMark F. Adams /*J
12008f6c3df8SBarry Smith     MatCoarsenType - String with the name of a PETSc matrix coarsen
1201b43b03e9SMark F. Adams 
1202b43b03e9SMark F. Adams    Level: beginner
12038f6c3df8SBarry Smith 
1204b43b03e9SMark F. Adams .seealso: MatCoarsenCreate(), MatCoarsen
1205b43b03e9SMark F. Adams J*/
120619fd82e9SBarry Smith typedef const char* MatCoarsenType;
1207b43b03e9SMark F. Adams #define MATCOARSENMIS  "mis"
12089057884aSMark F. Adams #define MATCOARSENHEM  "hem"
1209b43b03e9SMark F. Adams 
12100cbbd2e1SMark F. Adams /* linked list for aggregates */
121141b27cdeSMark F. Adams typedef struct _PetscCDIntNd{
121241b27cdeSMark F. Adams   struct _PetscCDIntNd *next;
12130cbbd2e1SMark F. Adams   PetscInt             gid;
121441b27cdeSMark F. Adams }PetscCDIntNd;
12150cbbd2e1SMark F. Adams 
12160cbbd2e1SMark F. Adams /* only used by node pool */
121741b27cdeSMark F. Adams typedef struct _PetscCDArrNd{
121841b27cdeSMark F. Adams   struct _PetscCDArrNd *next;
121941b27cdeSMark F. Adams   struct _PetscCDIntNd *array;
122041b27cdeSMark F. Adams }PetscCDArrNd;
12210cbbd2e1SMark F. Adams 
12220cbbd2e1SMark F. Adams typedef struct _PetscCoarsenData{
122382c86c8fSBarry Smith   PetscCDArrNd pool_list;  /* node pool */
122441b27cdeSMark F. Adams   PetscCDIntNd *new_node;
12250cbbd2e1SMark F. Adams   PetscInt     new_left;
12260cbbd2e1SMark F. Adams   PetscInt     chk_sz;
122741b27cdeSMark F. Adams   PetscCDIntNd *extra_nodes;
122882c86c8fSBarry Smith   PetscCDIntNd **array;  /* Array of lists */
12290cbbd2e1SMark F. Adams   PetscInt     size;
123082c86c8fSBarry Smith   Mat          mat;  /* cache a Mat for communication data */
12310cbbd2e1SMark F. Adams }PetscCoarsenData;
12320cbbd2e1SMark F. Adams 
1233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*);
123419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,MatCoarsenType);
1235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat);
1236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS);
1237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool);
1238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetVerbose(MatCoarsen,PetscInt);
1239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** );
1240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen);
1241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*);
1242b43b03e9SMark F. Adams 
1243bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],PetscErrorCode (*)(MatCoarsen));
1244b43b03e9SMark F. Adams 
1245014dd563SJed Brown PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled;
1246b43b03e9SMark F. Adams 
1247607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void);
1248b43b03e9SMark F. Adams 
1249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer);
1250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen);
125119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,MatCoarsenType*);
1252b43b03e9SMark F. Adams 
1253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1255591294e4SBarry Smith 
12560752156aSBarry Smith /*
12570a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1258d4fbbf0eSBarry Smith */
12591c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
12601c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
12611c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
12621c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
12631c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
12647c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
12657c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
12661c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
12671c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
12687c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
12697c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
12701c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
12711c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1272a714c06dSBarry Smith                MATOP_SOR=13,
12731c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
12741c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
12751c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
12761c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
12771c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
12781c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
12791c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
12801c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1281d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1282d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1283d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1284d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1285d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1286d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1287d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1288d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1289d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1290d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
12919d39f69dSJed Brown                /* MATOP_PLACEHOLDER_32=32, */
12929d39f69dSJed Brown                /* MATOP_PLACEHOLDER_33=33, */
1293d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1294d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1295d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1296d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1297d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1298d519adbfSMatthew Knepley                MATOP_AXPY=39,
1299d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1300d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1301d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1302d519adbfSMatthew Knepley                MATOP_COPY=43,
1303d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1304d519adbfSMatthew Knepley                MATOP_SCALE=45,
1305d519adbfSMatthew Knepley                MATOP_SHIFT=46,
130635153367SBarry Smith                MATOP_DIAGONAL_SET=47,
13079d39f69dSJed Brown                MATOP_ZERO_ROWS_COLUMNS=48,
13089d39f69dSJed Brown                MATOP_SET_RANDOM=49,
1309d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1310d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1311d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1312d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1313d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1314d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1315d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1316d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1317d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1318d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1319d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1320d519adbfSMatthew Knepley                MATOP_VIEW=61,
1321d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
13227bab7c10SHong Zhang                MATOP_MATMAT_MULT=63,
13237bab7c10SHong Zhang                MATOP_MATMAT_MULT_SYMBOLIC=64,
13247bab7c10SHong Zhang                MATOP_MATMAT_MULT_NUMERIC=65,
1325d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1326d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1327d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1328d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1329d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1330d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1331d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
13329d39f69dSJed Brown                /* MATOP_PLACEHOLDER_73=73, */
1333d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1334d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1335d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1336bda6c4cbSBarry Smith                MATOP_MULT_CONSTRAINED=77,
1337bda6c4cbSBarry Smith                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
13389d39f69dSJed Brown                MATOP_FIND_ZERO_DIAGONALS=79,
1339d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1340d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1341d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1342d519adbfSMatthew Knepley                MATOP_LOAD=83,
1343d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1344d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1345d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1346820469bcSHong Zhang                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1347d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1348d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1349d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1350d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1351d519adbfSMatthew Knepley                MATOP_PTAP=92,
1352d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1353d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1354bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT=95,
1355bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1356bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
13579d39f69dSJed Brown                /* MATOP_PLACEHOLDER_98=98, */
13589d39f69dSJed Brown                /* MATOP_PLACEHOLDER_99=99, */
13599d39f69dSJed Brown                /* MATOP_PLACEHOLDER_100=100, */
13609d39f69dSJed Brown                /* MATOP_PLACEHOLDER_101=101, */
1361d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
13629d39f69dSJed Brown                /* MATOP_PLACEHOLDER_103=103, */
1363d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1364d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1365bda6c4cbSBarry Smith                MATOP_IMAGINARY_PART=106,
1366bda6c4cbSBarry Smith                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1367bda6c4cbSBarry Smith                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1368bda6c4cbSBarry Smith                MATOP_MAT_SOLVE=109,
13690e8d04f7SBarry Smith                MATOP_GET_REDUNDANT_MATRIX=110,
1370d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
13710e8d04f7SBarry Smith                MATOP_GET_COLUMN_VECTOR=112,
1372d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
13730e8d04f7SBarry Smith                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
137489c6957cSBarry Smith                MATOP_CREATE=115,
137589c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
13760e8d04f7SBarry Smith                MATOP_GET_LOCAL_SUB_MATRIX=117,
13770e8d04f7SBarry Smith                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1378eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
13790e8d04f7SBarry Smith                MATOP_HERMITIAN_TRANSPOSE=120,
13800e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
13810e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
13820e8d04f7SBarry Smith                MATOP_GET_MULTI_PROC_BLOCK=123,
13839d39f69dSJed Brown                MATOP_FIND_NONZERO_ROWS=124,
13840e8d04f7SBarry Smith                MATOP_GET_COLUMN_NORMS=125,
13859d39f69dSJed Brown                MATOP_INVERT_BLOCK_DIAGONAL=126,
13869d39f69dSJed Brown                /* MATOP_PLACEHOLDER_127=127, */
13870e8d04f7SBarry Smith                MATOP_GET_SUB_MATRICES_PARALLE=128,
13882b8ad9a3SHong Zhang                MATOP_SET_VALUES_BATCH=129,
13890e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT=130,
13900e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
13910e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
13920e8d04f7SBarry Smith                MATOP_TRANSPOSE_COLORING_CREAT=133,
13930e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_SPT=134,
13940e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_DEN=135,
13950e8d04f7SBarry Smith                MATOP_RART=136,
13960e8d04f7SBarry Smith                MATOP_RART_SYMBOLIC=137,
13970e8d04f7SBarry Smith                MATOP_RART_NUMERIC=138,
1398e09a3074SHong Zhang                MATOP_SET_BLOCK_SIZES=139,
1399f9426fe0SMark Adams                MATOP_AYPX=140,
14001919a2e2SJed Brown                MATOP_RESIDUAL=141,
14011919a2e2SJed Brown                MATOP_FDCOLORING_SETUP= 142
1402fae171e0SBarry Smith              } MatOperation;
1403014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *);
1404014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1405014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1406014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1407112a2221SBarry Smith 
140890ace30eSBarry Smith /*
140990ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
141090ace30eSBarry Smith    stored in a universal format. By changing the format with
14117973ad04SBarry Smith    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
141290ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
141390ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1414f4403165SShri Abhyankar    read into matrices of the same type.
141590ace30eSBarry Smith */
141690ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
141790ace30eSBarry Smith 
1418014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1419014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1420014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
14211f4e1ec7SBarry Smith 
1422d9274352SBarry Smith /*S
1423d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1424d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1425d9274352SBarry Smith 
1426f7a9e4ceSBarry Smith    Level: advanced
1427d9274352SBarry Smith 
1428d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1429d9274352SBarry Smith 
14306e1639daSBarry Smith   Users manual sections:
14316e1639daSBarry Smith .   sec_singular
14326e1639daSBarry Smith 
1433d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1434d9274352SBarry Smith S*/
143574637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1436d9274352SBarry Smith 
1437014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1438014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1439014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1440d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1441014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1442014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1443014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1444014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1445014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1446014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1447014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1448014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
144974637425SBarry Smith 
1450014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1452014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
14543f1d51d7SBarry Smith 
1455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1458c4f061fbSSatish Balay 
1459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*);
1460b0a32e0cSBarry Smith 
1461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
146204f1ad80SBarry Smith 
1463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace);
1469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1477e884886fSBarry Smith 
14786370053bSBarry Smith /*S
14796370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
14806370053bSBarry Smith               Jacobian vector products
1481e884886fSBarry Smith 
14826370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
14836370053bSBarry Smith 
14846370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
14856370053bSBarry Smith 
14866370053bSBarry Smith     Level: developer
14876370053bSBarry Smith 
14886370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
14896370053bSBarry Smith S*/
1490e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1491e884886fSBarry Smith 
149276bdecfbSBarry Smith /*J
1493e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1494e884886fSBarry Smith 
1495e884886fSBarry Smith    Level: beginner
1496e884886fSBarry Smith 
1497e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
149876bdecfbSBarry Smith J*/
149919fd82e9SBarry Smith typedef const char* MatMFFDType;
1500e884886fSBarry Smith #define MATMFFD_DS  "ds"
1501e884886fSBarry Smith #define MATMFFD_WP  "wp"
1502e884886fSBarry Smith 
150319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1504bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));
1505e884886fSBarry Smith 
1506607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegisterAll(void);
1507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool );
1509e884886fSBarry Smith 
1510014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1511014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
15127dbadf16SMatthew Knepley 
151397969023SHong Zhang /*
151497969023SHong Zhang    PETSc interface to MUMPS
151597969023SHong Zhang */
151697969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1518b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
151997969023SHong Zhang #endif
152023a5497aSJed Brown 
1521d954db57SHong Zhang /*
1522d954db57SHong Zhang    PETSc interface to SUPERLU
1523d954db57SHong Zhang */
1524d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1526d954db57SHong Zhang #endif
1527d954db57SHong Zhang 
1528aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA
15293f9c0db1SPaul Mullowney /*E
1530e057df02SPaul Mullowney     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
15312692e278SPaul Mullowney     matrices.
1532e057df02SPaul Mullowney 
1533e057df02SPaul Mullowney     Not Collective
1534e057df02SPaul Mullowney 
15358468deeeSKarl Rupp +   MAT_CUSPARSE_CSR - Compressed Sparse Row
15362692e278SPaul Mullowney .   MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later).
15372692e278SPaul Mullowney -   MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).
1538e057df02SPaul Mullowney 
1539e057df02SPaul Mullowney     Level: intermediate
1540e057df02SPaul Mullowney 
1541e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1542e057df02SPaul Mullowney 
1543e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1544e057df02SPaul Mullowney E*/
1545e057df02SPaul Mullowney 
15463f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;
1547e057df02SPaul Mullowney 
1548e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1549e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
1550e057df02SPaul Mullowney 
15513f9c0db1SPaul Mullowney /*E
1552e057df02SPaul Mullowney     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
15532692e278SPaul Mullowney     matrices whose operation should use a particular storage format.
1554e057df02SPaul Mullowney 
1555e057df02SPaul Mullowney     Not Collective
1556e057df02SPaul Mullowney 
15578468deeeSKarl Rupp +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
15588468deeeSKarl Rupp .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
15598468deeeSKarl Rupp .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
15608468deeeSKarl Rupp -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices
1561e057df02SPaul Mullowney 
1562e057df02SPaul Mullowney     Level: intermediate
1563e057df02SPaul Mullowney 
1564e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1565e057df02SPaul Mullowney E*/
156636d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;
1567e057df02SPaul Mullowney 
156810e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
156910e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1570e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
15719ae82921SPaul Mullowney #endif
15729ae82921SPaul Mullowney 
157390c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
1574014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1575014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1576e057df02SPaul Mullowney 
15773f9c0db1SPaul Mullowney /*E
1578e057df02SPaul Mullowney     MatCUSPStorageFormat - indicates the storage format for CUSP (GPU)
157936d62e41SPaul Mullowney     matrices.
1580e057df02SPaul Mullowney 
1581e057df02SPaul Mullowney     Not Collective
1582e057df02SPaul Mullowney 
15838468deeeSKarl Rupp +   MAT_CUSP_CSR - Compressed Sparse Row
15848468deeeSKarl Rupp .   MAT_CUSP_DIA - Diagonal
15858468deeeSKarl Rupp -   MAT_CUSP_ELL - Ellpack
1586e057df02SPaul Mullowney 
1587e057df02SPaul Mullowney     Level: intermediate
1588e057df02SPaul Mullowney 
1589e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1590e057df02SPaul Mullowney 
1591e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
1592e057df02SPaul Mullowney E*/
15933f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat;
1594e057df02SPaul Mullowney 
1595e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1596e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[];
1597e057df02SPaul Mullowney 
15983f9c0db1SPaul Mullowney /*E
1599e057df02SPaul Mullowney     MatCUSPFormatOperation - indicates the operation of CUSP (GPU)
160036d62e41SPaul Mullowney     matrices whose operation should use a particular storage format.
1601e057df02SPaul Mullowney 
1602e057df02SPaul Mullowney     Not Collective
1603e057df02SPaul Mullowney 
16048468deeeSKarl Rupp +   MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
16058468deeeSKarl Rupp .   MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
16068468deeeSKarl Rupp .   MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
16078468deeeSKarl Rupp -   MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices
1608e057df02SPaul Mullowney 
1609e057df02SPaul Mullowney     Level: intermediate
1610e057df02SPaul Mullowney 
1611e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1612e057df02SPaul Mullowney 
1613e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat
1614e057df02SPaul Mullowney E*/
161536d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation;
1616e057df02SPaul Mullowney 
1617e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat);
1618e057df02SPaul Mullowney #endif
161990c53902SBarry Smith 
1620d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL)
1621d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1622d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1623d67ff14aSKarl Rupp #endif
1624d67ff14aSKarl Rupp 
162554efbe56SHong Zhang /*
162654efbe56SHong Zhang    PETSc interface to FFTW
162754efbe56SHong Zhang */
162854efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
1629014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1630014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
1631014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*);
163254efbe56SHong Zhang #endif
163354efbe56SHong Zhang 
1634db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL)
1635607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
1636db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
1637db31f6deSJed Brown #endif
1638db31f6deSJed Brown 
1639014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1640014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1641014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1642014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1643014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1644014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
164519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1646014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1647014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1648d8588912SDave May 
16494325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
1650e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);
16514325cce7SMatthew G Knepley 
165223a5497aSJed Brown #endif
1653