xref: /petsc/include/petscmat.h (revision bd04181c04923b6fd56517bde870143f8e9ec3da)
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;
139335efc43SPeter Brune PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID;
140014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID;
141014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
142014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID;
143014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID;
144014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID;
145014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID;
146c06d978dSMatthew Knepley 
147ceb03754SKris Buschelman /*E
148ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
149d6eff37eSBarry Smith      or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate
150d6eff37eSBarry Smith      that the input matrix is to be replaced with the converted matrix.
151ceb03754SKris Buschelman 
152ceb03754SKris Buschelman     Level: beginner
153ceb03754SKris Buschelman 
154ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
155ceb03754SKris Buschelman 
1560c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
157ceb03754SKris Buschelman E*/
158dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse;
159ceb03754SKris Buschelman 
1605494a064SHong Zhang /*E
1615494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
162829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1635494a064SHong Zhang 
1645494a064SHong Zhang     Level: beginner
1655494a064SHong Zhang 
166829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1675494a064SHong Zhang E*/
1685494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1695494a064SHong Zhang 
170607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void);
171c06d978dSMatthew Knepley 
172014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*);
173014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
17419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType);
175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat);
176ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
177607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatRegisterAll(void);
178bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat));
179014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]);
180014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]);
181014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]);
182014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]);
183f69a0ea3SMatthew Knepley 
184014dd563SJed Brown PETSC_EXTERN PetscBool         MatRegisterAllCalled;
185140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList;
186140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList;
187140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList;
188140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatCoarsenList;
18928988994SBarry Smith 
1903b224e63SBarry Smith /*E
1913b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
1923b224e63SBarry Smith 
1933b224e63SBarry Smith     Level: beginner
1943b224e63SBarry Smith 
1953b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1963b224e63SBarry Smith 
1973b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
1983b224e63SBarry Smith E*/
199214bc6a2SJed Brown typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN,SAME_PRECONDITIONER} MatStructure;
2003b224e63SBarry Smith 
201014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
202014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
206014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2078d7a6e47SBarry Smith 
208014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
211d21a29f3SJed Brown 
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
214d21a29f3SJed Brown 
215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
21738f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
218014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
2194cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]);
220dfb205c3SBarry Smith 
221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*);
224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
227c0cdd4a1SDahai Guo 
228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
232c0cdd4a1SDahai Guo 
233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2406d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
2426d7c1e57SBarry Smith 
24319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*);
244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
245dedccee8SHong Zhang 
246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*);
248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS);
249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
2501472f72bSBarry Smith 
251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]);
2521d6018f0SLisandro Dalcin 
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
25521c89e3eSBarry Smith 
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
261713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
26299cafbc1SBarry Smith 
2638ed539a5SBarry Smith /* ------------------------------------------------------------*/
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
267014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
268014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
26973a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);
27084cb2905SBarry Smith 
2712ef4de8bSBarry Smith /*S
2722ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
273be479b30SJed Brown         column of a matrix as indexed on an associated grid.
274be479b30SJed Brown 
275be479b30SJed Brown    Fortran usage is different, see MatSetValuesStencil() for details.
2762ef4de8bSBarry Smith 
2772ef4de8bSBarry Smith    Level: beginner
2782ef4de8bSBarry Smith 
2792ef4de8bSBarry Smith   Concepts: matrix; linear operator
2802ef4de8bSBarry Smith 
281be479b30SJed Brown .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil()
2822ef4de8bSBarry Smith S*/
283435da068SBarry Smith typedef struct {
284c1ac3661SBarry Smith   PetscInt k,j,i,c;
285435da068SBarry Smith } MatStencil;
2862ef4de8bSBarry Smith 
287014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
288014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
289014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
290435da068SBarry Smith 
291014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetColoring(Mat,ISColoring);
292014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*);
2933a7fca6bSBarry Smith 
294d91e6319SBarry Smith /*E
295d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
296d91e6319SBarry Smith      to continue to add values to it
297d91e6319SBarry Smith 
298d91e6319SBarry Smith     Level: beginner
299d91e6319SBarry Smith 
300d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
301d91e6319SBarry Smith E*/
3026d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
303014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
304014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
305014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);
3064f9c727eSBarry Smith 
3071d73ed98SBarry Smith 
30830de9b25SBarry Smith 
309d91e6319SBarry Smith /*E
310d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
311d91e6319SBarry Smith 
312d91e6319SBarry Smith     Level: beginner
313d91e6319SBarry Smith 
3140a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
315d91e6319SBarry Smith 
316382164b0SBarry Smith    Developer Notes: Entries that are negative need not be called collectively by all processes.
317382164b0SBarry Smith 
318d91e6319SBarry Smith .seealso: MatSetOption()
319d91e6319SBarry Smith E*/
320382164b0SBarry Smith typedef enum {MAT_OPTION_MIN = -8,
321382164b0SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR = -7,
322382164b0SBarry Smith               MAT_NO_OFF_PROC_ZERO_ROWS = -6,
323382164b0SBarry Smith               MAT_NO_OFF_PROC_ENTRIES = -5,
324382164b0SBarry Smith               MAT_UNUSED_NONZERO_LOCATION_ERR = -4,
325382164b0SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR = -3,
326382164b0SBarry Smith               MAT_ROW_ORIENTED = -2,
327382164b0SBarry Smith               MAT_NEW_NONZERO_LOCATIONS = -1,
328382164b0SBarry Smith               MAT_SYMMETRIC = 1,
329382164b0SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC = 2,
330382164b0SBarry Smith               MAT_NEW_DIAGONALS = 3,
331382164b0SBarry Smith               MAT_IGNORE_OFF_PROC_ENTRIES = 4,
332382164b0SBarry Smith               MAT_USE_HASH_TABLE = 5,
333382164b0SBarry Smith               MAT_KEEP_NONZERO_PATTERN = 6,
334382164b0SBarry Smith               MAT_IGNORE_ZERO_ENTRIES = 7,
335382164b0SBarry Smith               MAT_USE_INODES = 8,
336382164b0SBarry Smith               MAT_HERMITIAN = 9,
337382164b0SBarry Smith               MAT_SYMMETRY_ETERNAL = 10,
33811e456e1SBarry Smith               MAT_DUMMY = 11,
339382164b0SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR = 12,
340382164b0SBarry Smith               MAT_ERROR_LOWER_TRIANGULAR = 13,
341382164b0SBarry Smith               MAT_GETROW_UPPERTRIANGULAR = 14,
342382164b0SBarry Smith               MAT_SPD = 15,
343382164b0SBarry Smith               MAT_OPTION_MAX = 16} MatOption;
344e057df02SPaul Mullowney 
345014dd563SJed Brown PETSC_EXTERN const char *MatOptions[];
346014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool );
34719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);
34884cb2905SBarry Smith 
349014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
350014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
351014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
3578c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
3588c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
359*bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3608c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
3618c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
363014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
364014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
365014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
366014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNThreads(Mat,PetscInt);
367014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNThreads(Mat,PetscInt*);
3681620fd73SBarry Smith 
369014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
381f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);
382f5edf698SHong Zhang 
383d91e6319SBarry Smith /*E
384d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
385d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
386d91e6319SBarry Smith 
387d91e6319SBarry Smith     Level: beginner
388d91e6319SBarry Smith 
389d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
390d91e6319SBarry Smith 
39170dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
39270dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
39370dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
39470dcbbb9SBarry Smith 
395d91e6319SBarry Smith .seealso: MatDuplicate()
396d91e6319SBarry Smith E*/
39770dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
3982e8a6d31SBarry Smith 
39919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
400014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);
40194a9d846SBarry Smith 
402cb5b572fSBarry Smith 
403014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
404014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
405014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
406014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
407014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
409014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
410014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
411014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
4127b80b807SBarry Smith 
4131a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4141a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
4151a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4161a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
417d4fbbf0eSBarry Smith 
418d91e6319SBarry Smith /*S
419d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
420d91e6319SBarry Smith 
421d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
422d91e6319SBarry Smith 
423d91e6319SBarry Smith    Level: intermediate
424d91e6319SBarry Smith 
425d91e6319SBarry Smith   Concepts: matrix^nonzero information
426d91e6319SBarry Smith 
427d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
428d91e6319SBarry Smith S*/
4294e220ebcSLois Curfman McInnes typedef struct {
430b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
431b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
432b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
433b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
434b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
435b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
436b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4374e220ebcSLois Curfman McInnes } MatInfo;
4384e220ebcSLois Curfman McInnes 
439d9274352SBarry Smith /*E
440d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
441d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
442d9274352SBarry Smith 
443d9274352SBarry Smith     Level: beginner
444d9274352SBarry Smith 
445d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
446d9274352SBarry Smith 
447d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
448d9274352SBarry Smith E*/
4497b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
450014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
452014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat *);
460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);
462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool *);
463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
4677b80b807SBarry Smith 
468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal *);
469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *);
470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
4777b80b807SBarry Smith 
478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
484566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);
4855ef9f2a5SBarry Smith 
486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);
4947b80b807SBarry Smith 
495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,Mat*);
497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJNumeric(MPI_Comm,Mat,PetscInt,Mat);
498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5058efafbd8SBarry Smith 
506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5077b80b807SBarry Smith 
508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat);
51122440eb1SKris Buschelman 
5127bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
5137bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultSymbolic(Mat,Mat,Mat,PetscReal,Mat*);
5147bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultNumeric(Mat,Mat,Mat,Mat);
5157bab7c10SHong Zhang 
516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat);
519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat);
52222440eb1SKris Buschelman 
523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat);
526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat);
529bc011b1eSHong Zhang 
530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);
5321c741599SHong Zhang 
533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);
5357b80b807SBarry Smith 
536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
546052efed2SBarry Smith 
547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
54990f02eecSBarry Smith 
550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecs(Mat,Vec*,Vec*);
554b2bf6370SHong Zhang PETSC_EXTERN PetscErrorCode MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
5573a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
558bd481603SBarry Smith 
559bd481603SBarry Smith /*MC
5601620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5611620fd73SBarry Smith 
5621620fd73SBarry Smith    Not collective
5631620fd73SBarry Smith 
5641620fd73SBarry Smith    Input Parameters:
5651620fd73SBarry Smith +  m - the matrix
5661620fd73SBarry Smith .  row - the row location of the entry
5671620fd73SBarry Smith .  col - the column location of the entry
5681620fd73SBarry Smith .  value - the value to insert
5691620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5701620fd73SBarry Smith 
5711620fd73SBarry Smith    Notes:
5721620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
5731620fd73SBarry Smith    values simultaneously if possible.
5741620fd73SBarry Smith 
5751620fd73SBarry Smith    Level: beginner
5761620fd73SBarry Smith 
5771620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
5781620fd73SBarry Smith M*/
5791620fd73SBarry 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);}
5801620fd73SBarry Smith 
5811620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
5821620fd73SBarry Smith 
5831620fd73SBarry 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);}
5841620fd73SBarry Smith 
5851620fd73SBarry Smith /*MC
5860d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
587bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
588bd481603SBarry Smith 
589bd481603SBarry Smith    Synopsis:
590f2ba6396SBarry Smith    #include "petscmat.h"
591c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
592bd481603SBarry Smith 
593bd481603SBarry Smith    Collective on MPI_Comm
594bd481603SBarry Smith 
595bd481603SBarry Smith    Input Parameters:
596bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
597859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
598859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
599bd481603SBarry Smith 
600bd481603SBarry Smith    Output Parameters:
601bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
602bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
603bd481603SBarry Smith 
604bd481603SBarry Smith 
605bd481603SBarry Smith    Level: intermediate
606bd481603SBarry Smith 
607bd481603SBarry Smith    Notes:
6080598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
609bd481603SBarry Smith 
6101d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
611bd481603SBarry Smith 
612bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
613bd481603SBarry Smith 
6141620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6151620fd73SBarry Smith 
616bd481603SBarry Smith   Concepts: preallocation^Matrix
617bd481603SBarry Smith 
618bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
619bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
620bd481603SBarry Smith M*/
621c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6227c922b88SBarry Smith { \
623a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
6241795a4d1SJed Brown   _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \
6251795a4d1SJed Brown   __start = 0; __end = __start;                                         \
626c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
627a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6287c922b88SBarry Smith 
629bd481603SBarry Smith /*MC
630bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
631bd481603SBarry Smith        inserted using a local number of the rows and columns
632bd481603SBarry Smith 
633bd481603SBarry Smith    Synopsis:
634f2ba6396SBarry Smith    #include "petscmat.h"
635c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
636bd481603SBarry Smith 
637bd481603SBarry Smith    Not Collective
638bd481603SBarry Smith 
639bd481603SBarry Smith    Input Parameters:
640784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
641bd481603SBarry Smith .  nrows - the number of rows indicated
6421d73ed98SBarry Smith .  rows - the indices of the rows
643784ac674SJed Brown .  cmap - the column mapping from local to global numbering
644bd481603SBarry Smith .  ncols - the number of columns in the matrix
645bd481603SBarry Smith .  cols - the columns indicated
646bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
647bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
648bd481603SBarry Smith 
649bd481603SBarry Smith 
650bd481603SBarry Smith    Level: intermediate
651bd481603SBarry Smith 
652bd481603SBarry Smith    Notes:
6530598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
654bd481603SBarry Smith 
6551d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
656bd481603SBarry Smith 
657bd481603SBarry Smith   Concepts: preallocation^Matrix
658bd481603SBarry Smith 
659bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
660bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
661bd481603SBarry Smith M*/
662784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
663c4f061fbSSatish Balay {\
664c1ac3661SBarry Smith   PetscInt __l;\
665784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
666784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
667c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
668ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
669c4f061fbSSatish Balay   }\
670c4f061fbSSatish Balay }
671c4f061fbSSatish Balay 
672bd481603SBarry Smith /*MC
673bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
674bd481603SBarry Smith        inserted using a local number of the rows and columns
675bd481603SBarry Smith 
676bd481603SBarry Smith    Synopsis:
677f2ba6396SBarry Smith    #include "petscmat.h"
678c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
679bd481603SBarry Smith 
680bd481603SBarry Smith    Not Collective
681bd481603SBarry Smith 
682bd481603SBarry Smith    Input Parameters:
683bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
684bd481603SBarry Smith .  nrows - the number of rows indicated
6851d73ed98SBarry Smith .  rows - the indices of the rows
686bd481603SBarry Smith .  ncols - the number of columns in the matrix
687bd481603SBarry Smith .  cols - the columns indicated
688bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
689bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
690bd481603SBarry Smith 
691bd481603SBarry Smith 
692bd481603SBarry Smith    Level: intermediate
693bd481603SBarry Smith 
694bd481603SBarry Smith    Notes:
6950598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
696bd481603SBarry Smith 
697bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
698bd481603SBarry Smith 
699bd481603SBarry Smith   Concepts: preallocation^Matrix
700bd481603SBarry Smith 
701bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
702bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
703bd481603SBarry Smith M*/
704d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
705d3d32019SBarry Smith {\
706c1ac3661SBarry Smith   PetscInt __l;\
707d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
708d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
709d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
710d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
711d3d32019SBarry Smith   }\
712d3d32019SBarry Smith }
713d3d32019SBarry Smith 
714bd481603SBarry Smith /*MC
715bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
716bd481603SBarry Smith        inserted using a local number of the rows and columns
717bd481603SBarry Smith 
718bd481603SBarry Smith    Synopsis:
719f2ba6396SBarry Smith    #include "petscmat.h"
720c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
721bd481603SBarry Smith 
722bd481603SBarry Smith    Not Collective
723bd481603SBarry Smith 
724bd481603SBarry Smith    Input Parameters:
72564075487SBarry Smith +  row - the row
726bd481603SBarry Smith .  ncols - the number of columns in the matrix
727a50f8bf6SHong Zhang -  cols - the columns indicated
728a50f8bf6SHong Zhang 
729a50f8bf6SHong Zhang    Output Parameters:
730a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
731bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
732bd481603SBarry Smith 
733bd481603SBarry Smith 
734bd481603SBarry Smith    Level: intermediate
735bd481603SBarry Smith 
736bd481603SBarry Smith    Notes:
7370598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
738bd481603SBarry Smith 
739bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
740bd481603SBarry Smith 
7411620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7421620fd73SBarry Smith 
743bd481603SBarry Smith   Concepts: preallocation^Matrix
744bd481603SBarry Smith 
745bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
746bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
747bd481603SBarry Smith M*/
748c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
749c1ac3661SBarry Smith { PetscInt __i; \
750e32f2f54SBarry 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);\
751e32f2f54SBarry 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);\
7527c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
75364075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
7547cc688d7SBarry Smith     else dnz[row - __rstart]++;\
7557c922b88SBarry Smith   }\
7567c922b88SBarry Smith }
7577c922b88SBarry Smith 
758bd481603SBarry Smith /*MC
759bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
760bd481603SBarry Smith        inserted using a local number of the rows and columns
761bd481603SBarry Smith 
762bd481603SBarry Smith    Synopsis:
763f2ba6396SBarry Smith    #include "petscmat.h"
764c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
765bd481603SBarry Smith 
766bd481603SBarry Smith    Not Collective
767bd481603SBarry Smith 
768bd481603SBarry Smith    Input Parameters:
769bd481603SBarry Smith +  nrows - the number of rows indicated
7701d73ed98SBarry Smith .  rows - the indices of the rows
771bd481603SBarry Smith .  ncols - the number of columns in the matrix
772bd481603SBarry Smith .  cols - the columns indicated
773bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
774bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
775bd481603SBarry Smith 
776bd481603SBarry Smith 
777bd481603SBarry Smith    Level: intermediate
778bd481603SBarry Smith 
779bd481603SBarry Smith    Notes:
7800598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
781bd481603SBarry Smith 
782bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
783bd481603SBarry Smith 
7841620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7851620fd73SBarry Smith 
786bd481603SBarry Smith   Concepts: preallocation^Matrix
787bd481603SBarry Smith 
788bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
789bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
790bd481603SBarry Smith M*/
791d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
792c1ac3661SBarry Smith { PetscInt __i; \
793d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
794d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
795d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
796d3d32019SBarry Smith   }\
797d3d32019SBarry Smith }
798d3d32019SBarry Smith 
799bd481603SBarry Smith /*MC
80016371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
80116371a99SBarry Smith 
80216371a99SBarry Smith    Synopsis:
803f2ba6396SBarry Smith    #include "petscmat.h"
80416371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
80516371a99SBarry Smith 
80616371a99SBarry Smith    Not Collective
80716371a99SBarry Smith 
80816371a99SBarry Smith    Input Parameters:
80916371a99SBarry Smith .  A - matrix
81016371a99SBarry Smith .  row - row where values exist (must be local to this process)
81116371a99SBarry Smith .  ncols - number of columns
81216371a99SBarry Smith .  cols - columns with nonzeros
81316371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
81416371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
81516371a99SBarry Smith 
81616371a99SBarry Smith 
81716371a99SBarry Smith    Level: intermediate
81816371a99SBarry Smith 
81916371a99SBarry Smith    Notes:
8200598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
82116371a99SBarry Smith 
82216371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
82316371a99SBarry Smith 
82416371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
82516371a99SBarry Smith 
82616371a99SBarry Smith   Concepts: preallocation^Matrix
82716371a99SBarry Smith 
82816371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
829eabe889fSLisandro Dalcin           MatPreallocateSymmetricSetLocal()
83016371a99SBarry Smith M*/
8310298fd71SBarry 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);}
83216371a99SBarry Smith 
83316371a99SBarry Smith 
83416371a99SBarry Smith /*MC
8350d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
836bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
837bd481603SBarry Smith 
838bd481603SBarry Smith    Synopsis:
839f2ba6396SBarry Smith    #include "petscmat.h"
840c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
841bd481603SBarry Smith 
842bd481603SBarry Smith    Collective on MPI_Comm
843bd481603SBarry Smith 
844bd481603SBarry Smith    Input Parameters:
84516371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
846bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
847bd481603SBarry Smith 
848bd481603SBarry Smith 
849bd481603SBarry Smith    Level: intermediate
850bd481603SBarry Smith 
851bd481603SBarry Smith    Notes:
8520598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
853bd481603SBarry Smith 
854bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
855bd481603SBarry Smith 
8561620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
8571620fd73SBarry Smith 
858bd481603SBarry Smith   Concepts: preallocation^Matrix
859bd481603SBarry Smith 
860bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
861eabe889fSLisandro Dalcin           MatPreallocateSymmetricSetLocal()
862bd481603SBarry Smith M*/
863a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
8647c922b88SBarry Smith 
865bd481603SBarry Smith 
866bd481603SBarry Smith 
8677b80b807SBarry Smith /* Routines unique to particular data structures */
868014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
869698d4c6aSKris Buschelman 
870014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
871014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
872698d4c6aSKris Buschelman 
873014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
874014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
875014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
876014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
877014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
878014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
8797b80b807SBarry Smith 
880a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
881a96a251dSBarry Smith 
882014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
883014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
884014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
885273d9f13SBarry Smith 
886014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
887014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
888014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
889014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
890014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
891014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
892014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
893014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
894014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
895014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
8969230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
8979230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
898014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
899273d9f13SBarry Smith 
900014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt);
901014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
9021b807ce4Svictorle 
903014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
904014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
9052e8a6d31SBarry Smith 
906014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
9073a7fca6bSBarry Smith 
908014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
9097b80b807SBarry Smith /*
9107b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
91194b7f48cSBarry Smith   done through the KSP and PC interfaces.
9127b80b807SBarry Smith */
9137b80b807SBarry Smith 
91476bdecfbSBarry Smith /*J
9158f6c3df8SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering
916d9274352SBarry Smith 
917d9274352SBarry Smith    Level: beginner
918d9274352SBarry Smith 
919d9274352SBarry Smith .seealso: MatGetOrdering()
92076bdecfbSBarry Smith J*/
92119fd82e9SBarry Smith typedef const char* MatOrderingType;
9222692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
9232692d6eeSBarry Smith #define MATORDERINGND          "nd"
9242692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
9252692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
9262692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
9272692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
928510184a7SJed Brown #define MATORDERINGWBM         "wbm"
929fc1bef75SJed Brown #define MATORDERINGSPECTRAL    "spectral"
930312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
931b12f92e5SBarry Smith 
93219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
933140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
934bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
935607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void);
936014dd563SJed Brown PETSC_EXTERN PetscBool         MatOrderingRegisterAllCalled;
937140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList;
938d4fbbf0eSBarry Smith 
939014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
940fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);
941a2ce50c7SBarry Smith 
942d91e6319SBarry Smith /*S
943d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
944d90ac83dSHong Zhang 
945d90ac83dSHong Zhang    Level: beginner
946d90ac83dSHong Zhang 
947d90ac83dSHong Zhang S*/
948d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
9496a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[];
9505e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];
951d90ac83dSHong Zhang 
952d90ac83dSHong Zhang /*S
9532401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
954b00f7748SHong Zhang 
95561cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
95661cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
957b00f7748SHong Zhang 
95815e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
959b00f7748SHong Zhang 
96061cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
96161cad860SBarry Smith 
962b00f7748SHong Zhang    Level: developer
963b00f7748SHong Zhang 
964d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
965d7d82daaSBarry Smith           MatFactorInfoInitialize()
966b00f7748SHong Zhang 
967b00f7748SHong Zhang S*/
968b00f7748SHong Zhang typedef struct {
96915e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
97085317021SBarry Smith   PetscReal     usedt;
97115e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
972b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
97315e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
97467e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
975348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
976bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
977bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
97815e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
979f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
980f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
98115e8a5b3SHong Zhang } MatFactorInfo;
982ffa6d0a5SLois Curfman McInnes 
983014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
984014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
985014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
986014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
987014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
988014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
989014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
990014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
991014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
992014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
993014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
994014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
995014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
996014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
997014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
998014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
999014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1000014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1001014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
10028ed539a5SBarry Smith 
1003014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
1004bb5a7306SBarry Smith 
1005d91e6319SBarry Smith /*E
1006d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1007bb1eb677SSatish Balay 
1008d91e6319SBarry Smith     Level: beginner
1009d91e6319SBarry Smith 
1010d9274352SBarry Smith    May be bitwise ORd together
1011d9274352SBarry Smith 
1012d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1013d91e6319SBarry Smith 
10144e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10154e7234bfSBarry Smith 
101641f059aeSBarry Smith .seealso: MatSOR()
1017d91e6319SBarry Smith E*/
1018ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1019ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1020ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
102184cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1022014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10238ed539a5SBarry Smith 
1024d4fbbf0eSBarry Smith /*
1025639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1026639f9d9dSBarry Smith */
1027b12f92e5SBarry Smith 
10287086a01eSPeter Brune /*S
10297086a01eSPeter Brune      MatColoring - Object for managing the coloring of matrices.
10307086a01eSPeter Brune 
10317086a01eSPeter Brune    Level: beginner
10327086a01eSPeter Brune 
10337086a01eSPeter Brune   Concepts: matrix, coloring
10347086a01eSPeter Brune 
10357086a01eSPeter Brune .seealso:  MatFDColoringCreate() ISColoring MatFDColoring
10367086a01eSPeter Brune S*/
10377086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring;
103876bdecfbSBarry Smith /*J
10398f6c3df8SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring
1040d9274352SBarry Smith 
1041d9274352SBarry Smith    Level: beginner
1042d9274352SBarry Smith 
10437086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring
104476bdecfbSBarry Smith J*/
10457086a01eSPeter Brune 
104619fd82e9SBarry Smith typedef const  char*           MatColoringType;
10475567d71aSPeter Brune #define MATCOLORINGJP      "jp"
10482aa85f04SPeter Brune #define MATCOLORINGMIS     "mis"
10492692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
10502692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
10512692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
10522692d6eeSBarry Smith #define MATCOLORINGID      "id"
1053b12f92e5SBarry Smith 
1054335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1055335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1056335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1057335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1058335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1059335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1060335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1061335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1062335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1063335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1064607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void);
1065335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1066335efc43SPeter Brune PETSC_EXTERN PetscBool MatColoringRegisterAllCalled;
1067014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1068639f9d9dSBarry Smith 
1069d9274352SBarry Smith /*S
1070d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1071d9274352SBarry Smith         and coloring
1072639f9d9dSBarry Smith 
1073d9274352SBarry Smith    Level: beginner
1074d9274352SBarry Smith 
1075d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1076d9274352SBarry Smith 
1077d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1078d9274352SBarry Smith S*/
1079e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1080639f9d9dSBarry Smith 
1081014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1082014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1083014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1084014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1085014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1086014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1087014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1088014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1089014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1090014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1091f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1092f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1093f86b9fbaSHong Zhang 
1094b1683b59SHong Zhang 
1095b1683b59SHong Zhang /*S
1096b9af6bddSHong Zhang      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1097b1683b59SHong Zhang 
1098b1683b59SHong Zhang    Level: beginner
1099b1683b59SHong Zhang 
1100b1683b59SHong Zhang   Concepts: coloring, sparse matrix product
1101b1683b59SHong Zhang 
1102b9af6bddSHong Zhang .seealso:  MatTransposeColoringCreate()
1103b1683b59SHong Zhang S*/
1104b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1105b1683b59SHong Zhang 
1106014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1107014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1108014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1109014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1110b1683b59SHong Zhang 
1111639f9d9dSBarry Smith /*
11120752156aSBarry Smith     These routines are for partitioning matrices: currently used only
11133eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
11140752156aSBarry Smith */
1115ca161407SBarry Smith 
1116d9274352SBarry Smith /*S
1117d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1118d9274352SBarry Smith 
1119d9274352SBarry Smith    Level: beginner
1120d9274352SBarry Smith 
1121d9274352SBarry Smith   Concepts: partitioning
1122d9274352SBarry Smith 
1123743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1124d9274352SBarry Smith S*/
112591e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1126d9274352SBarry Smith 
112776bdecfbSBarry Smith /*J
11288f6c3df8SBarry Smith     MatPartitioningType - String with the name of a PETSc matrix partitioning
1129d9274352SBarry Smith 
1130d9274352SBarry Smith    Level: beginner
11310d04baf8SBarry Smith dm
1132b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
113376bdecfbSBarry Smith J*/
113419fd82e9SBarry Smith typedef const char* MatPartitioningType;
11352692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
11362692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
11372692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
11382692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
11392692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
114050d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
114171306c60Sjroman 
1142ca161407SBarry Smith 
1143014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
114419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1145014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1146014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1147014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1148014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1149014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1150014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
11512aabb6bbSBarry Smith 
1152bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
115330de9b25SBarry Smith 
1154014dd563SJed Brown PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled;
1155f1144a30SSatish Balay 
1156607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void);
11572bad1931SBarry Smith 
1158014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1159014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
116019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1161ca161407SBarry Smith 
1162014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1163014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
11640752156aSBarry Smith 
1165b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
11666a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1167b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
11686a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[];
1169b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
11706a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[];
1171b6956c12SJose Roman 
1172014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1173014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1174014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1178014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1179014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1180014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1181014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1182014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
118371306c60Sjroman 
118471306c60Sjroman #define MP_PARTY_OPT "opt"
118571306c60Sjroman #define MP_PARTY_LIN "lin"
118671306c60Sjroman #define MP_PARTY_SCA "sca"
118771306c60Sjroman #define MP_PARTY_RAN "ran"
118871306c60Sjroman #define MP_PARTY_GBF "gbf"
118971306c60Sjroman #define MP_PARTY_GCF "gcf"
119071306c60Sjroman #define MP_PARTY_BUB "bub"
119171306c60Sjroman #define MP_PARTY_DEF "def"
1192014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
119371306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
119471306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
119571306c60Sjroman #define MP_PARTY_NONE "no"
1196014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1197014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1198014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1199014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
120071306c60Sjroman 
120150d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
12026a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1203e0f1cffaSJose Roman 
1204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1206014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
120871306c60Sjroman 
1209b43b03e9SMark F. Adams /*
1210b43b03e9SMark F. Adams     These routines are for coarsening matrices:
1211b43b03e9SMark F. Adams */
1212b43b03e9SMark F. Adams 
1213b43b03e9SMark F. Adams /*S
1214b43b03e9SMark F. Adams      MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix)
1215b43b03e9SMark F. Adams 
1216b43b03e9SMark F. Adams    Level: beginner
1217b43b03e9SMark F. Adams 
1218b43b03e9SMark F. Adams   Concepts: coarsen
1219b43b03e9SMark F. Adams 
1220b43b03e9SMark F. Adams .seealso:  MatCoarsenCreate), MatCoarsenType
1221b43b03e9SMark F. Adams S*/
1222b43b03e9SMark F. Adams typedef struct _p_MatCoarsen* MatCoarsen;
1223b43b03e9SMark F. Adams 
1224b43b03e9SMark F. Adams /*J
12258f6c3df8SBarry Smith     MatCoarsenType - String with the name of a PETSc matrix coarsen
1226b43b03e9SMark F. Adams 
1227b43b03e9SMark F. Adams    Level: beginner
12288f6c3df8SBarry Smith 
1229b43b03e9SMark F. Adams .seealso: MatCoarsenCreate(), MatCoarsen
1230b43b03e9SMark F. Adams J*/
123119fd82e9SBarry Smith typedef const char* MatCoarsenType;
1232b43b03e9SMark F. Adams #define MATCOARSENMIS  "mis"
12339057884aSMark F. Adams #define MATCOARSENHEM  "hem"
1234b43b03e9SMark F. Adams 
12350cbbd2e1SMark F. Adams /* linked list for aggregates */
123641b27cdeSMark F. Adams typedef struct _PetscCDIntNd{
123741b27cdeSMark F. Adams   struct _PetscCDIntNd *next;
12380cbbd2e1SMark F. Adams   PetscInt             gid;
123941b27cdeSMark F. Adams }PetscCDIntNd;
12400cbbd2e1SMark F. Adams 
12410cbbd2e1SMark F. Adams /* only used by node pool */
124241b27cdeSMark F. Adams typedef struct _PetscCDArrNd{
124341b27cdeSMark F. Adams   struct _PetscCDArrNd *next;
124441b27cdeSMark F. Adams   struct _PetscCDIntNd *array;
124541b27cdeSMark F. Adams }PetscCDArrNd;
12460cbbd2e1SMark F. Adams 
12470cbbd2e1SMark F. Adams typedef struct _PetscCoarsenData{
124882c86c8fSBarry Smith   PetscCDArrNd pool_list;  /* node pool */
124941b27cdeSMark F. Adams   PetscCDIntNd *new_node;
12500cbbd2e1SMark F. Adams   PetscInt     new_left;
12510cbbd2e1SMark F. Adams   PetscInt     chk_sz;
125241b27cdeSMark F. Adams   PetscCDIntNd *extra_nodes;
125382c86c8fSBarry Smith   PetscCDIntNd **array;  /* Array of lists */
12540cbbd2e1SMark F. Adams   PetscInt     size;
125582c86c8fSBarry Smith   Mat          mat;  /* cache a Mat for communication data */
12560cbbd2e1SMark F. Adams }PetscCoarsenData;
12570cbbd2e1SMark F. Adams 
1258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*);
125919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,MatCoarsenType);
1260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat);
1261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS);
1262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool);
1263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetVerbose(MatCoarsen,PetscInt);
1264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** );
1265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen);
1266014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*);
1267b43b03e9SMark F. Adams 
1268bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],PetscErrorCode (*)(MatCoarsen));
1269b43b03e9SMark F. Adams 
1270014dd563SJed Brown PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled;
1271b43b03e9SMark F. Adams 
1272607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void);
1273b43b03e9SMark F. Adams 
1274014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer);
1275014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen);
127619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,MatCoarsenType*);
1277ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatCoarsenViewFromOptions(MatCoarsen A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
1278b43b03e9SMark F. Adams 
1279014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1280014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1281591294e4SBarry Smith 
12820752156aSBarry Smith /*
12830a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1284d4fbbf0eSBarry Smith */
12851c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
12861c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
12871c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
12881c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
12891c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
12907c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
12917c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
12921c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
12931c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
12947c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
12957c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
12961c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
12971c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1298a714c06dSBarry Smith                MATOP_SOR=13,
12991c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13001c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13011c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13021c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13031c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13041c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13051c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13061c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1307d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1308d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1309d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1310d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1311d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1312d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1313d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1314d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1315d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1316d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
13179d39f69dSJed Brown                /* MATOP_PLACEHOLDER_32=32, */
13189d39f69dSJed Brown                /* MATOP_PLACEHOLDER_33=33, */
1319d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1320d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1321d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1322d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1323d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1324d519adbfSMatthew Knepley                MATOP_AXPY=39,
1325d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1326d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1327d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1328d519adbfSMatthew Knepley                MATOP_COPY=43,
1329d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1330d519adbfSMatthew Knepley                MATOP_SCALE=45,
1331d519adbfSMatthew Knepley                MATOP_SHIFT=46,
133235153367SBarry Smith                MATOP_DIAGONAL_SET=47,
13339d39f69dSJed Brown                MATOP_ZERO_ROWS_COLUMNS=48,
13349d39f69dSJed Brown                MATOP_SET_RANDOM=49,
1335d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1336d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1337d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1338d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1339d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1340d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1341d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1342d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1343d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1344d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1345d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1346d519adbfSMatthew Knepley                MATOP_VIEW=61,
1347d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
13487bab7c10SHong Zhang                MATOP_MATMAT_MULT=63,
13497bab7c10SHong Zhang                MATOP_MATMAT_MULT_SYMBOLIC=64,
13507bab7c10SHong Zhang                MATOP_MATMAT_MULT_NUMERIC=65,
1351d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1352d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1353d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1354d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1355d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1356d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1357d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
13589d39f69dSJed Brown                /* MATOP_PLACEHOLDER_73=73, */
1359d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1360d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1361d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1362bda6c4cbSBarry Smith                MATOP_MULT_CONSTRAINED=77,
1363bda6c4cbSBarry Smith                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
13649d39f69dSJed Brown                MATOP_FIND_ZERO_DIAGONALS=79,
1365d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1366d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1367d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1368d519adbfSMatthew Knepley                MATOP_LOAD=83,
1369d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1370d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1371d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1372820469bcSHong Zhang                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1373d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1374d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1375d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1376d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1377d519adbfSMatthew Knepley                MATOP_PTAP=92,
1378d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1379d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1380bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT=95,
1381bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1382bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
13839d39f69dSJed Brown                /* MATOP_PLACEHOLDER_98=98, */
13849d39f69dSJed Brown                /* MATOP_PLACEHOLDER_99=99, */
13859d39f69dSJed Brown                /* MATOP_PLACEHOLDER_100=100, */
13869d39f69dSJed Brown                /* MATOP_PLACEHOLDER_101=101, */
1387d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
13889d39f69dSJed Brown                /* MATOP_PLACEHOLDER_103=103, */
1389d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1390d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1391bda6c4cbSBarry Smith                MATOP_IMAGINARY_PART=106,
1392bda6c4cbSBarry Smith                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1393bda6c4cbSBarry Smith                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1394bda6c4cbSBarry Smith                MATOP_MAT_SOLVE=109,
13950e8d04f7SBarry Smith                MATOP_GET_REDUNDANT_MATRIX=110,
1396d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
13970e8d04f7SBarry Smith                MATOP_GET_COLUMN_VECTOR=112,
1398d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
13990e8d04f7SBarry Smith                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
140089c6957cSBarry Smith                MATOP_CREATE=115,
140189c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
14020e8d04f7SBarry Smith                MATOP_GET_LOCAL_SUB_MATRIX=117,
14030e8d04f7SBarry Smith                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1404eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
14050e8d04f7SBarry Smith                MATOP_HERMITIAN_TRANSPOSE=120,
14060e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
14070e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
14080e8d04f7SBarry Smith                MATOP_GET_MULTI_PROC_BLOCK=123,
14099d39f69dSJed Brown                MATOP_FIND_NONZERO_ROWS=124,
14100e8d04f7SBarry Smith                MATOP_GET_COLUMN_NORMS=125,
14119d39f69dSJed Brown                MATOP_INVERT_BLOCK_DIAGONAL=126,
14129d39f69dSJed Brown                /* MATOP_PLACEHOLDER_127=127, */
14130e8d04f7SBarry Smith                MATOP_GET_SUB_MATRICES_PARALLE=128,
14142b8ad9a3SHong Zhang                MATOP_SET_VALUES_BATCH=129,
14150e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT=130,
14160e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
14170e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
14180e8d04f7SBarry Smith                MATOP_TRANSPOSE_COLORING_CREAT=133,
14190e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_SPT=134,
14200e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_DEN=135,
14210e8d04f7SBarry Smith                MATOP_RART=136,
14220e8d04f7SBarry Smith                MATOP_RART_SYMBOLIC=137,
14230e8d04f7SBarry Smith                MATOP_RART_NUMERIC=138,
1424e09a3074SHong Zhang                MATOP_SET_BLOCK_SIZES=139,
1425f9426fe0SMark Adams                MATOP_AYPX=140,
14261919a2e2SJed Brown                MATOP_RESIDUAL=141,
14271919a2e2SJed Brown                MATOP_FDCOLORING_SETUP= 142
1428fae171e0SBarry Smith              } MatOperation;
1429014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *);
1430014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1431014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1432014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1433112a2221SBarry Smith 
143490ace30eSBarry Smith /*
143590ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
143690ace30eSBarry Smith    stored in a universal format. By changing the format with
14377973ad04SBarry Smith    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
143890ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
143990ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1440f4403165SShri Abhyankar    read into matrices of the same type.
144190ace30eSBarry Smith */
144290ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
144390ace30eSBarry Smith 
1444014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1445014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1446014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
14471f4e1ec7SBarry Smith 
1448d9274352SBarry Smith /*S
1449d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1450d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1451d9274352SBarry Smith 
1452f7a9e4ceSBarry Smith    Level: advanced
1453d9274352SBarry Smith 
1454d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1455d9274352SBarry Smith 
14566e1639daSBarry Smith   Users manual sections:
14576e1639daSBarry Smith .   sec_singular
14586e1639daSBarry Smith 
1459d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1460d9274352SBarry Smith S*/
146174637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1462d9274352SBarry Smith 
1463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1466d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
147574637425SBarry Smith 
1476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
14803f1d51d7SBarry Smith 
1481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1484c4f061fbSSatish Balay 
1485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*);
1486b0a32e0cSBarry Smith 
1487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
148804f1ad80SBarry Smith 
1489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace);
1495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1503e884886fSBarry Smith 
15046370053bSBarry Smith /*S
15056370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
15066370053bSBarry Smith               Jacobian vector products
1507e884886fSBarry Smith 
15086370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
15096370053bSBarry Smith 
15106370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
15116370053bSBarry Smith 
15126370053bSBarry Smith     Level: developer
15136370053bSBarry Smith 
15146370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
15156370053bSBarry Smith S*/
1516e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1517e884886fSBarry Smith 
151876bdecfbSBarry Smith /*J
1519e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1520e884886fSBarry Smith 
1521e884886fSBarry Smith    Level: beginner
1522e884886fSBarry Smith 
1523e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
152476bdecfbSBarry Smith J*/
152519fd82e9SBarry Smith typedef const char* MatMFFDType;
1526e884886fSBarry Smith #define MATMFFD_DS  "ds"
1527e884886fSBarry Smith #define MATMFFD_WP  "wp"
1528e884886fSBarry Smith 
152919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1530bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));
1531e884886fSBarry Smith 
1532607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegisterAll(void);
1533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool );
1535e884886fSBarry Smith 
1536014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1537014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
15387dbadf16SMatthew Knepley 
153997969023SHong Zhang /*
154097969023SHong Zhang    PETSc interface to MUMPS
154197969023SHong Zhang */
154297969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1544b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
154597969023SHong Zhang #endif
154623a5497aSJed Brown 
1547d954db57SHong Zhang /*
1548d954db57SHong Zhang    PETSc interface to SUPERLU
1549d954db57SHong Zhang */
1550d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1552d954db57SHong Zhang #endif
1553d954db57SHong Zhang 
1554aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA
15553f9c0db1SPaul Mullowney /*E
1556e057df02SPaul Mullowney     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
15572692e278SPaul Mullowney     matrices.
1558e057df02SPaul Mullowney 
1559e057df02SPaul Mullowney     Not Collective
1560e057df02SPaul Mullowney 
15618468deeeSKarl Rupp +   MAT_CUSPARSE_CSR - Compressed Sparse Row
15622692e278SPaul Mullowney .   MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later).
15632692e278SPaul Mullowney -   MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).
1564e057df02SPaul Mullowney 
1565e057df02SPaul Mullowney     Level: intermediate
1566e057df02SPaul Mullowney 
1567e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1568e057df02SPaul Mullowney 
1569e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1570e057df02SPaul Mullowney E*/
1571e057df02SPaul Mullowney 
15723f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;
1573e057df02SPaul Mullowney 
1574e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1575e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
1576e057df02SPaul Mullowney 
15773f9c0db1SPaul Mullowney /*E
1578e057df02SPaul Mullowney     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
15792692e278SPaul Mullowney     matrices whose operation should use a particular storage format.
1580e057df02SPaul Mullowney 
1581e057df02SPaul Mullowney     Not Collective
1582e057df02SPaul Mullowney 
15838468deeeSKarl Rupp +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
15848468deeeSKarl Rupp .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
15858468deeeSKarl Rupp .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
15868468deeeSKarl Rupp -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices
1587e057df02SPaul Mullowney 
1588e057df02SPaul Mullowney     Level: intermediate
1589e057df02SPaul Mullowney 
1590e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1591e057df02SPaul Mullowney E*/
159236d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;
1593e057df02SPaul Mullowney 
159410e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
159510e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1596e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
15979ae82921SPaul Mullowney #endif
15989ae82921SPaul Mullowney 
159990c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
1600014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1601014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1602e057df02SPaul Mullowney 
16033f9c0db1SPaul Mullowney /*E
1604e057df02SPaul Mullowney     MatCUSPStorageFormat - indicates the storage format for CUSP (GPU)
160536d62e41SPaul Mullowney     matrices.
1606e057df02SPaul Mullowney 
1607e057df02SPaul Mullowney     Not Collective
1608e057df02SPaul Mullowney 
16098468deeeSKarl Rupp +   MAT_CUSP_CSR - Compressed Sparse Row
16108468deeeSKarl Rupp .   MAT_CUSP_DIA - Diagonal
16118468deeeSKarl Rupp -   MAT_CUSP_ELL - Ellpack
1612e057df02SPaul Mullowney 
1613e057df02SPaul Mullowney     Level: intermediate
1614e057df02SPaul Mullowney 
1615e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1616e057df02SPaul Mullowney 
1617e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
1618e057df02SPaul Mullowney E*/
16193f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat;
1620e057df02SPaul Mullowney 
1621e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1622e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[];
1623e057df02SPaul Mullowney 
16243f9c0db1SPaul Mullowney /*E
1625e057df02SPaul Mullowney     MatCUSPFormatOperation - indicates the operation of CUSP (GPU)
162636d62e41SPaul Mullowney     matrices whose operation should use a particular storage format.
1627e057df02SPaul Mullowney 
1628e057df02SPaul Mullowney     Not Collective
1629e057df02SPaul Mullowney 
16308468deeeSKarl Rupp +   MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
16318468deeeSKarl Rupp .   MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
16328468deeeSKarl Rupp .   MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
16338468deeeSKarl Rupp -   MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices
1634e057df02SPaul Mullowney 
1635e057df02SPaul Mullowney     Level: intermediate
1636e057df02SPaul Mullowney 
1637e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1638e057df02SPaul Mullowney 
1639e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat
1640e057df02SPaul Mullowney E*/
164136d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation;
1642e057df02SPaul Mullowney 
1643e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat);
1644e057df02SPaul Mullowney #endif
164590c53902SBarry Smith 
1646d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL)
1647d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1648d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1649d67ff14aSKarl Rupp #endif
1650d67ff14aSKarl Rupp 
165154efbe56SHong Zhang /*
165254efbe56SHong Zhang    PETSc interface to FFTW
165354efbe56SHong Zhang */
165454efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
1655014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1656014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
1657014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*);
165854efbe56SHong Zhang #endif
165954efbe56SHong Zhang 
1660db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL)
1661607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
1662db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
1663db31f6deSJed Brown #endif
1664db31f6deSJed Brown 
1665014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1666014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1667014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1668014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1669014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1670014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
167119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1672014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1673014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1674d8588912SDave May 
16754325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
1676e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);
16774325cce7SMatthew G Knepley 
167823a5497aSJed Brown #endif
1679