xref: /petsc/include/petscmat.h (revision d0195637390c19153758d8270d49d9757db8dc56)
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
9d9274352SBarry Smith      Mat - Abstract PETSc matrix object
102eac72dbSBarry Smith 
11d91e6319SBarry Smith    Level: beginner
12d91e6319SBarry Smith 
13d9274352SBarry Smith   Concepts: matrix; linear operator
14d9274352SBarry Smith 
15d9274352SBarry Smith .seealso:  MatCreate(), MatType, MatSetType()
16d9274352SBarry Smith S*/
17d9274352SBarry Smith typedef struct _p_Mat*           Mat;
18d9274352SBarry Smith 
1976bdecfbSBarry Smith /*J
20d9274352SBarry Smith     MatType - String with the name of a PETSc matrix or the creation function
21d9274352SBarry Smith        with an optional dynamic library name, for example
22d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mymatcreate()
23d9274352SBarry Smith 
24d9274352SBarry Smith    Level: beginner
25d9274352SBarry Smith 
26c7393fdbSBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage
2776bdecfbSBarry Smith J*/
2819fd82e9SBarry Smith typedef const char* MatType;
29273d9f13SBarry Smith #define MATSAME            "same"
305a11e1b2SBarry Smith #define MATMAIJ            "maij"
31273d9f13SBarry Smith #define MATSEQMAIJ         "seqmaij"
32273d9f13SBarry Smith #define MATMPIMAIJ         "mpimaij"
33273d9f13SBarry Smith #define MATIS              "is"
345a11e1b2SBarry Smith #define MATAIJ             "aij"
35273d9f13SBarry Smith #define MATSEQAIJ          "seqaij"
367d6a0e61SBarry Smith #define MATSEQAIJPTHREAD   "seqaijpthread"
37bf2c1783SBarry Smith #define MATAIJPTHREAD      "aijpthread"
38273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
395a11e1b2SBarry Smith #define MATAIJCRL          "aijcrl"
405a11e1b2SBarry Smith #define MATSEQAIJCRL       "seqaijcrl"
415a11e1b2SBarry Smith #define MATMPIAIJCRL       "mpiaijcrl"
428154be41SBarry Smith #define MATAIJCUSP         "aijcusp"
438154be41SBarry Smith #define MATSEQAIJCUSP      "seqaijcusp"
448154be41SBarry Smith #define MATMPIAIJCUSP      "mpiaijcusp"
459ae82921SPaul Mullowney #define MATAIJCUSPARSE     "aijcusparse"
469ae82921SPaul Mullowney #define MATSEQAIJCUSPARSE  "seqaijcusparse"
479ae82921SPaul Mullowney #define MATMPIAIJCUSPARSE  "mpiaijcusparse"
485a11e1b2SBarry Smith #define MATAIJPERM         "aijperm"
495a11e1b2SBarry Smith #define MATSEQAIJPERM      "seqaijperm"
505a11e1b2SBarry Smith #define MATMPIAIJPERM      "mpiaijperm"
51273d9f13SBarry Smith #define MATSHELL           "shell"
525a11e1b2SBarry Smith #define MATDENSE           "dense"
53209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
54273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
55db31f6deSJed Brown #define MATELEMENTAL       "elemental"
565a11e1b2SBarry Smith #define MATBAIJ            "baij"
57273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
58273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
59273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
605a11e1b2SBarry Smith #define MATSBAIJ           "sbaij"
61273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
62273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
63c0cdd4a1SDahai Guo #define MATSEQBSTRM        "seqbstrm"
64c0cdd4a1SDahai Guo #define MATMPIBSTRM        "mpibstrm"
65c0cdd4a1SDahai Guo #define MATBSTRM           "bstrm"
66aa5a9175SDahai Guo #define MATSEQSBSTRM       "seqsbstrm"
67aa5a9175SDahai Guo #define MATMPISBSTRM       "mpisbstrm"
68aa5a9175SDahai Guo #define MATSBSTRM          "sbstrm"
69cebc7f6cSBarry Smith #define MATDAAD            "daad"
70cebc7f6cSBarry Smith #define MATMFFD            "mffd"
71c8a8475eSBarry Smith #define MATNORMAL          "normal"
72ab92ecdeSBarry Smith #define MATLRC             "lrc"
732a6744ebSBarry Smith #define MATSCATTER         "scatter"
74421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
75793850ffSBarry Smith #define MATCOMPOSITE       "composite"
761f8c7532SHong Zhang #define MATFFT             "fft"
771f8c7532SHong Zhang #define MATFFTW            "fftw"
78e133240eSMatthew G Knepley #define MATSEQCUFFT        "seqcufft"
79557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
8072ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
811d6018f0SLisandro Dalcin #define MATPYTHON          "python"
82f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
83a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT    "hypresstruct"
84ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
852c0dbf93SJed Brown #define MATLOCALREF        "localref"
86d8588912SDave May #define MATNEST            "nest"
87773941d6SBarry Smith 
8876bdecfbSBarry Smith /*J
89c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
909989ab13SBarry Smith 
919989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
929989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
939989ab13SBarry Smith 
949989ab13SBarry Smith 
959989ab13SBarry Smith    Level: beginner
969989ab13SBarry Smith 
975c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
9876bdecfbSBarry Smith J*/
99c7393fdbSBarry Smith #define MatSolverPackage char*
1002692d6eeSBarry Smith #define MATSOLVERSUPERLU      "superlu"
1012692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist"
1022692d6eeSBarry Smith #define MATSOLVERUMFPACK      "umfpack"
10320db9a53SJed Brown #define MATSOLVERCHOLMOD      "cholmod"
1042692d6eeSBarry Smith #define MATSOLVERESSL         "essl"
1052692d6eeSBarry Smith #define MATSOLVERLUSOL        "lusol"
1062692d6eeSBarry Smith #define MATSOLVERMUMPS        "mumps"
1072692d6eeSBarry Smith #define MATSOLVERPASTIX       "pastix"
1082692d6eeSBarry Smith #define MATSOLVERMATLAB       "matlab"
1092692d6eeSBarry Smith #define MATSOLVERPETSC        "petsc"
1102692d6eeSBarry Smith #define MATSOLVERBAS          "bas"
1119ae82921SPaul Mullowney #define MATSOLVERCUSPARSE     "cusparse"
112aa5a9175SDahai Guo #define MATSOLVERBSTRM        "bstrm"
113aa5a9175SDahai Guo #define MATSOLVERSBSTRM       "sbstrm"
11415767789SHong Zhang #define MATSOLVERELEMENTAL    "elemental"
11517f1a0eaSHong Zhang #define MATSOLVERCLIQUE       "clique"
116c0cdd4a1SDahai Guo 
117b24902e0SBarry Smith /*E
118b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
119b24902e0SBarry Smith 
120b24902e0SBarry Smith     Level: beginner
121b24902e0SBarry Smith 
122b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
123b24902e0SBarry Smith 
124c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
125b24902e0SBarry Smith E*/
126599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
127014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[];
128e92e720dSBarry Smith 
129014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
130014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
131014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
132014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*);
1339989ab13SBarry Smith 
134c06d978dSMatthew Knepley /* Logging support */
1350700a824SBarry Smith #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
136014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID;
137014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID;
138014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
139014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID;
140014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID;
141014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID;
142014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID;
143c06d978dSMatthew Knepley 
144ceb03754SKris Buschelman /*E
145ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
146d6eff37eSBarry Smith      or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate
147d6eff37eSBarry Smith      that the input matrix is to be replaced with the converted matrix.
148ceb03754SKris Buschelman 
149ceb03754SKris Buschelman     Level: beginner
150ceb03754SKris Buschelman 
151ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
152ceb03754SKris Buschelman 
1530c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
154ceb03754SKris Buschelman E*/
155dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse;
156ceb03754SKris Buschelman 
1575494a064SHong Zhang /*E
1585494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
159829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1605494a064SHong Zhang 
1615494a064SHong Zhang     Level: beginner
1625494a064SHong Zhang 
163829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1645494a064SHong Zhang E*/
1655494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1665494a064SHong Zhang 
167607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void);
168c06d978dSMatthew Knepley 
169014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*);
170014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
17119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType);
172014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat);
173146574abSBarry Smith PETSC_EXTERN PetscErrorCode MatViewFromOptions(Mat,const char[],const char[]);
174607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatRegisterAll(void);
175bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat));
176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]);
177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]);
178014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]);
179014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]);
180f69a0ea3SMatthew Knepley 
181014dd563SJed Brown PETSC_EXTERN PetscBool         MatRegisterAllCalled;
182140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList;
183140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList;
184140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList;
185140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatCoarsenList;
18628988994SBarry Smith 
1873b224e63SBarry Smith /*E
1883b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
1893b224e63SBarry Smith 
1903b224e63SBarry Smith     Level: beginner
1913b224e63SBarry Smith 
1923b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1933b224e63SBarry Smith 
1943b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
1953b224e63SBarry Smith E*/
196214bc6a2SJed Brown typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN,SAME_PRECONDITIONER} MatStructure;
1973b224e63SBarry Smith 
198014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
199014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
200014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
201014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
202014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2048d7a6e47SBarry Smith 
205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
206014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
208d21a29f3SJed Brown 
209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
211d21a29f3SJed Brown 
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*);
216dfb205c3SBarry Smith 
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
218014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
219014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*);
220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
223c0cdd4a1SDahai Guo 
224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
228c0cdd4a1SDahai Guo 
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2366d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
2386d7c1e57SBarry Smith 
23919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*);
240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
241dedccee8SHong Zhang 
242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*);
244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS);
245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
2461472f72bSBarry Smith 
247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]);
2481d6018f0SLisandro Dalcin 
249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
25121c89e3eSBarry Smith 
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
257713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
25899cafbc1SBarry Smith 
2598ed539a5SBarry Smith /* ------------------------------------------------------------*/
260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
26573a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);
26684cb2905SBarry Smith 
2672ef4de8bSBarry Smith /*S
2682ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
269be479b30SJed Brown         column of a matrix as indexed on an associated grid.
270be479b30SJed Brown 
271be479b30SJed Brown    Fortran usage is different, see MatSetValuesStencil() for details.
2722ef4de8bSBarry Smith 
2732ef4de8bSBarry Smith    Level: beginner
2742ef4de8bSBarry Smith 
2752ef4de8bSBarry Smith   Concepts: matrix; linear operator
2762ef4de8bSBarry Smith 
277be479b30SJed Brown .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil()
2782ef4de8bSBarry Smith S*/
279435da068SBarry Smith typedef struct {
280c1ac3661SBarry Smith   PetscInt k,j,i,c;
281435da068SBarry Smith } MatStencil;
2822ef4de8bSBarry Smith 
283014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
284014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
285014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
286435da068SBarry Smith 
287014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetColoring(Mat,ISColoring);
288014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*);
2893a7fca6bSBarry Smith 
290d91e6319SBarry Smith /*E
291d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
292d91e6319SBarry Smith      to continue to add values to it
293d91e6319SBarry Smith 
294d91e6319SBarry Smith     Level: beginner
295d91e6319SBarry Smith 
296d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
297d91e6319SBarry Smith E*/
2986d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
299014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
300014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
301014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);
3024f9c727eSBarry Smith 
3031d73ed98SBarry Smith 
30430de9b25SBarry Smith 
305d91e6319SBarry Smith /*E
306d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
307d91e6319SBarry Smith 
308d91e6319SBarry Smith     Level: beginner
309d91e6319SBarry Smith 
3100a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
311d91e6319SBarry Smith 
312382164b0SBarry Smith    Developer Notes: Entries that are negative need not be called collectively by all processes.
313382164b0SBarry Smith 
314d91e6319SBarry Smith .seealso: MatSetOption()
315d91e6319SBarry Smith E*/
316382164b0SBarry Smith typedef enum {MAT_OPTION_MIN = -8,
317382164b0SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR = -7,
318382164b0SBarry Smith               MAT_NO_OFF_PROC_ZERO_ROWS = -6,
319382164b0SBarry Smith               MAT_NO_OFF_PROC_ENTRIES = -5,
320382164b0SBarry Smith               MAT_UNUSED_NONZERO_LOCATION_ERR = -4,
321382164b0SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR = -3,
322382164b0SBarry Smith               MAT_ROW_ORIENTED = -2,
323382164b0SBarry Smith               MAT_NEW_NONZERO_LOCATIONS = -1,
324382164b0SBarry Smith               MAT_SYMMETRIC = 1,
325382164b0SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC = 2,
326382164b0SBarry Smith               MAT_NEW_DIAGONALS = 3,
327382164b0SBarry Smith               MAT_IGNORE_OFF_PROC_ENTRIES = 4,
328382164b0SBarry Smith               MAT_USE_HASH_TABLE = 5,
329382164b0SBarry Smith               MAT_KEEP_NONZERO_PATTERN = 6,
330382164b0SBarry Smith               MAT_IGNORE_ZERO_ENTRIES = 7,
331382164b0SBarry Smith               MAT_USE_INODES = 8,
332382164b0SBarry Smith               MAT_HERMITIAN = 9,
333382164b0SBarry Smith               MAT_SYMMETRY_ETERNAL = 10,
334382164b0SBarry Smith               MAT_CHECK_COMPRESSED_ROW = 11,
335382164b0SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR = 12,
336382164b0SBarry Smith               MAT_ERROR_LOWER_TRIANGULAR = 13,
337382164b0SBarry Smith               MAT_GETROW_UPPERTRIANGULAR = 14,
338382164b0SBarry Smith               MAT_SPD = 15,
339382164b0SBarry Smith               MAT_OPTION_MAX = 16} MatOption;
340e057df02SPaul Mullowney 
341014dd563SJed Brown PETSC_EXTERN const char *MatOptions[];
342014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool );
34319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);
34484cb2905SBarry Smith 
345014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
346014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
347014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
348014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
349014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
350014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
351014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
3538c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
3548c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
3558c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
3568c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNThreads(Mat,PetscInt);
362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNThreads(Mat,PetscInt*);
3631620fd73SBarry Smith 
364014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
365014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
366014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
367014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
368014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
369014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
376f5edf698SHong Zhang 
377d91e6319SBarry Smith /*E
378d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
379d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
380d91e6319SBarry Smith 
381d91e6319SBarry Smith     Level: beginner
382d91e6319SBarry Smith 
383d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
384d91e6319SBarry Smith 
38570dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
38670dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
38770dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
38870dcbbb9SBarry Smith 
389d91e6319SBarry Smith .seealso: MatDuplicate()
390d91e6319SBarry Smith E*/
39170dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
3922e8a6d31SBarry Smith 
39319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
394014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);
39594a9d846SBarry Smith 
396cb5b572fSBarry Smith 
397014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
398014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
399014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
400014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
401014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
402014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
403014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
404014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
405014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
4067b80b807SBarry Smith 
4071a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4081a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
4091a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4101a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
411d4fbbf0eSBarry Smith 
412d91e6319SBarry Smith /*S
413d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
414d91e6319SBarry Smith 
415d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
416d91e6319SBarry Smith 
417d91e6319SBarry Smith    Level: intermediate
418d91e6319SBarry Smith 
419d91e6319SBarry Smith   Concepts: matrix^nonzero information
420d91e6319SBarry Smith 
421d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
422d91e6319SBarry Smith S*/
4234e220ebcSLois Curfman McInnes typedef struct {
424b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
425b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
426b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
427b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
428b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
429b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
430b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4314e220ebcSLois Curfman McInnes } MatInfo;
4324e220ebcSLois Curfman McInnes 
433d9274352SBarry Smith /*E
434d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
435d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
436d9274352SBarry Smith 
437d9274352SBarry Smith     Level: beginner
438d9274352SBarry Smith 
439d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
440d9274352SBarry Smith 
441d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
442d9274352SBarry Smith E*/
4437b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
444014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
445014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
446014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
447014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
448014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
449014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
450014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
452014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat *);
454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);
456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool *);
457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
4617b80b807SBarry Smith 
462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal *);
463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *);
464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
4717b80b807SBarry Smith 
472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
478566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);
4795ef9f2a5SBarry Smith 
480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);
4887b80b807SBarry Smith 
489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,Mat*);
491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJNumeric(MPI_Comm,Mat,PetscInt,Mat);
492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
4998efafbd8SBarry Smith 
500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5017b80b807SBarry Smith 
502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat);
50522440eb1SKris Buschelman 
5067bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
5077bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultSymbolic(Mat,Mat,Mat,PetscReal,Mat*);
5087bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultNumeric(Mat,Mat,Mat,Mat);
5097bab7c10SHong Zhang 
510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat);
513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat);
51622440eb1SKris Buschelman 
517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat);
520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat);
523bc011b1eSHong Zhang 
524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);
5261c741599SHong Zhang 
527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);
5297b80b807SBarry Smith 
530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
540052efed2SBarry Smith 
541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
54390f02eecSBarry Smith 
544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecs(Mat,Vec*,Vec*);
548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
551bd481603SBarry Smith 
552bd481603SBarry Smith /*MC
5531620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5541620fd73SBarry Smith 
5551620fd73SBarry Smith    Not collective
5561620fd73SBarry Smith 
5571620fd73SBarry Smith    Input Parameters:
5581620fd73SBarry Smith +  m - the matrix
5591620fd73SBarry Smith .  row - the row location of the entry
5601620fd73SBarry Smith .  col - the column location of the entry
5611620fd73SBarry Smith .  value - the value to insert
5621620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5631620fd73SBarry Smith 
5641620fd73SBarry Smith    Notes:
5651620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
5661620fd73SBarry Smith    values simultaneously if possible.
5671620fd73SBarry Smith 
5681620fd73SBarry Smith    Level: beginner
5691620fd73SBarry Smith 
5701620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
5711620fd73SBarry Smith M*/
5721620fd73SBarry 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);}
5731620fd73SBarry Smith 
5741620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
5751620fd73SBarry Smith 
5761620fd73SBarry 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);}
5771620fd73SBarry Smith 
5781620fd73SBarry Smith /*MC
5790d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
580bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
581bd481603SBarry Smith 
582bd481603SBarry Smith    Synopsis:
583f2ba6396SBarry Smith    #include "petscmat.h"
584c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
585bd481603SBarry Smith 
586bd481603SBarry Smith    Collective on MPI_Comm
587bd481603SBarry Smith 
588bd481603SBarry Smith    Input Parameters:
589bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
590859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
591859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
592bd481603SBarry Smith 
593bd481603SBarry Smith    Output Parameters:
594bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
595bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
596bd481603SBarry Smith 
597bd481603SBarry Smith 
598bd481603SBarry Smith    Level: intermediate
599bd481603SBarry Smith 
600bd481603SBarry Smith    Notes:
6010598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
602bd481603SBarry Smith 
6031d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
604bd481603SBarry Smith 
605bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
606bd481603SBarry Smith 
6071620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6081620fd73SBarry Smith 
609bd481603SBarry Smith   Concepts: preallocation^Matrix
610bd481603SBarry Smith 
611bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
612bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
613bd481603SBarry Smith M*/
614c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6157c922b88SBarry Smith { \
616a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
617a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
618a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
619eabe889fSLisandro Dalcin   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr); __start = 0; __end = __start; \
620c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
621a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6227c922b88SBarry Smith 
623bd481603SBarry Smith /*MC
624bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
625bd481603SBarry Smith        inserted using a local number of the rows and columns
626bd481603SBarry Smith 
627bd481603SBarry Smith    Synopsis:
628f2ba6396SBarry Smith    #include "petscmat.h"
629c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
630bd481603SBarry Smith 
631bd481603SBarry Smith    Not Collective
632bd481603SBarry Smith 
633bd481603SBarry Smith    Input Parameters:
634784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
635bd481603SBarry Smith .  nrows - the number of rows indicated
6361d73ed98SBarry Smith .  rows - the indices of the rows
637784ac674SJed Brown .  cmap - the column mapping from local to global numbering
638bd481603SBarry Smith .  ncols - the number of columns in the matrix
639bd481603SBarry Smith .  cols - the columns indicated
640bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
641bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
642bd481603SBarry Smith 
643bd481603SBarry Smith 
644bd481603SBarry Smith    Level: intermediate
645bd481603SBarry Smith 
646bd481603SBarry Smith    Notes:
6470598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
648bd481603SBarry Smith 
6491d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
650bd481603SBarry Smith 
651bd481603SBarry Smith   Concepts: preallocation^Matrix
652bd481603SBarry Smith 
653bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
654bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
655bd481603SBarry Smith M*/
656784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
657c4f061fbSSatish Balay {\
658c1ac3661SBarry Smith   PetscInt __l;\
659784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
660784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
661c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
662ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
663c4f061fbSSatish Balay   }\
664c4f061fbSSatish Balay }
665c4f061fbSSatish Balay 
666bd481603SBarry Smith /*MC
667bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
668bd481603SBarry Smith        inserted using a local number of the rows and columns
669bd481603SBarry Smith 
670bd481603SBarry Smith    Synopsis:
671f2ba6396SBarry Smith    #include "petscmat.h"
672c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
673bd481603SBarry Smith 
674bd481603SBarry Smith    Not Collective
675bd481603SBarry Smith 
676bd481603SBarry Smith    Input Parameters:
677bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
678bd481603SBarry Smith .  nrows - the number of rows indicated
6791d73ed98SBarry Smith .  rows - the indices of the rows
680bd481603SBarry Smith .  ncols - the number of columns in the matrix
681bd481603SBarry Smith .  cols - the columns indicated
682bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
683bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
684bd481603SBarry Smith 
685bd481603SBarry Smith 
686bd481603SBarry Smith    Level: intermediate
687bd481603SBarry Smith 
688bd481603SBarry Smith    Notes:
6890598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
690bd481603SBarry Smith 
691bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
692bd481603SBarry Smith 
693bd481603SBarry Smith   Concepts: preallocation^Matrix
694bd481603SBarry Smith 
695bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
696bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
697bd481603SBarry Smith M*/
698d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
699d3d32019SBarry Smith {\
700c1ac3661SBarry Smith   PetscInt __l;\
701d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
702d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
703d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
704d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
705d3d32019SBarry Smith   }\
706d3d32019SBarry Smith }
707d3d32019SBarry Smith 
708bd481603SBarry Smith /*MC
709bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
710bd481603SBarry Smith        inserted using a local number of the rows and columns
711bd481603SBarry Smith 
712bd481603SBarry Smith    Synopsis:
713f2ba6396SBarry Smith    #include "petscmat.h"
714c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
715bd481603SBarry Smith 
716bd481603SBarry Smith    Not Collective
717bd481603SBarry Smith 
718bd481603SBarry Smith    Input Parameters:
71964075487SBarry Smith +  row - the row
720bd481603SBarry Smith .  ncols - the number of columns in the matrix
721a50f8bf6SHong Zhang -  cols - the columns indicated
722a50f8bf6SHong Zhang 
723a50f8bf6SHong Zhang    Output Parameters:
724a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
725bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
726bd481603SBarry Smith 
727bd481603SBarry Smith 
728bd481603SBarry Smith    Level: intermediate
729bd481603SBarry Smith 
730bd481603SBarry Smith    Notes:
7310598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
732bd481603SBarry Smith 
733bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
734bd481603SBarry Smith 
7351620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7361620fd73SBarry Smith 
737bd481603SBarry Smith   Concepts: preallocation^Matrix
738bd481603SBarry Smith 
739bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
740bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
741bd481603SBarry Smith M*/
742c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
743c1ac3661SBarry Smith { PetscInt __i; \
744e32f2f54SBarry 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);\
745e32f2f54SBarry 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);\
7467c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
74764075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
7487cc688d7SBarry Smith     else dnz[row - __rstart]++;\
7497c922b88SBarry Smith   }\
7507c922b88SBarry Smith }
7517c922b88SBarry Smith 
752bd481603SBarry Smith /*MC
753bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
754bd481603SBarry Smith        inserted using a local number of the rows and columns
755bd481603SBarry Smith 
756bd481603SBarry Smith    Synopsis:
757f2ba6396SBarry Smith    #include "petscmat.h"
758c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
759bd481603SBarry Smith 
760bd481603SBarry Smith    Not Collective
761bd481603SBarry Smith 
762bd481603SBarry Smith    Input Parameters:
763bd481603SBarry Smith +  nrows - the number of rows indicated
7641d73ed98SBarry Smith .  rows - the indices of the rows
765bd481603SBarry Smith .  ncols - the number of columns in the matrix
766bd481603SBarry Smith .  cols - the columns indicated
767bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
768bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
769bd481603SBarry Smith 
770bd481603SBarry Smith 
771bd481603SBarry Smith    Level: intermediate
772bd481603SBarry Smith 
773bd481603SBarry Smith    Notes:
7740598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
775bd481603SBarry Smith 
776bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
777bd481603SBarry Smith 
7781620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7791620fd73SBarry Smith 
780bd481603SBarry Smith   Concepts: preallocation^Matrix
781bd481603SBarry Smith 
782bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
783bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
784bd481603SBarry Smith M*/
785d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
786c1ac3661SBarry Smith { PetscInt __i; \
787d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
788d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
789d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
790d3d32019SBarry Smith   }\
791d3d32019SBarry Smith }
792d3d32019SBarry Smith 
793bd481603SBarry Smith /*MC
79416371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
79516371a99SBarry Smith 
79616371a99SBarry Smith    Synopsis:
797f2ba6396SBarry Smith    #include "petscmat.h"
79816371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
79916371a99SBarry Smith 
80016371a99SBarry Smith    Not Collective
80116371a99SBarry Smith 
80216371a99SBarry Smith    Input Parameters:
80316371a99SBarry Smith .  A - matrix
80416371a99SBarry Smith .  row - row where values exist (must be local to this process)
80516371a99SBarry Smith .  ncols - number of columns
80616371a99SBarry Smith .  cols - columns with nonzeros
80716371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
80816371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
80916371a99SBarry Smith 
81016371a99SBarry Smith 
81116371a99SBarry Smith    Level: intermediate
81216371a99SBarry Smith 
81316371a99SBarry Smith    Notes:
8140598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
81516371a99SBarry Smith 
81616371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
81716371a99SBarry Smith 
81816371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
81916371a99SBarry Smith 
82016371a99SBarry Smith   Concepts: preallocation^Matrix
82116371a99SBarry Smith 
82216371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
823eabe889fSLisandro Dalcin           MatPreallocateSymmetricSetLocal()
82416371a99SBarry Smith M*/
8250298fd71SBarry 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);}
82616371a99SBarry Smith 
82716371a99SBarry Smith 
82816371a99SBarry Smith /*MC
8290d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
830bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
831bd481603SBarry Smith 
832bd481603SBarry Smith    Synopsis:
833f2ba6396SBarry Smith    #include "petscmat.h"
834c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
835bd481603SBarry Smith 
836bd481603SBarry Smith    Collective on MPI_Comm
837bd481603SBarry Smith 
838bd481603SBarry Smith    Input Parameters:
83916371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
840bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
841bd481603SBarry Smith 
842bd481603SBarry Smith 
843bd481603SBarry Smith    Level: intermediate
844bd481603SBarry Smith 
845bd481603SBarry Smith    Notes:
8460598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
847bd481603SBarry Smith 
848bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
849bd481603SBarry Smith 
8501620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
8511620fd73SBarry Smith 
852bd481603SBarry Smith   Concepts: preallocation^Matrix
853bd481603SBarry Smith 
854bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
855eabe889fSLisandro Dalcin           MatPreallocateSymmetricSetLocal()
856bd481603SBarry Smith M*/
857a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
8587c922b88SBarry Smith 
859bd481603SBarry Smith 
860bd481603SBarry Smith 
8617b80b807SBarry Smith /* Routines unique to particular data structures */
862014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
863698d4c6aSKris Buschelman 
864014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
865014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
866698d4c6aSKris Buschelman 
867014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
868014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
869014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
870014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
871014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
872014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
8737b80b807SBarry Smith 
874a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
875a96a251dSBarry Smith 
876014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
877014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
878014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
879273d9f13SBarry Smith 
880014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
881014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
882014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
883014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
884014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
885014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
886014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
887014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
888014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
889014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
8909230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
8919230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
892014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
893273d9f13SBarry Smith 
894014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt);
895014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
8961b807ce4Svictorle 
897014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
898014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
8992e8a6d31SBarry Smith 
900014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
9013a7fca6bSBarry Smith 
902014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
9037b80b807SBarry Smith /*
9047b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
90594b7f48cSBarry Smith   done through the KSP and PC interfaces.
9067b80b807SBarry Smith */
9077b80b807SBarry Smith 
90876bdecfbSBarry Smith /*J
909d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
910d9274352SBarry Smith        with an optional dynamic library name, for example
911d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
912d9274352SBarry Smith 
913d9274352SBarry Smith    Level: beginner
914d9274352SBarry Smith 
915e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
916e5a9bf91SBarry Smith 
917d9274352SBarry Smith .seealso: MatGetOrdering()
91876bdecfbSBarry Smith J*/
91919fd82e9SBarry Smith typedef const char* MatOrderingType;
9202692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
9212692d6eeSBarry Smith #define MATORDERINGND          "nd"
9222692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
9232692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
9242692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
9252692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
926312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
927b12f92e5SBarry Smith 
92819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
929140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
930bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
931607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void);
932014dd563SJed Brown PETSC_EXTERN PetscBool         MatOrderingRegisterAllCalled;
933140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList;
934d4fbbf0eSBarry Smith 
935014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
936a2ce50c7SBarry Smith 
937d91e6319SBarry Smith /*S
938d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
939d90ac83dSHong Zhang 
940d90ac83dSHong Zhang    Level: beginner
941d90ac83dSHong Zhang 
942d90ac83dSHong Zhang S*/
943d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
9446a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[];
945d90ac83dSHong Zhang 
946d90ac83dSHong Zhang /*S
9472401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
948b00f7748SHong Zhang 
94961cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
95061cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
951b00f7748SHong Zhang 
95215e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
953b00f7748SHong Zhang 
95461cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
95561cad860SBarry Smith 
956b00f7748SHong Zhang    Level: developer
957b00f7748SHong Zhang 
958d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
959d7d82daaSBarry Smith           MatFactorInfoInitialize()
960b00f7748SHong Zhang 
961b00f7748SHong Zhang S*/
962b00f7748SHong Zhang typedef struct {
96315e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
96485317021SBarry Smith   PetscReal     usedt;
96515e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
966b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
96715e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
96867e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
969348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
970bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
971bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
97215e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
973f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
974f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
97515e8a5b3SHong Zhang } MatFactorInfo;
976ffa6d0a5SLois Curfman McInnes 
977014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
978014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
979014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
980014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
981014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
982014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
983014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
984014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
985014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
986014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
987014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
988014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
989014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
990014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
991014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
992014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
993014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
994014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
995014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
9968ed539a5SBarry Smith 
997014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
998bb5a7306SBarry Smith 
999d91e6319SBarry Smith /*E
1000d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1001bb1eb677SSatish Balay 
1002d91e6319SBarry Smith     Level: beginner
1003d91e6319SBarry Smith 
1004d9274352SBarry Smith    May be bitwise ORd together
1005d9274352SBarry Smith 
1006d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1007d91e6319SBarry Smith 
10084e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10094e7234bfSBarry Smith 
101041f059aeSBarry Smith .seealso: MatSOR()
1011d91e6319SBarry Smith E*/
1012ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1013ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1014ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
101584cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1016014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10178ed539a5SBarry Smith 
1018d4fbbf0eSBarry Smith /*
1019639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1020639f9d9dSBarry Smith */
1021b12f92e5SBarry Smith 
102276bdecfbSBarry Smith /*J
1023d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1024d9274352SBarry Smith        with an optional dynamic library name, for example
1025d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1026d9274352SBarry Smith 
1027d9274352SBarry Smith    Level: beginner
1028d9274352SBarry Smith 
1029d9274352SBarry Smith .seealso: MatGetColoring()
103076bdecfbSBarry Smith J*/
103119fd82e9SBarry Smith typedef const char* MatColoringType;
10322692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
10332692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
10342692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
10352692d6eeSBarry Smith #define MATCOLORINGID      "id"
1036b12f92e5SBarry Smith 
103719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetColoring(Mat,MatColoringType,ISColoring*);
1038bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
103930de9b25SBarry Smith 
1040014dd563SJed Brown PETSC_EXTERN PetscBool MatColoringRegisterAllCalled;
1041f1144a30SSatish Balay 
1042607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void);
1043014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1044639f9d9dSBarry Smith 
1045d9274352SBarry Smith /*S
1046d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1047d9274352SBarry Smith         and coloring
1048639f9d9dSBarry Smith 
1049d9274352SBarry Smith    Level: beginner
1050d9274352SBarry Smith 
1051d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1052d9274352SBarry Smith 
1053d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1054d9274352SBarry Smith S*/
1055e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1056639f9d9dSBarry Smith 
1057014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1058014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1059014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1060014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1061014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1062014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1063014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1064014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1065014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1066014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1067b1683b59SHong Zhang 
1068b1683b59SHong Zhang /*S
1069b9af6bddSHong Zhang      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1070b1683b59SHong Zhang 
1071b1683b59SHong Zhang    Level: beginner
1072b1683b59SHong Zhang 
1073b1683b59SHong Zhang   Concepts: coloring, sparse matrix product
1074b1683b59SHong Zhang 
1075b9af6bddSHong Zhang .seealso:  MatTransposeColoringCreate()
1076b1683b59SHong Zhang S*/
1077b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1078b1683b59SHong Zhang 
1079014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1080014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1081014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1082014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1083b1683b59SHong Zhang 
1084639f9d9dSBarry Smith /*
10850752156aSBarry Smith     These routines are for partitioning matrices: currently used only
10863eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
10870752156aSBarry Smith */
1088ca161407SBarry Smith 
1089d9274352SBarry Smith /*S
1090d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1091d9274352SBarry Smith 
1092d9274352SBarry Smith    Level: beginner
1093d9274352SBarry Smith 
1094d9274352SBarry Smith   Concepts: partitioning
1095d9274352SBarry Smith 
1096743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1097d9274352SBarry Smith S*/
109891e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1099d9274352SBarry Smith 
110076bdecfbSBarry Smith /*J
11015bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1102d9274352SBarry Smith        with an optional dynamic library name, for example
1103d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1104d9274352SBarry Smith 
1105d9274352SBarry Smith    Level: beginner
11060d04baf8SBarry Smith dm
1107b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
110876bdecfbSBarry Smith J*/
110919fd82e9SBarry Smith typedef const char* MatPartitioningType;
11102692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
11112692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
11122692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
11132692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
11142692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
111550d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
111671306c60Sjroman 
1117ca161407SBarry Smith 
1118014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
111919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1120014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1121014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1122014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1123014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1124014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1125014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
11262aabb6bbSBarry Smith 
1127bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
112830de9b25SBarry Smith 
1129014dd563SJed Brown PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled;
1130f1144a30SSatish Balay 
1131607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void);
11322bad1931SBarry Smith 
1133014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1134014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
113519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1136ca161407SBarry Smith 
1137014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1138014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
11390752156aSBarry Smith 
1140b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
11416a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1142b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
11436a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[];
1144b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
11456a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[];
1146b6956c12SJose Roman 
1147014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1148014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1149014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1150014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1151014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1152014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1153014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1154014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1155014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1156014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1157014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
115871306c60Sjroman 
115971306c60Sjroman #define MP_PARTY_OPT "opt"
116071306c60Sjroman #define MP_PARTY_LIN "lin"
116171306c60Sjroman #define MP_PARTY_SCA "sca"
116271306c60Sjroman #define MP_PARTY_RAN "ran"
116371306c60Sjroman #define MP_PARTY_GBF "gbf"
116471306c60Sjroman #define MP_PARTY_GCF "gcf"
116571306c60Sjroman #define MP_PARTY_BUB "bub"
116671306c60Sjroman #define MP_PARTY_DEF "def"
1167014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
116871306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
116971306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
117071306c60Sjroman #define MP_PARTY_NONE "no"
1171014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1172014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1173014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1174014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
117571306c60Sjroman 
117650d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
11776a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1178e0f1cffaSJose Roman 
1179014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1180014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1181014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1182014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
118371306c60Sjroman 
1184b43b03e9SMark F. Adams /*
1185b43b03e9SMark F. Adams     These routines are for coarsening matrices:
1186b43b03e9SMark F. Adams */
1187b43b03e9SMark F. Adams 
1188b43b03e9SMark F. Adams /*S
1189b43b03e9SMark F. Adams      MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix)
1190b43b03e9SMark F. Adams 
1191b43b03e9SMark F. Adams    Level: beginner
1192b43b03e9SMark F. Adams 
1193b43b03e9SMark F. Adams   Concepts: coarsen
1194b43b03e9SMark F. Adams 
1195b43b03e9SMark F. Adams .seealso:  MatCoarsenCreate), MatCoarsenType
1196b43b03e9SMark F. Adams S*/
1197b43b03e9SMark F. Adams typedef struct _p_MatCoarsen* MatCoarsen;
1198b43b03e9SMark F. Adams 
1199b43b03e9SMark F. Adams /*J
1200b43b03e9SMark F. Adams     MatCoarsenType - String with the name of a PETSc matrix coarsen or the creation function
1201b43b03e9SMark F. Adams        with an optional dynamic library name, for example
1202b43b03e9SMark F. Adams        http://www.mcs.anl.gov/petsc/lib.a:coarsencreate()
1203b43b03e9SMark F. Adams 
1204b43b03e9SMark F. Adams    Level: beginner
1205b43b03e9SMark F. Adams dm
1206b43b03e9SMark F. Adams .seealso: MatCoarsenCreate(), MatCoarsen
1207b43b03e9SMark F. Adams J*/
120819fd82e9SBarry Smith typedef const char* MatCoarsenType;
1209b43b03e9SMark F. Adams #define MATCOARSENMIS  "mis"
12109057884aSMark F. Adams #define MATCOARSENHEM  "hem"
1211b43b03e9SMark F. Adams 
12120cbbd2e1SMark F. Adams /* linked list for aggregates */
121341b27cdeSMark F. Adams typedef struct _PetscCDIntNd{
121441b27cdeSMark F. Adams   struct _PetscCDIntNd *next;
12150cbbd2e1SMark F. Adams   PetscInt             gid;
121641b27cdeSMark F. Adams }PetscCDIntNd;
12170cbbd2e1SMark F. Adams 
12180cbbd2e1SMark F. Adams /* only used by node pool */
121941b27cdeSMark F. Adams typedef struct _PetscCDArrNd{
122041b27cdeSMark F. Adams   struct _PetscCDArrNd *next;
122141b27cdeSMark F. Adams   struct _PetscCDIntNd *array;
122241b27cdeSMark F. Adams }PetscCDArrNd;
12230cbbd2e1SMark F. Adams 
12240cbbd2e1SMark F. Adams typedef struct _PetscCoarsenData{
122582c86c8fSBarry Smith   PetscCDArrNd pool_list;  /* node pool */
122641b27cdeSMark F. Adams   PetscCDIntNd *new_node;
12270cbbd2e1SMark F. Adams   PetscInt     new_left;
12280cbbd2e1SMark F. Adams   PetscInt     chk_sz;
122941b27cdeSMark F. Adams   PetscCDIntNd *extra_nodes;
123082c86c8fSBarry Smith   PetscCDIntNd **array;  /* Array of lists */
12310cbbd2e1SMark F. Adams   PetscInt     size;
123282c86c8fSBarry Smith   Mat          mat;  /* cache a Mat for communication data */
12330cbbd2e1SMark F. Adams }PetscCoarsenData;
12340cbbd2e1SMark F. Adams 
1235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*);
123619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,MatCoarsenType);
1237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat);
1238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS);
1239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool);
1240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetVerbose(MatCoarsen,PetscInt);
1241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** );
1242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen);
1243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*);
1244b43b03e9SMark F. Adams 
1245bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],PetscErrorCode (*)(MatCoarsen));
1246b43b03e9SMark F. Adams 
1247014dd563SJed Brown PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled;
1248b43b03e9SMark F. Adams 
1249607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void);
1250b43b03e9SMark F. Adams 
1251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer);
1252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen);
125319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,MatCoarsenType*);
1254b43b03e9SMark F. Adams 
1255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1257591294e4SBarry Smith 
12580752156aSBarry Smith /*
12590a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1260d4fbbf0eSBarry Smith */
12611c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
12621c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
12631c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
12641c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
12651c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
12667c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
12677c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
12681c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
12691c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
12707c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
12717c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
12721c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
12731c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1274a714c06dSBarry Smith                MATOP_SOR=13,
12751c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
12761c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
12771c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
12781c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
12791c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
12801c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
12811c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
12821c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1283d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1284d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1285d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1286d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1287d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1288d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1289d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1290d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1291d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1292d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1293d519adbfSMatthew Knepley                MATOP_GET_ARRAY=32,
1294d519adbfSMatthew Knepley                MATOP_RESTORE_ARRAY=33,
1295d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1296d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1297d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1298d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1299d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1300d519adbfSMatthew Knepley                MATOP_AXPY=39,
1301d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1302d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1303d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1304d519adbfSMatthew Knepley                MATOP_COPY=43,
1305d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1306d519adbfSMatthew Knepley                MATOP_SCALE=45,
1307d519adbfSMatthew Knepley                MATOP_SHIFT=46,
130835153367SBarry Smith                MATOP_DIAGONAL_SET=47,
1309d519adbfSMatthew Knepley                MATOP_ILUDT_FACTOR=48,
1310d519adbfSMatthew Knepley                MATOP_SET_BLOCK_SIZE=49,
1311d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1312d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1313d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1314d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1315d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1316d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1317d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1318d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1319d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1320d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1321d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1322d519adbfSMatthew Knepley                MATOP_VIEW=61,
1323d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
13247bab7c10SHong Zhang                MATOP_MATMAT_MULT=63,
13257bab7c10SHong Zhang                MATOP_MATMAT_MULT_SYMBOLIC=64,
13267bab7c10SHong Zhang                MATOP_MATMAT_MULT_NUMERIC=65,
1327d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1328d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1329d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1330d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1331d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1332d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1333d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
1334bcaeba4dSBarry Smith                MATOP_PLACEHOLDER=73,
1335d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1336d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1337d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1338bda6c4cbSBarry Smith                MATOP_MULT_CONSTRAINED=77,
1339bda6c4cbSBarry Smith                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
1340d519adbfSMatthew Knepley                MATOP_PERMUTE_SPARSIFY=79,
1341d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1342d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1343d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1344d519adbfSMatthew Knepley                MATOP_LOAD=83,
1345d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1346d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1347d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1348820469bcSHong Zhang                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1349d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1350d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1351d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1352d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1353d519adbfSMatthew Knepley                MATOP_PTAP=92,
1354d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1355d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1356bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT=95,
1357bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1358bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
1359820469bcSHong Zhang                MATOP_DUMMY98=98,
1360820469bcSHong Zhang                MATOP_DUMMY99=99,
1361820469bcSHong Zhang                MATOP_DUMMY100=100,
1362820469bcSHong Zhang                MATOP_DUMMY101=101,
1363d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
1364d519adbfSMatthew Knepley                MATOP_SET_SIZES=103,
1365d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1366d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1367bda6c4cbSBarry Smith                MATOP_IMAGINARY_PART=106,
1368bda6c4cbSBarry Smith                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1369bda6c4cbSBarry Smith                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1370bda6c4cbSBarry Smith                MATOP_MAT_SOLVE=109,
13710e8d04f7SBarry Smith                MATOP_GET_REDUNDANT_MATRIX=110,
1372d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
13730e8d04f7SBarry Smith                MATOP_GET_COLUMN_VECTOR=112,
1374d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
13750e8d04f7SBarry Smith                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
137689c6957cSBarry Smith                MATOP_CREATE=115,
137789c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
13780e8d04f7SBarry Smith                MATOP_GET_LOCAL_SUB_MATRIX=117,
13790e8d04f7SBarry Smith                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1380eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
13810e8d04f7SBarry Smith                MATOP_HERMITIAN_TRANSPOSE=120,
13820e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
13830e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
13840e8d04f7SBarry Smith                MATOP_GET_MULTI_PROC_BLOCK=123,
13850e8d04f7SBarry Smith                MATOP_GET_COLUMN_NORMS=125,
13860e8d04f7SBarry Smith                MATOP_GET_SUB_MATRICES_PARALLE=128,
13872b8ad9a3SHong Zhang                MATOP_SET_VALUES_BATCH=129,
13880e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT=130,
13890e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
13900e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
13910e8d04f7SBarry Smith                MATOP_TRANSPOSE_COLORING_CREAT=133,
13920e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_SPT=134,
13930e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_DEN=135,
13940e8d04f7SBarry Smith                MATOP_RART=136,
13950e8d04f7SBarry Smith                MATOP_RART_SYMBOLIC=137,
13960e8d04f7SBarry Smith                MATOP_RART_NUMERIC=138,
1397e09a3074SHong Zhang                MATOP_SET_BLOCK_SIZES=139,
1398e09a3074SHong Zhang                MATOP_AYPX=140
1399fae171e0SBarry Smith              } MatOperation;
1400014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *);
1401014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1402014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1403014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1404112a2221SBarry Smith 
140590ace30eSBarry Smith /*
140690ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
140790ace30eSBarry Smith    stored in a universal format. By changing the format with
14087973ad04SBarry Smith    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
140990ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
141090ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1411f4403165SShri Abhyankar    read into matrices of the same type.
141290ace30eSBarry Smith */
141390ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
141490ace30eSBarry Smith 
1415014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1416014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1417014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
14181f4e1ec7SBarry Smith 
1419d9274352SBarry Smith /*S
1420d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1421d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1422d9274352SBarry Smith 
1423f7a9e4ceSBarry Smith    Level: advanced
1424d9274352SBarry Smith 
1425d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1426d9274352SBarry Smith 
14276e1639daSBarry Smith   Users manual sections:
14286e1639daSBarry Smith .   sec_singular
14296e1639daSBarry Smith 
1430d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1431d9274352SBarry Smith S*/
143274637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1433d9274352SBarry Smith 
1434014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1435014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1436014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1437*d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1438014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1439014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1440014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1441014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1442014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1443014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1444014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1445014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
144674637425SBarry Smith 
1447014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1448014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1449014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1450014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
14513f1d51d7SBarry Smith 
1452014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1455c4f061fbSSatish Balay 
1456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*);
1457b0a32e0cSBarry Smith 
1458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
145904f1ad80SBarry Smith 
1460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace);
1466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1474e884886fSBarry Smith 
14756370053bSBarry Smith /*S
14766370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
14776370053bSBarry Smith               Jacobian vector products
1478e884886fSBarry Smith 
14796370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
14806370053bSBarry Smith 
14816370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
14826370053bSBarry Smith 
14836370053bSBarry Smith     Level: developer
14846370053bSBarry Smith 
14856370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
14866370053bSBarry Smith S*/
1487e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1488e884886fSBarry Smith 
148976bdecfbSBarry Smith /*J
1490e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1491e884886fSBarry Smith 
1492e884886fSBarry Smith    Level: beginner
1493e884886fSBarry Smith 
1494e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
149576bdecfbSBarry Smith J*/
149619fd82e9SBarry Smith typedef const char* MatMFFDType;
1497e884886fSBarry Smith #define MATMFFD_DS  "ds"
1498e884886fSBarry Smith #define MATMFFD_WP  "wp"
1499e884886fSBarry Smith 
150019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1501bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));
1502e884886fSBarry Smith 
1503607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegisterAll(void);
1504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool );
1506e884886fSBarry Smith 
1507014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1508014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
15097dbadf16SMatthew Knepley 
151097969023SHong Zhang /*
151197969023SHong Zhang    PETSc interface to MUMPS
151297969023SHong Zhang */
151397969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1515b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
151697969023SHong Zhang #endif
151723a5497aSJed Brown 
1518d954db57SHong Zhang /*
1519d954db57SHong Zhang    PETSc interface to SUPERLU
1520d954db57SHong Zhang */
1521d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1523d954db57SHong Zhang #endif
1524d954db57SHong Zhang 
15259ae82921SPaul Mullowney #if defined PETSC_HAVE_TXPETSCGPU
1526e057df02SPaul Mullowney 
15273f9c0db1SPaul Mullowney /*E
1528e057df02SPaul Mullowney     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
15298468deeeSKarl Rupp     matrices. Requires the txpetscgpu package to use. Configure with
15308468deeeSKarl Rupp     --download-txpetscgpu to build/install petsc with the txpetscgpu library.
1531e057df02SPaul Mullowney 
1532e057df02SPaul Mullowney     Not Collective
1533e057df02SPaul Mullowney 
15348468deeeSKarl Rupp +   MAT_CUSPARSE_CSR - Compressed Sparse Row
15358468deeeSKarl Rupp .   MAT_CUSPARSE_ELL - Ellpack
15368468deeeSKarl Rupp -   MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format.
1537e057df02SPaul Mullowney 
1538e057df02SPaul Mullowney     Level: intermediate
1539e057df02SPaul Mullowney 
1540e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1541e057df02SPaul Mullowney 
1542e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1543e057df02SPaul Mullowney E*/
1544e057df02SPaul Mullowney 
15453f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;
1546e057df02SPaul Mullowney 
1547e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1548e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
1549e057df02SPaul Mullowney 
15503f9c0db1SPaul Mullowney /*E
1551e057df02SPaul Mullowney     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
15528468deeeSKarl Rupp     matrices whose operation should use a particular storage format. Requires
15538468deeeSKarl Rupp     the txpetscgpu package to use. Configure with --download-txpetscgpu to
15548468deeeSKarl Rupp     build/install petsc with the txpetscgpu library.
1555e057df02SPaul Mullowney 
1556e057df02SPaul Mullowney     Not Collective
1557e057df02SPaul Mullowney 
15588468deeeSKarl Rupp +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
15598468deeeSKarl Rupp .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
15608468deeeSKarl Rupp .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
15618468deeeSKarl Rupp .   MAT_CUSPARSE_SOLVE - sets the storage format for the triangular factors in the serial (single GPU) MatSolve
15628468deeeSKarl Rupp -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices
1563e057df02SPaul Mullowney 
1564e057df02SPaul Mullowney     Level: intermediate
1565e057df02SPaul Mullowney 
1566e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1567e057df02SPaul Mullowney E*/
15683f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_SOLVE, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;
1569e057df02SPaul Mullowney 
157010e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
157110e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1572e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
1573e057df02SPaul Mullowney 
15749ae82921SPaul Mullowney #endif
15759ae82921SPaul Mullowney 
157690c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
1577014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1578014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1579e057df02SPaul Mullowney 
15803f9c0db1SPaul Mullowney /*E
1581e057df02SPaul Mullowney     MatCUSPStorageFormat - indicates the storage format for CUSP (GPU)
15828468deeeSKarl Rupp     matrices. Requires the txpetscgpu package to use. Configure with
15838468deeeSKarl Rupp     --download-txpetscgpu to build/install petsc with the txpetscgpu library.
1584e057df02SPaul Mullowney 
1585e057df02SPaul Mullowney     Not Collective
1586e057df02SPaul Mullowney 
15878468deeeSKarl Rupp +   MAT_CUSP_CSR - Compressed Sparse Row
15888468deeeSKarl Rupp .   MAT_CUSP_DIA - Diagonal
15898468deeeSKarl Rupp -   MAT_CUSP_ELL - Ellpack
1590e057df02SPaul Mullowney 
1591e057df02SPaul Mullowney     Level: intermediate
1592e057df02SPaul Mullowney 
1593e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1594e057df02SPaul Mullowney 
1595e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
1596e057df02SPaul Mullowney E*/
15973f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat;
1598e057df02SPaul Mullowney 
1599e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1600e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[];
1601e057df02SPaul Mullowney 
16023f9c0db1SPaul Mullowney /*E
1603e057df02SPaul Mullowney     MatCUSPFormatOperation - indicates the operation of CUSP (GPU)
16048468deeeSKarl Rupp     matrices whose operation should use a particular storage format. Requires
16058468deeeSKarl Rupp     the txpetscgpu package to use. Configure with --download-txpetscgpu to
16068468deeeSKarl Rupp     build/install petsc with the txpetscgpu library.
1607e057df02SPaul Mullowney 
1608e057df02SPaul Mullowney     Not Collective
1609e057df02SPaul Mullowney 
16108468deeeSKarl Rupp +   MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
16118468deeeSKarl Rupp .   MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
16128468deeeSKarl Rupp .   MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
16138468deeeSKarl Rupp .   MAT_CUSP_SOLVE - sets the storage format for the triangular factors in the serial (single GPU) MatSolve
16148468deeeSKarl Rupp -   MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices
1615e057df02SPaul Mullowney 
1616e057df02SPaul Mullowney     Level: intermediate
1617e057df02SPaul Mullowney 
1618e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1619e057df02SPaul Mullowney 
1620e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat
1621e057df02SPaul Mullowney E*/
16223f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_SOLVE, MAT_CUSP_ALL} MatCUSPFormatOperation;
1623e057df02SPaul Mullowney 
1624e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat);
1625e057df02SPaul Mullowney #endif
162690c53902SBarry Smith 
162754efbe56SHong Zhang /*
162854efbe56SHong Zhang    PETSc interface to FFTW
162954efbe56SHong Zhang */
163054efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
1631014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1632014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
1633014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*);
163454efbe56SHong Zhang #endif
163554efbe56SHong Zhang 
1636db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL)
1637607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
1638db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
1639db31f6deSJed Brown #endif
1640db31f6deSJed Brown 
1641014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1642014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1643014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1644014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1645014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1646014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
164719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1648014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1649014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1650d8588912SDave May 
16514325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
16524325cce7SMatthew G Knepley 
165323a5497aSJed Brown #endif
1654