xref: /petsc/include/petscmat.h (revision 2d033e1fe096f11e11b9bec3b697499293687b30)
12eac72dbSBarry Smith /*
22eac72dbSBarry Smith      Include file for the matrix component of PETSc
32eac72dbSBarry Smith */
40a835dfdSSatish Balay #ifndef __PETSCMAT_H
50a835dfdSSatish Balay #define __PETSCMAT_H
62c8e378dSBarry Smith #include <petscvec.h>
72eac72dbSBarry Smith 
8d9274352SBarry Smith /*S
98f6c3df8SBarry Smith      Mat - Abstract PETSc matrix object used to manage all linear operators in PETSc, even those without
108f6c3df8SBarry Smith            an explicit sparse representation (such as matrix-free operators)
112eac72dbSBarry Smith 
12d91e6319SBarry Smith    Level: beginner
13d91e6319SBarry Smith 
14d9274352SBarry Smith   Concepts: matrix; linear operator
15d9274352SBarry Smith 
168f6c3df8SBarry Smith .seealso:  MatCreate(), MatType, MatSetType(), MatDestroy()
17d9274352SBarry Smith S*/
18d9274352SBarry Smith typedef struct _p_Mat*           Mat;
19d9274352SBarry Smith 
2076bdecfbSBarry Smith /*J
218f6c3df8SBarry Smith     MatType - String with the name of a PETSc matrix type
22d9274352SBarry Smith 
23d9274352SBarry Smith    Level: beginner
24d9274352SBarry Smith 
258f6c3df8SBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage, MatRegister()
2676bdecfbSBarry Smith J*/
2719fd82e9SBarry Smith typedef const char* MatType;
28273d9f13SBarry Smith #define MATSAME            "same"
295a11e1b2SBarry Smith #define MATMAIJ            "maij"
30273d9f13SBarry Smith #define MATSEQMAIJ         "seqmaij"
31273d9f13SBarry Smith #define MATMPIMAIJ         "mpimaij"
32273d9f13SBarry Smith #define MATIS              "is"
335a11e1b2SBarry Smith #define MATAIJ             "aij"
34273d9f13SBarry Smith #define MATSEQAIJ          "seqaij"
35273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
365a11e1b2SBarry Smith #define MATAIJCRL          "aijcrl"
375a11e1b2SBarry Smith #define MATSEQAIJCRL       "seqaijcrl"
385a11e1b2SBarry Smith #define MATMPIAIJCRL       "mpiaijcrl"
398154be41SBarry Smith #define MATAIJCUSP         "aijcusp"
408154be41SBarry Smith #define MATSEQAIJCUSP      "seqaijcusp"
418154be41SBarry Smith #define MATMPIAIJCUSP      "mpiaijcusp"
429ae82921SPaul Mullowney #define MATAIJCUSPARSE     "aijcusparse"
439ae82921SPaul Mullowney #define MATSEQAIJCUSPARSE  "seqaijcusparse"
449ae82921SPaul Mullowney #define MATMPIAIJCUSPARSE  "mpiaijcusparse"
45d67ff14aSKarl Rupp #define MATAIJVIENNACL     "aijviennacl"
46d67ff14aSKarl Rupp #define MATSEQAIJVIENNACL  "seqaijviennacl"
47d67ff14aSKarl Rupp #define MATMPIAIJVIENNACL  "mpiaijviennacl"
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"
63cebc7f6cSBarry Smith #define MATDAAD            "daad"
64cebc7f6cSBarry Smith #define MATMFFD            "mffd"
65c8a8475eSBarry Smith #define MATNORMAL          "normal"
66c5e4d11fSDmitry Karpeev #define MATNORMALHERMITIAN "normalh"
67ab92ecdeSBarry Smith #define MATLRC             "lrc"
682a6744ebSBarry Smith #define MATSCATTER         "scatter"
69421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
70793850ffSBarry Smith #define MATCOMPOSITE       "composite"
711f8c7532SHong Zhang #define MATFFT             "fft"
721f8c7532SHong Zhang #define MATFFTW            "fftw"
73e133240eSMatthew G Knepley #define MATSEQCUFFT        "seqcufft"
74557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
7572ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
761d6018f0SLisandro Dalcin #define MATPYTHON          "python"
7763c07aadSStefano Zampini #define MATHYPRE           "hypre"
78f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
79a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT    "hypresstruct"
80ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
812c0dbf93SJed Brown #define MATLOCALREF        "localref"
82d8588912SDave May #define MATNEST            "nest"
83c094ef40SMatthew G. Knepley #define MATPREALLOCATOR    "preallocator"
84a3b2e22bSHong Zhang #define MATDUMMY           "dummy"
85773941d6SBarry Smith 
8676bdecfbSBarry Smith /*J
87c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
889989ab13SBarry Smith 
89c2b89b5dSBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu_dist" the parallel SuperLU_DIST package etc
909989ab13SBarry Smith 
919989ab13SBarry Smith    Level: beginner
929989ab13SBarry Smith 
935c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
9476bdecfbSBarry Smith J*/
95c7393fdbSBarry Smith #define MatSolverPackage char*
962692d6eeSBarry Smith #define MATSOLVERSUPERLU          "superlu"
972692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST     "superlu_dist"
9808f5efcfSPieter Ghysels #define MATSOLVERSTRUMPACK        "strumpack"
992692d6eeSBarry Smith #define MATSOLVERUMFPACK          "umfpack"
10020db9a53SJed Brown #define MATSOLVERCHOLMOD          "cholmod"
101d89f5e7aSBarry Smith #define MATSOLVERKLU              "klu"
102418810c4SBarry Smith #define MATSOLVERSPARSEELEMENTAL  "sparseelemental"
103d89f5e7aSBarry Smith #define MATSOLVERELEMENTAL        "elemental"
1042692d6eeSBarry Smith #define MATSOLVERESSL             "essl"
1052692d6eeSBarry Smith #define MATSOLVERLUSOL            "lusol"
1062692d6eeSBarry Smith #define MATSOLVERMUMPS            "mumps"
107d615f992SSatish Balay #define MATSOLVERMKL_PARDISO      "mkl_pardiso"
108d305a81bSVasiliy Kozyrev #define MATSOLVERMKL_CPARDISO     "mkl_cpardiso"
1092692d6eeSBarry Smith #define MATSOLVERPASTIX           "pastix"
1102692d6eeSBarry Smith #define MATSOLVERMATLAB           "matlab"
1112692d6eeSBarry Smith #define MATSOLVERPETSC            "petsc"
1122692d6eeSBarry Smith #define MATSOLVERBAS              "bas"
1139ae82921SPaul Mullowney #define MATSOLVERCUSPARSE         "cusparse"
114c0cdd4a1SDahai Guo 
115b24902e0SBarry Smith /*E
116b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
117b24902e0SBarry Smith 
118b24902e0SBarry Smith     Level: beginner
119b24902e0SBarry Smith 
120af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
121b24902e0SBarry Smith 
122c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
123b24902e0SBarry Smith E*/
124599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
125014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[];
126e92e720dSBarry Smith 
127014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
128014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
129014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
130014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*);
13118713533SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister(const MatSolverPackage,const MatType,MatFactorType,PetscErrorCode(*)(Mat,MatFactorType,Mat*));
13242c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageGet(const MatSolverPackage,const MatType,MatFactorType,PetscBool*,PetscBool*,PetscErrorCode (**)(Mat,MatFactorType,Mat*));
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;
137335efc43SPeter Brune PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID;
138014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID;
139014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
140014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID;
141014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID;
142014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID;
143014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID;
144c06d978dSMatthew Knepley 
145ceb03754SKris Buschelman /*E
1467dae84e0SHong Zhang     MatReuse - Indicates if matrices obtained from a previous call to MatCreateSubMatrices(), MatCreateSubMatrix(), MatConvert() or several other functions
147cf37664fSBarry Smith      are to be reused to store the new matrix values.
148cf37664fSBarry Smith 
149cf37664fSBarry Smith $  MAT_INITIAL_MATRIX - create a new matrix
150cf37664fSBarry Smith $  MAT_REUSE_MATRIX - reuse the matrix created with a previous call that used MAT_INITIAL_MATRIX
151cf37664fSBarry Smith $  MAT_INPLACE_MATRIX - replace the first input matrix with the new matrix (not applicable to all functions)
152cf37664fSBarry Smith $  MAT_IGNORE_MATRIX - do not create a new matrix or reuse a give matrix, just ignore that matrix argument (not applicable to all functions)
153ceb03754SKris Buschelman 
154ceb03754SKris Buschelman     Level: beginner
155ceb03754SKris Buschelman 
156af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
157ceb03754SKris Buschelman 
1587dae84e0SHong Zhang .seealso: MatCreateSubMatrices(), MatCreateSubMatrix(), MatDestroyMatrices(), MatConvert()
159ceb03754SKris Buschelman E*/
160511c6705SHong Zhang typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse;
161ceb03754SKris Buschelman 
1625494a064SHong Zhang /*E
1637dae84e0SHong Zhang     MatCreateSubMatrixOption - Indicates if matrices obtained from a call to MatCreateSubMatrices()
164829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1655494a064SHong Zhang 
1665494a064SHong Zhang     Level: beginner
1675494a064SHong Zhang 
168829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1695494a064SHong Zhang E*/
1707dae84e0SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatCreateSubMatrixOption;
1715494a064SHong Zhang 
172607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void);
173c06d978dSMatthew Knepley 
174014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*);
175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
17619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType);
177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat);
178685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
179bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat));
180014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]);
181014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]);
182014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]);
183014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]);
18484d44b13SHong Zhang PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat,PetscBool);
185f69a0ea3SMatthew Knepley 
186140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList;
187140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList;
188140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList;
18928988994SBarry Smith 
1903b224e63SBarry Smith /*E
191aa6c7ce3SBarry Smith     MatStructure - Indicates if two matrices have the same nonzero structure
1923b224e63SBarry Smith 
1933b224e63SBarry Smith     Level: beginner
1943b224e63SBarry Smith 
195af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1963b224e63SBarry Smith 
197aa6c7ce3SBarry Smith .seealso: MatCopy(), MatAXPY()
1983b224e63SBarry Smith E*/
199aa6c7ce3SBarry Smith typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure;
2003b224e63SBarry Smith 
201014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
202014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
206014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2073b00a383SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSeqAIJ(MPI_Comm,Mat,Mat,const PetscInt[],Mat*);
2083b00a383SHong Zhang 
209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
212d21a29f3SJed Brown 
213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
215d21a29f3SJed Brown 
216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
21838f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
219014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
2204cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]);
221dfb205c3SBarry Smith 
222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
224c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*);
225986c4d72SJose E. Roman PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Vec,Mat,Mat*);
226267ca84cSJose E. Roman PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat,Mat*,Mat*,Vec*,Mat*);
227c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*);
228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
230c0cdd4a1SDahai Guo 
231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2386d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
2406d7c1e57SBarry Smith 
24119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*);
242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
243dedccee8SHong Zhang 
244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
2458060fb66Sstefano_zampini PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat,Mat*);
246d0de2241SAndrew Spott PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*);
24754e05e6cSHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrixVirtual(Mat,IS,IS,Mat*);
24854e05e6cSHong Zhang PETSC_EXTERN PetscErrorCode MatSubMatrixVirtualUpdate(Mat,Mat,IS,IS);
249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
2501472f72bSBarry Smith 
251978814f1SStefano Zampini #if defined(PETSC_HAVE_HYPRE)
252d975228cSstefano_zampini PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
253978814f1SStefano Zampini #endif
254978814f1SStefano Zampini 
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]);
2561d6018f0SLisandro Dalcin 
257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
259e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*);
26021c89e3eSBarry Smith 
261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
266713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
26799cafbc1SBarry Smith 
2688ed539a5SBarry Smith /* ------------------------------------------------------------*/
269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
271014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
272014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
273014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
27473a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);
27584cb2905SBarry Smith 
2762ef4de8bSBarry Smith /*S
2772ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
27862ef0227SBarry Smith         column of a matrix as indexed on an associated grid. These are arguments to MatSetStencil() and MatSetBlockStencil()
27962ef0227SBarry Smith 
28062ef0227SBarry Smith    The i,j, and k represent the logical coordinates over the entire grid (for 2 and 1 dimensional problems the k and j entries are ignored).
28162ef0227SBarry Smith    The c represents the the degrees of freedom at each grid point (the dof argument to DMDASetDOF()). If dof is 1 then this entry is ignored.
28262ef0227SBarry Smith 
28362ef0227SBarry Smith    For stencil access to vectors see DMDAVecGetArray(), DMDAVecGetArrayF90().
284be479b30SJed Brown 
285be479b30SJed Brown    Fortran usage is different, see MatSetValuesStencil() for details.
2862ef4de8bSBarry Smith 
2872ef4de8bSBarry Smith    Level: beginner
2882ef4de8bSBarry Smith 
2892ef4de8bSBarry Smith   Concepts: matrix; linear operator
2902ef4de8bSBarry Smith 
29162ef0227SBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90()
2922ef4de8bSBarry Smith S*/
293435da068SBarry Smith typedef struct {
294c1ac3661SBarry Smith   PetscInt k,j,i,c;
295435da068SBarry Smith } MatStencil;
2962ef4de8bSBarry Smith 
297014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
298014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
299014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
300435da068SBarry Smith 
301d91e6319SBarry Smith /*E
302d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
30362ef0227SBarry Smith      to continue to add or insert values to it
304d91e6319SBarry Smith 
305d91e6319SBarry Smith     Level: beginner
306d91e6319SBarry Smith 
307d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
308d91e6319SBarry Smith E*/
3096d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
310014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
311014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
312014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);
3134f9c727eSBarry Smith 
3141d73ed98SBarry Smith 
31530de9b25SBarry Smith 
316d91e6319SBarry Smith /*E
317d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
318d91e6319SBarry Smith 
319d91e6319SBarry Smith     Level: beginner
320d91e6319SBarry Smith 
321af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
3220f8fb01aSBarry Smith    Any additions/changes here must also be made in src/mat/interface/dlregismat.c in MatOptions[]
323d91e6319SBarry Smith 
324382164b0SBarry Smith    Developer Notes: Entries that are negative need not be called collectively by all processes.
325382164b0SBarry Smith 
326d91e6319SBarry Smith .seealso: MatSetOption()
327d91e6319SBarry Smith E*/
328c5e4d11fSDmitry Karpeev typedef enum {MAT_OPTION_MIN = -3,
329c5e4d11fSDmitry Karpeev               MAT_UNUSED_NONZERO_LOCATION_ERR = -2,
33092d4d306SBarry Smith               MAT_ROW_ORIENTED = -1,
331382164b0SBarry Smith               MAT_SYMMETRIC = 1,
332382164b0SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC = 2,
333382164b0SBarry Smith               MAT_NEW_DIAGONALS = 3,
334382164b0SBarry Smith               MAT_IGNORE_OFF_PROC_ENTRIES = 4,
335382164b0SBarry Smith               MAT_USE_HASH_TABLE = 5,
336382164b0SBarry Smith               MAT_KEEP_NONZERO_PATTERN = 6,
337382164b0SBarry Smith               MAT_IGNORE_ZERO_ENTRIES = 7,
338382164b0SBarry Smith               MAT_USE_INODES = 8,
339382164b0SBarry Smith               MAT_HERMITIAN = 9,
340382164b0SBarry Smith               MAT_SYMMETRY_ETERNAL = 10,
341c5e4d11fSDmitry Karpeev               MAT_NEW_NONZERO_LOCATION_ERR = 11,
342382164b0SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR = 12,
343382164b0SBarry Smith               MAT_ERROR_LOWER_TRIANGULAR = 13,
344382164b0SBarry Smith               MAT_GETROW_UPPERTRIANGULAR = 14,
345382164b0SBarry Smith               MAT_SPD = 15,
34692d4d306SBarry Smith               MAT_NO_OFF_PROC_ZERO_ROWS = 16,
34792d4d306SBarry Smith               MAT_NO_OFF_PROC_ENTRIES = 17,
34892d4d306SBarry Smith               MAT_NEW_NONZERO_LOCATIONS = 18,
349c5e4d11fSDmitry Karpeev               MAT_NEW_NONZERO_ALLOCATION_ERR = 19,
350c5e4d11fSDmitry Karpeev               MAT_SUBSET_OFF_PROC_ENTRIES = 20,
351c10200c1SHong Zhang               MAT_SUBMAT_SINGLEIS = 21,
352957cac9fSHong Zhang               MAT_STRUCTURE_ONLY = 22,
353957cac9fSHong Zhang               MAT_OPTION_MAX = 23} MatOption;
354e057df02SPaul Mullowney 
3550f8fb01aSBarry Smith PETSC_EXTERN const char *const *MatOptions;
356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool);
3577d68702bSBarry Smith PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*);
35819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);
35984cb2905SBarry Smith 
360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
363014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
364014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
365014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
3668c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
3678c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
36821e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*);
369bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3704099cc6bSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetType(Mat,MatType);
3714099cc6bSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRegister(const char[],PetscErrorCode (*)(Mat,const MatType,MatReuse,Mat *));
3724099cc6bSBarry Smith PETSC_EXTERN PetscFunctionList MatSeqAIJList;
3738397e458SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat,PetscScalar *[]);
3748397e458SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat,PetscScalar *[]);
3758c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
3768c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
377d3042a70SBarry Smith PETSC_EXTERN PetscErrorCode MatDensePlaceArray(Mat,const PetscScalar[]);
378d3042a70SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseResetArray(Mat);
379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
38333d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat);
3841620fd73SBarry Smith 
385014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
386014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
387014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
388014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
389014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
390014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
391014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
392014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
393014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
394014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
395014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
396014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
397bdc285e1SStefano Zampini PETSC_EXTERN PetscErrorCode MatMatSolveTranspose(Mat,Mat,Mat);
398f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);
399f5edf698SHong Zhang 
400d91e6319SBarry Smith /*E
401d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
402d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
403d91e6319SBarry Smith 
404d91e6319SBarry Smith     Level: beginner
405d91e6319SBarry Smith 
406af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
407d91e6319SBarry Smith 
40870dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
40970dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
41070dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
41170dcbbb9SBarry Smith 
412d91e6319SBarry Smith .seealso: MatDuplicate()
413d91e6319SBarry Smith E*/
41470dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
4152e8a6d31SBarry Smith 
41619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
417014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);
41894a9d846SBarry Smith 
419cb5b572fSBarry Smith 
420014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
421014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
422014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
423014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
424014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
425014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
426014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
427014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
428014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
4297b80b807SBarry Smith 
4301a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4311a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
4321a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4331a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
434d4fbbf0eSBarry Smith 
435d91e6319SBarry Smith /*S
436d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
437d91e6319SBarry Smith 
438d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
439d91e6319SBarry Smith 
440d91e6319SBarry Smith    Level: intermediate
441d91e6319SBarry Smith 
442d91e6319SBarry Smith   Concepts: matrix^nonzero information
443d91e6319SBarry Smith 
444d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
445d91e6319SBarry Smith S*/
4464e220ebcSLois Curfman McInnes typedef struct {
447b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
448b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
449b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
450b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
451b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
452b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
453b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4544e220ebcSLois Curfman McInnes } MatInfo;
4554e220ebcSLois Curfman McInnes 
456d9274352SBarry Smith /*E
457d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
458d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
459d9274352SBarry Smith 
460d9274352SBarry Smith     Level: beginner
461d9274352SBarry Smith 
462af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
463d9274352SBarry Smith 
464d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
465d9274352SBarry Smith E*/
4667b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*);
477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);
479a52676ecSHong Zhang 
480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*);
481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*);
482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*);
483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*);
484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*);
485a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
486a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
4877b80b807SBarry Smith 
488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*);
489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*);
490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
4977b80b807SBarry Smith 
498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
504566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);
5055ef9f2a5SBarry Smith 
5067dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
507cd2534ddSHong Zhang PETSC_DEPRECATED("Use MatCreateSubMatrices()") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatrices(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) {return MatCreateSubMatrices(mat,n,irow,icol,scall,submat);}
5087dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
509cd2534ddSHong Zhang PETSC_DEPRECATED("Use MatCreateSubMatricesMPI()") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatricesMPI(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) {return MatCreateSubMatricesMPI(mat,n,irow,icol,scall,submat);}
510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
511df750dc8SHong Zhang PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]);
5127dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *);
513cd2534ddSHong Zhang PETSC_DEPRECATED("Use MatCreateSubMatrix()") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatrix(Mat mat,IS isrow,IS iscol,MatReuse cll,Mat *newmat) {return MatCreateSubMatrix(mat,isrow,iscol,cll,newmat);}
514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);
5187b80b807SBarry Smith 
519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5268efafbd8SBarry Smith 
527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
528aa1e27eaSFande Kong PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov);
529d2b2a242SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool);
5307b80b807SBarry Smith 
531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat);
53422440eb1SKris Buschelman 
5357bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
5362df6c5c3SPatrick Farrell PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
5377bab7c10SHong Zhang 
538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat);
541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat);
54422440eb1SKris Buschelman 
545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat);
548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat);
551bc011b1eSHong Zhang 
552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);
5541c741599SHong Zhang 
555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);
5577b80b807SBarry Smith 
558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
560a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*);
561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
565014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
567052efed2SBarry Smith 
568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
569014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
57090f02eecSBarry Smith 
571014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
572014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
573014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
5742a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*);
57584cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);}
57653cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
577014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
578014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
5793a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
5809c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
581bd481603SBarry Smith 
582bd481603SBarry Smith /*MC
5831620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5841620fd73SBarry Smith 
5851620fd73SBarry Smith    Not collective
5861620fd73SBarry Smith 
587a9834a7dSSatish Balay    Synopsis:
588a9834a7dSSatish Balay      #include <petscmat.h>
589a9834a7dSSatish Balay      PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode)
590a9834a7dSSatish Balay 
5911620fd73SBarry Smith    Input Parameters:
5921620fd73SBarry Smith +  m - the matrix
5931620fd73SBarry Smith .  row - the row location of the entry
5941620fd73SBarry Smith .  col - the column location of the entry
5951620fd73SBarry Smith .  value - the value to insert
5961620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5971620fd73SBarry Smith 
5981620fd73SBarry Smith    Notes:
5991620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6001620fd73SBarry Smith    values simultaneously if possible.
6011620fd73SBarry Smith 
6021620fd73SBarry Smith    Level: beginner
6031620fd73SBarry Smith 
6041620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6051620fd73SBarry Smith M*/
6061620fd73SBarry 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);}
6071620fd73SBarry Smith 
6081620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6091620fd73SBarry Smith 
6101620fd73SBarry 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);}
6111620fd73SBarry Smith 
6121620fd73SBarry Smith /*MC
6130d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
614bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
615bd481603SBarry Smith 
616bd481603SBarry Smith    Synopsis:
617aaa7dc30SBarry Smith    #include <petscmat.h>
618c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
619bd481603SBarry Smith 
620bd481603SBarry Smith    Collective on MPI_Comm
621bd481603SBarry Smith 
622bd481603SBarry Smith    Input Parameters:
623bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
624859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
625859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
626bd481603SBarry Smith 
627bd481603SBarry Smith    Output Parameters:
628bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
629bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
630bd481603SBarry Smith 
631bd481603SBarry Smith    Level: intermediate
632bd481603SBarry Smith 
633bd481603SBarry Smith    Notes:
634a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
635bd481603SBarry Smith 
6361d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
637bd481603SBarry Smith 
638bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
639bd481603SBarry Smith 
6401620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6411620fd73SBarry Smith 
642bd481603SBarry Smith   Concepts: preallocation^Matrix
643bd481603SBarry Smith 
644d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
645d6e23781SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock()
646bd481603SBarry Smith M*/
647c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6487c922b88SBarry Smith { \
649a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
6501795a4d1SJed Brown   _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \
6511795a4d1SJed Brown   __start = 0; __end = __start;                                         \
652c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
653a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6547c922b88SBarry Smith 
655bd481603SBarry Smith /*MC
656bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
657bd481603SBarry Smith        inserted using a local number of the rows and columns
658bd481603SBarry Smith 
659bd481603SBarry Smith    Synopsis:
660aaa7dc30SBarry Smith    #include <petscmat.h>
661c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
662bd481603SBarry Smith 
663bd481603SBarry Smith    Not Collective
664bd481603SBarry Smith 
665bd481603SBarry Smith    Input Parameters:
666784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
667bd481603SBarry Smith .  nrows - the number of rows indicated
6681d73ed98SBarry Smith .  rows - the indices of the rows
669784ac674SJed Brown .  cmap - the column mapping from local to global numbering
670bd481603SBarry Smith .  ncols - the number of columns in the matrix
671bd481603SBarry Smith .  cols - the columns indicated
672bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
673bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
674bd481603SBarry Smith 
675bd481603SBarry Smith    Level: intermediate
676bd481603SBarry Smith 
677bd481603SBarry Smith    Notes:
678a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
679bd481603SBarry Smith 
6801d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
681bd481603SBarry Smith 
682bd481603SBarry Smith   Concepts: preallocation^Matrix
683bd481603SBarry Smith 
684d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
685c1154cd5SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups()
686bd481603SBarry Smith M*/
687784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
688c4f061fbSSatish Balay {\
689c1ac3661SBarry Smith   PetscInt __l;\
690784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
691784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
692c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
693ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
694c4f061fbSSatish Balay   }\
695c4f061fbSSatish Balay }
696c4f061fbSSatish Balay 
697bd481603SBarry Smith /*MC
698c1154cd5SBarry Smith    MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be
699c1154cd5SBarry Smith        inserted using a local number of the rows and columns. This version removes any duplicate columns in cols
700c1154cd5SBarry Smith 
701c1154cd5SBarry Smith    Synopsis:
702c1154cd5SBarry Smith    #include <petscmat.h>
703c1154cd5SBarry Smith    PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
704c1154cd5SBarry Smith 
705c1154cd5SBarry Smith    Not Collective
706c1154cd5SBarry Smith 
707c1154cd5SBarry Smith    Input Parameters:
708c1154cd5SBarry Smith +  map - the row mapping from local numbering to global numbering
709c1154cd5SBarry Smith .  nrows - the number of rows indicated
710c1154cd5SBarry Smith .  rows - the indices of the rows (these values are mapped to the global values)
711c1154cd5SBarry Smith .  cmap - the column mapping from local to global numbering
712c1154cd5SBarry Smith .  ncols - the number of columns in the matrix   (this value will be changed if duplicate columns are found)
713c1154cd5SBarry Smith .  cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed)
714c1154cd5SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
715c1154cd5SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
716c1154cd5SBarry Smith 
717c1154cd5SBarry Smith    Level: intermediate
718c1154cd5SBarry Smith 
719c1154cd5SBarry Smith    Notes:
720c1154cd5SBarry Smith     See Users-Manual: ch_performance for more details.
721c1154cd5SBarry Smith 
722c1154cd5SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
723c1154cd5SBarry Smith 
724c1154cd5SBarry Smith   Concepts: preallocation^Matrix
725c1154cd5SBarry Smith 
726c1154cd5SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
727c1154cd5SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
728c1154cd5SBarry Smith M*/
729c1154cd5SBarry Smith #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
730c1154cd5SBarry Smith {\
731c1154cd5SBarry Smith   PetscInt __l;\
732c1154cd5SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
733c1154cd5SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
734c1154cd5SBarry Smith   _4_ierr = PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\
735c1154cd5SBarry Smith   for (__l=0;__l<nrows;__l++) {\
736c1154cd5SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
737c1154cd5SBarry Smith   }\
738c1154cd5SBarry Smith }
739c1154cd5SBarry Smith 
740c1154cd5SBarry Smith /*MC
741d6e23781SBarry Smith    MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
742bd481603SBarry Smith        inserted using a local number of the rows and columns
743bd481603SBarry Smith 
744bd481603SBarry Smith    Synopsis:
745aaa7dc30SBarry Smith    #include <petscmat.h>
746d6e23781SBarry Smith    PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
747d6e23781SBarry Smith 
748d6e23781SBarry Smith    Not Collective
749d6e23781SBarry Smith 
750d6e23781SBarry Smith    Input Parameters:
751d6e23781SBarry Smith +  map - the row mapping from local numbering to global numbering
752d6e23781SBarry Smith .  nrows - the number of rows indicated
753d6e23781SBarry Smith .  rows - the indices of the rows
754d6e23781SBarry Smith .  cmap - the column mapping from local to global numbering
755d6e23781SBarry Smith .  ncols - the number of columns in the matrix
756d6e23781SBarry Smith .  cols - the columns indicated
757d6e23781SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
758d6e23781SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
759d6e23781SBarry Smith 
760d6e23781SBarry Smith    Level: intermediate
761d6e23781SBarry Smith 
762d6e23781SBarry Smith    Notes:
763d6e23781SBarry Smith     See Users-Manual: ch_performance for more details.
764d6e23781SBarry Smith 
765d6e23781SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
766d6e23781SBarry Smith 
767d6e23781SBarry Smith   Concepts: preallocation^Matrix
768d6e23781SBarry Smith 
769d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
770d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
771d6e23781SBarry Smith M*/
772d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
773d6e23781SBarry Smith {\
774d6e23781SBarry Smith   PetscInt __l;\
775d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
776d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
777d6e23781SBarry Smith   for (__l=0;__l<nrows;__l++) {\
778d6e23781SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
779d6e23781SBarry Smith   }\
780d6e23781SBarry Smith }
781d6e23781SBarry Smith 
782d6e23781SBarry Smith /*MC
783d6e23781SBarry Smith    MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
784d6e23781SBarry Smith        inserted using a local number of the rows and columns
785d6e23781SBarry Smith 
786d6e23781SBarry Smith    Synopsis:
787d6e23781SBarry Smith    #include <petscmat.h>
788d6e23781SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
789bd481603SBarry Smith 
790bd481603SBarry Smith    Not Collective
791bd481603SBarry Smith 
792bd481603SBarry Smith    Input Parameters:
793bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
794bd481603SBarry Smith .  nrows - the number of rows indicated
7951d73ed98SBarry Smith .  rows - the indices of the rows
796bd481603SBarry Smith .  ncols - the number of columns in the matrix
797bd481603SBarry Smith .  cols - the columns indicated
798bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
799bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
800bd481603SBarry Smith 
801bd481603SBarry Smith    Level: intermediate
802bd481603SBarry Smith 
803bd481603SBarry Smith    Notes:
804a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
805bd481603SBarry Smith 
806bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
807bd481603SBarry Smith 
808bd481603SBarry Smith   Concepts: preallocation^Matrix
809bd481603SBarry Smith 
810d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(),
811d6e23781SBarry Smith           MatPreallocateInitialize(),  MatPreallocateSetLocal()
812bd481603SBarry Smith M*/
813d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\
814d3d32019SBarry Smith {\
815c1ac3661SBarry Smith   PetscInt __l;\
816d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
817d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
818d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
819d6e23781SBarry Smith     _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
820d3d32019SBarry Smith   }\
821d3d32019SBarry Smith }
822bd481603SBarry Smith /*MC
823bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
824bd481603SBarry Smith        inserted using a local number of the rows and columns
825bd481603SBarry Smith 
826bd481603SBarry Smith    Synopsis:
827aaa7dc30SBarry Smith    #include <petscmat.h>
828c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
829bd481603SBarry Smith 
830bd481603SBarry Smith    Not Collective
831bd481603SBarry Smith 
832bd481603SBarry Smith    Input Parameters:
83364075487SBarry Smith +  row - the row
834bd481603SBarry Smith .  ncols - the number of columns in the matrix
835a50f8bf6SHong Zhang -  cols - the columns indicated
836a50f8bf6SHong Zhang 
837a50f8bf6SHong Zhang    Output Parameters:
838a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
839bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
840bd481603SBarry Smith 
841bd481603SBarry Smith    Level: intermediate
842bd481603SBarry Smith 
843bd481603SBarry Smith    Notes:
844a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
845bd481603SBarry Smith 
846bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
847bd481603SBarry Smith 
8481620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8491620fd73SBarry Smith 
850bd481603SBarry Smith   Concepts: preallocation^Matrix
851bd481603SBarry Smith 
852d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
853d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSetLocal()
854bd481603SBarry Smith M*/
855c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
856c1ac3661SBarry Smith { PetscInt __i; \
857e32f2f54SBarry 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);\
858e32f2f54SBarry 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);\
8597c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
86064075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8617cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8627c922b88SBarry Smith   }\
8637c922b88SBarry Smith }
8647c922b88SBarry Smith 
865bd481603SBarry Smith /*MC
866d6e23781SBarry Smith    MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
867bd481603SBarry Smith        inserted using a local number of the rows and columns
868bd481603SBarry Smith 
869bd481603SBarry Smith    Synopsis:
870aaa7dc30SBarry Smith    #include <petscmat.h>
871d6e23781SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
872bd481603SBarry Smith 
873bd481603SBarry Smith    Not Collective
874bd481603SBarry Smith 
875bd481603SBarry Smith    Input Parameters:
876bd481603SBarry Smith +  nrows - the number of rows indicated
8771d73ed98SBarry Smith .  rows - the indices of the rows
878bd481603SBarry Smith .  ncols - the number of columns in the matrix
879bd481603SBarry Smith .  cols - the columns indicated
880bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
881bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
882bd481603SBarry Smith 
883bd481603SBarry Smith    Level: intermediate
884bd481603SBarry Smith 
885bd481603SBarry Smith    Notes:
886a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
887bd481603SBarry Smith 
888bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
889bd481603SBarry Smith 
8901620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8911620fd73SBarry Smith 
892bd481603SBarry Smith   Concepts: preallocation^Matrix
893bd481603SBarry Smith 
894d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(),  MatPreallocateInitialize(),
895d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
896bd481603SBarry Smith M*/
897d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\
898c1ac3661SBarry Smith { PetscInt __i; \
899d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
900d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
901d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
902d3d32019SBarry Smith   }\
903d3d32019SBarry Smith }
904d3d32019SBarry Smith 
905bd481603SBarry Smith /*MC
90616371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
90716371a99SBarry Smith 
90816371a99SBarry Smith    Synopsis:
909aaa7dc30SBarry Smith    #include <petscmat.h>
91016371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
91116371a99SBarry Smith 
91216371a99SBarry Smith    Not Collective
91316371a99SBarry Smith 
91416371a99SBarry Smith    Input Parameters:
91516371a99SBarry Smith .  A - matrix
91616371a99SBarry Smith .  row - row where values exist (must be local to this process)
91716371a99SBarry Smith .  ncols - number of columns
91816371a99SBarry Smith .  cols - columns with nonzeros
91916371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
92016371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
92116371a99SBarry Smith 
92216371a99SBarry Smith    Level: intermediate
92316371a99SBarry Smith 
92416371a99SBarry Smith    Notes:
925a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
92616371a99SBarry Smith 
92716371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
92816371a99SBarry Smith 
92916371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
93016371a99SBarry Smith 
93116371a99SBarry Smith   Concepts: preallocation^Matrix
93216371a99SBarry Smith 
933d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
934d6e23781SBarry Smith           MatPreallocateSymmetricSetLocalBlock()
93516371a99SBarry Smith M*/
9360298fd71SBarry 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);}
93716371a99SBarry Smith 
93816371a99SBarry Smith 
93916371a99SBarry Smith /*MC
9400d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
941bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
942bd481603SBarry Smith 
943bd481603SBarry Smith    Synopsis:
944aaa7dc30SBarry Smith    #include <petscmat.h>
945c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
946bd481603SBarry Smith 
947bd481603SBarry Smith    Collective on MPI_Comm
948bd481603SBarry Smith 
949bd481603SBarry Smith    Input Parameters:
95016371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
951bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
952bd481603SBarry Smith 
953bd481603SBarry Smith    Level: intermediate
954bd481603SBarry Smith 
955bd481603SBarry Smith    Notes:
956a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
957bd481603SBarry Smith 
958bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
959bd481603SBarry Smith 
9601620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9611620fd73SBarry Smith 
962bd481603SBarry Smith   Concepts: preallocation^Matrix
963bd481603SBarry Smith 
964d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
965d6e23781SBarry Smith           MatPreallocateSymmetricSetLocalBlock()
966bd481603SBarry Smith M*/
967a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9687c922b88SBarry Smith 
9697b80b807SBarry Smith /* Routines unique to particular data structures */
970014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
971698d4c6aSKris Buschelman 
972014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
973014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
974698d4c6aSKris Buschelman 
975014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
976014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
977014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
978014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
979014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
980014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
9817b80b807SBarry Smith 
982a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
983a96a251dSBarry Smith 
984014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
985014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
986014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
987273d9f13SBarry Smith 
988014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
989014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
990014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
991014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
992014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
993014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
994014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
995014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
996014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
997014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
9989230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
9999230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1000014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
1001273d9f13SBarry Smith 
10022e1947a5SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1003cf0a3239SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetUpSF(Mat);
1004014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt);
1005014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
10061b807ce4Svictorle 
1007014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
1008014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
10092e8a6d31SBarry Smith 
1010014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
10113a7fca6bSBarry Smith 
1012014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
10137cfe41e4SPatrick Farrell PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*);
10147b80b807SBarry Smith /*
10157b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
101694b7f48cSBarry Smith   done through the KSP and PC interfaces.
10177b80b807SBarry Smith */
10187b80b807SBarry Smith 
101976bdecfbSBarry Smith /*J
10208f6c3df8SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering
1021d9274352SBarry Smith 
1022d9274352SBarry Smith    Level: beginner
1023d9274352SBarry Smith 
1024d9274352SBarry Smith .seealso: MatGetOrdering()
102576bdecfbSBarry Smith J*/
102619fd82e9SBarry Smith typedef const char* MatOrderingType;
10272692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
10282692d6eeSBarry Smith #define MATORDERINGND          "nd"
10292692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
10302692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
10312692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
10322692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
1033510184a7SJed Brown #define MATORDERINGWBM         "wbm"
1034fc1bef75SJed Brown #define MATORDERINGSPECTRAL    "spectral"
1035312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
1036b12f92e5SBarry Smith 
103719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
1038140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
1039bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
1040140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList;
1041d4fbbf0eSBarry Smith 
1042014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1043fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);
1044a2ce50c7SBarry Smith 
1045d91e6319SBarry Smith /*S
1046d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
1047d90ac83dSHong Zhang 
1048d90ac83dSHong Zhang    Level: beginner
1049d90ac83dSHong Zhang 
1050d90ac83dSHong Zhang S*/
1051d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
10526a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[];
10535e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];
1054d90ac83dSHong Zhang 
1055d90ac83dSHong Zhang /*S
10569aa7eafdSHong Zhang     MatFactorError - indicates what type of error in matrix factor
10579aa7eafdSHong Zhang 
10589aa7eafdSHong Zhang     Level: beginner
10599aa7eafdSHong Zhang 
106000e125f8SBarry Smith     Developer Notes: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
106100e125f8SBarry Smith 
106200e125f8SBarry Smith .seealso: MatGetFactor()
10639aa7eafdSHong Zhang S*/
10645cd7cf9dSHong Zhang typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError;
10659aa7eafdSHong Zhang 
106600e125f8SBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*);
1067b8b68cfdSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat);
10687b6c816cSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*);
106900e125f8SBarry Smith 
10709aa7eafdSHong Zhang /*S
1071422a814eSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization
1072b00f7748SHong Zhang 
107361cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
107461cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1075b00f7748SHong Zhang 
107615e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1077b00f7748SHong Zhang 
107861cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
107961cad860SBarry Smith 
1080b00f7748SHong Zhang    Level: developer
1081b00f7748SHong Zhang 
1082d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1083d7d82daaSBarry Smith           MatFactorInfoInitialize()
1084b00f7748SHong Zhang 
1085b00f7748SHong Zhang S*/
1086b00f7748SHong Zhang typedef struct {
108715e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
108885317021SBarry Smith   PetscReal     usedt;
108915e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1090b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
109115e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
109267e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1093348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1094bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1095bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
109615e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1097f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1098f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
109915e8a5b3SHong Zhang } MatFactorInfo;
1100ffa6d0a5SLois Curfman McInnes 
1101014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
11029a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
11039a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11049a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
11059a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
11069a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
11079a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11089a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11099a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11109a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
11119a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1112014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1113014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1114014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1115014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1116014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1117014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1118014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1119014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
11207c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
11217c2f51b8SStefano Zampini 
11227c2f51b8SStefano Zampini typedef enum {MAT_FACTOR_SCHUR_UNFACTORED, MAT_FACTOR_SCHUR_FACTORED, MAT_FACTOR_SCHUR_INVERTED} MatFactorSchurStatus;
11238e7ba810SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS);
11247c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
11257c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*,MatFactorSchurStatus);
11265a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat);
11277c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
11285a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec);
11295a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec);
11306dba178dSStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat);
1131bb5a7306SBarry Smith 
1132d91e6319SBarry Smith /*E
1133d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1134bb1eb677SSatish Balay 
1135d91e6319SBarry Smith     Level: beginner
1136d91e6319SBarry Smith 
1137d9274352SBarry Smith    May be bitwise ORd together
1138d9274352SBarry Smith 
1139af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1140d91e6319SBarry Smith 
11414e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11424e7234bfSBarry Smith 
114341f059aeSBarry Smith .seealso: MatSOR()
1144d91e6319SBarry Smith E*/
1145ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1146ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1147ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
114884cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1149014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11508ed539a5SBarry Smith 
1151d4fbbf0eSBarry Smith /*
1152639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1153639f9d9dSBarry Smith */
1154b12f92e5SBarry Smith 
11557086a01eSPeter Brune /*S
11567086a01eSPeter Brune      MatColoring - Object for managing the coloring of matrices.
11577086a01eSPeter Brune 
11587086a01eSPeter Brune    Level: beginner
11597086a01eSPeter Brune 
11607086a01eSPeter Brune   Concepts: matrix, coloring
11617086a01eSPeter Brune 
11627086a01eSPeter Brune .seealso:  MatFDColoringCreate() ISColoring MatFDColoring
11637086a01eSPeter Brune S*/
11647086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring;
116576bdecfbSBarry Smith /*J
11668f6c3df8SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring
1167d9274352SBarry Smith 
1168d9274352SBarry Smith    Level: beginner
1169d9274352SBarry Smith 
11707086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring
117176bdecfbSBarry Smith J*/
11727086a01eSPeter Brune 
117319fd82e9SBarry Smith typedef const  char*           MatColoringType;
11745567d71aSPeter Brune #define MATCOLORINGJP      "jp"
1175b9acb617SPeter Brune #define MATCOLORINGPOWER   "power"
11762692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
11772692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
11782692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
11792692d6eeSBarry Smith #define MATCOLORINGID      "id"
11808121bdceSPeter Brune #define MATCOLORINGGREEDY  "greedy"
1181b12f92e5SBarry Smith 
11828ac6da0aSPeter Brune /*E
11838ac6da0aSPeter Brune    MatColoringWeightType - Type of weight scheme
11848ac6da0aSPeter Brune 
11858ac6da0aSPeter Brune     Not Collective
11868ac6da0aSPeter Brune 
11878ac6da0aSPeter Brune +   MAT_COLORING_RANDOM  - Random weights
11888ac6da0aSPeter Brune .   MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
11898ac6da0aSPeter Brune -   MAT_COLORING_LF      - Last-first weighting.
11908ac6da0aSPeter Brune 
11918ac6da0aSPeter Brune     Level: intermediate
11928ac6da0aSPeter Brune 
1193af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
11948ac6da0aSPeter Brune 
11958ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
11968ac6da0aSPeter Brune E*/
1197301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;
11988ac6da0aSPeter Brune 
1199335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1200d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*);
1201335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1202335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1203335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1204335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1205335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1206335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1207335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1208335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1209335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1210335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
12128ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType);
1213c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*);
12148ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm);
1215cb394dc5SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringTestValid(MatColoring,ISColoring);
1216639f9d9dSBarry Smith 
1217d9274352SBarry Smith /*S
1218d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1219d9274352SBarry Smith         and coloring
1220639f9d9dSBarry Smith 
1221d9274352SBarry Smith    Level: beginner
1222d9274352SBarry Smith 
1223d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1224d9274352SBarry Smith 
1225d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1226d9274352SBarry Smith S*/
1227e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1228639f9d9dSBarry Smith 
1229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1236d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *);
1237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1238edaa7c33SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,const PetscInt*[]);
1239f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1240f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1241f86b9fbaSHong Zhang 
1242b1683b59SHong Zhang 
1243b1683b59SHong Zhang /*S
1244b9af6bddSHong Zhang      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1245b1683b59SHong Zhang 
1246b1683b59SHong Zhang    Level: beginner
1247b1683b59SHong Zhang 
1248b1683b59SHong Zhang   Concepts: coloring, sparse matrix product
1249b1683b59SHong Zhang 
1250b9af6bddSHong Zhang .seealso:  MatTransposeColoringCreate()
1251b1683b59SHong Zhang S*/
1252b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1253b1683b59SHong Zhang 
1254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1258b1683b59SHong Zhang 
1259639f9d9dSBarry Smith /*
12600752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12613eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12620752156aSBarry Smith */
1263ca161407SBarry Smith 
1264d9274352SBarry Smith /*S
1265d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1266d9274352SBarry Smith 
1267d9274352SBarry Smith    Level: beginner
1268d9274352SBarry Smith 
1269d9274352SBarry Smith   Concepts: partitioning
1270d9274352SBarry Smith 
1271743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1272d9274352SBarry Smith S*/
127391e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1274d9274352SBarry Smith 
127576bdecfbSBarry Smith /*J
12768f6c3df8SBarry Smith     MatPartitioningType - String with the name of a PETSc matrix partitioning
1277d9274352SBarry Smith 
1278d9274352SBarry Smith    Level: beginner
12790d04baf8SBarry Smith dm
1280b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
128176bdecfbSBarry Smith J*/
128219fd82e9SBarry Smith typedef const char* MatPartitioningType;
12832692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
1284aa1e27eaSFande Kong #define MATPARTITIONINGAVERAGE   "average"
12852692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
12862692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
12872692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
12882692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
128950d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
129088d2ac2bSFande Kong #define MATPARTITIONINGHIERARCH  "hierarch"
129171306c60Sjroman 
1292ca161407SBarry Smith 
1293014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
129419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1295014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1296014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1297014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1298014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1299014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1300014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
13012aabb6bbSBarry Smith 
1302bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
130330de9b25SBarry Smith 
1304f1144a30SSatish Balay 
13052bad1931SBarry Smith 
1306014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1307014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
130819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1309ca161407SBarry Smith 
131027538973SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1311014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1312014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
13130752156aSBarry Smith 
1314b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
13156a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1316b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
13176a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[];
1318b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
13196a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[];
1320b6956c12SJose Roman 
1321014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1322014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1323014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1324014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1325014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1326014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1327014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1328014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1329014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1330014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1331014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
133271306c60Sjroman 
133371306c60Sjroman #define MP_PARTY_OPT "opt"
133471306c60Sjroman #define MP_PARTY_LIN "lin"
133571306c60Sjroman #define MP_PARTY_SCA "sca"
133671306c60Sjroman #define MP_PARTY_RAN "ran"
133771306c60Sjroman #define MP_PARTY_GBF "gbf"
133871306c60Sjroman #define MP_PARTY_GCF "gcf"
133971306c60Sjroman #define MP_PARTY_BUB "bub"
134071306c60Sjroman #define MP_PARTY_DEF "def"
1341014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
134271306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
134371306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
134471306c60Sjroman #define MP_PARTY_NONE "no"
1345014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1346014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1347014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1348014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
134971306c60Sjroman 
135050d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
13516a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1352e0f1cffaSJose Roman 
1353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
135771306c60Sjroman 
1358b43b03e9SMark F. Adams /*
13596294fa2bSFande Kong  * hierarchical partitioning
13606294fa2bSFande Kong  */
13614e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*);
13624e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*);
13634e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt);
13644e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);
13656294fa2bSFande Kong 
1366014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1367014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1368591294e4SBarry Smith 
13690752156aSBarry Smith /*
1370af0996ceSBarry Smith     If you add entries here you must also add them to petsc/finclude/petscmat.h
1371d4fbbf0eSBarry Smith */
13721c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13731c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13741c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13751c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13761c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13777c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13787c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13791c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13801c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13817c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13827c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13831c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13841c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1385a714c06dSBarry Smith                MATOP_SOR=13,
13861c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13871c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13881c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13891c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13901c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13911c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13921c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13931c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1394d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1395d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1396d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1397d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1398d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1399d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1400d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1401d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1402d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1403d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1404a5b7ff6bSBarry Smith                MATOP_GET_DIAGONAL_BLOCK=32,
14059d39f69dSJed Brown                /* MATOP_PLACEHOLDER_33=33, */
1406d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1407d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1408d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1409d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1410d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1411d519adbfSMatthew Knepley                MATOP_AXPY=39,
14127dae84e0SHong Zhang                MATOP_CREATE_SUBMATRICES=40,
1413d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1414d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1415d519adbfSMatthew Knepley                MATOP_COPY=43,
1416d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1417d519adbfSMatthew Knepley                MATOP_SCALE=45,
1418d519adbfSMatthew Knepley                MATOP_SHIFT=46,
141935153367SBarry Smith                MATOP_DIAGONAL_SET=47,
14209d39f69dSJed Brown                MATOP_ZERO_ROWS_COLUMNS=48,
14219d39f69dSJed Brown                MATOP_SET_RANDOM=49,
1422d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1423d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1424d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1425d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1426d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1427d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1428d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1429d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1430d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
14317dae84e0SHong Zhang                MATOP_CREATE_SUBMATRIX=59,
1432d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1433d519adbfSMatthew Knepley                MATOP_VIEW=61,
1434d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
14357bab7c10SHong Zhang                MATOP_MATMAT_MULT=63,
14367bab7c10SHong Zhang                MATOP_MATMAT_MULT_SYMBOLIC=64,
14377bab7c10SHong Zhang                MATOP_MATMAT_MULT_NUMERIC=65,
1438d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1439d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1440d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1441d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1442d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1443d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1444d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
14459d39f69dSJed Brown                /* MATOP_PLACEHOLDER_73=73, */
1446d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1447d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1448d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1449bda6c4cbSBarry Smith                MATOP_MULT_CONSTRAINED=77,
1450bda6c4cbSBarry Smith                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
14519d39f69dSJed Brown                MATOP_FIND_ZERO_DIAGONALS=79,
1452d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1453d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1454d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1455d519adbfSMatthew Knepley                MATOP_LOAD=83,
1456d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1457d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1458d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1459820469bcSHong Zhang                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1460d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1461d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1462d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1463d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1464d519adbfSMatthew Knepley                MATOP_PTAP=92,
1465d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1466d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1467bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT=95,
1468bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1469bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
14709d39f69dSJed Brown                /* MATOP_PLACEHOLDER_98=98, */
14719d39f69dSJed Brown                /* MATOP_PLACEHOLDER_99=99, */
14729d39f69dSJed Brown                /* MATOP_PLACEHOLDER_100=100, */
14739d39f69dSJed Brown                /* MATOP_PLACEHOLDER_101=101, */
1474d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
14759d39f69dSJed Brown                /* MATOP_PLACEHOLDER_103=103, */
1476d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1477d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1478bda6c4cbSBarry Smith                MATOP_IMAGINARY_PART=106,
1479bda6c4cbSBarry Smith                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1480bda6c4cbSBarry Smith                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1481bda6c4cbSBarry Smith                MATOP_MAT_SOLVE=109,
148229eadf9eSStefano Zampini                MATOP_MAT_SOLVE_TRANSPOSE=110,
1483d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
14840e8d04f7SBarry Smith                MATOP_GET_COLUMN_VECTOR=112,
1485d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
14860e8d04f7SBarry Smith                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
148789c6957cSBarry Smith                MATOP_CREATE=115,
148889c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
14890e8d04f7SBarry Smith                MATOP_GET_LOCAL_SUB_MATRIX=117,
14900e8d04f7SBarry Smith                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1491eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
14920e8d04f7SBarry Smith                MATOP_HERMITIAN_TRANSPOSE=120,
14930e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
14940e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
14950e8d04f7SBarry Smith                MATOP_GET_MULTI_PROC_BLOCK=123,
14969d39f69dSJed Brown                MATOP_FIND_NONZERO_ROWS=124,
14970e8d04f7SBarry Smith                MATOP_GET_COLUMN_NORMS=125,
14989d39f69dSJed Brown                MATOP_INVERT_BLOCK_DIAGONAL=126,
14999d39f69dSJed Brown                /* MATOP_PLACEHOLDER_127=127, */
15000e8d04f7SBarry Smith                MATOP_GET_SUB_MATRICES_PARALLE=128,
15012b8ad9a3SHong Zhang                MATOP_SET_VALUES_BATCH=129,
15020e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT=130,
15030e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
1504e9b602ebSSatish Balay                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
15050e8d04f7SBarry Smith                MATOP_TRANSPOSE_COLORING_CREAT=133,
15060e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_SPT=134,
15070e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_DEN=135,
15080e8d04f7SBarry Smith                MATOP_RART=136,
15090e8d04f7SBarry Smith                MATOP_RART_SYMBOLIC=137,
15100e8d04f7SBarry Smith                MATOP_RART_NUMERIC=138,
1511e09a3074SHong Zhang                MATOP_SET_BLOCK_SIZES=139,
1512f9426fe0SMark Adams                MATOP_AYPX=140,
15131919a2e2SJed Brown                MATOP_RESIDUAL=141,
15149c8f2541SHong Zhang                MATOP_FDCOLORING_SETUP=142,
1515*2d033e1fSHong Zhang                MATOP_MPICONCATENATESEQ=144,
1516*2d033e1fSHong Zhang                MATOP_DESTROYSUBMATRICES=145
1517fae171e0SBarry Smith              } MatOperation;
1518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *);
1519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1522112a2221SBarry Smith 
152390ace30eSBarry Smith /*
152490ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
152590ace30eSBarry Smith    stored in a universal format. By changing the format with
15266a9046bcSBarry Smith    PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
152790ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
152890ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1529f4403165SShri Abhyankar    read into matrices of the same type.
153090ace30eSBarry Smith */
153190ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
153290ace30eSBarry Smith 
1533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
15353b3b1effSJed Brown PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat,Mat*);
1536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
1537b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*);
15381f4e1ec7SBarry Smith 
1539d9274352SBarry Smith /*S
1540d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1541d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1542d9274352SBarry Smith 
1543f7a9e4ceSBarry Smith    Level: advanced
1544d9274352SBarry Smith 
1545d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1546d9274352SBarry Smith 
15476e1639daSBarry Smith   Users manual sections:
15486e1639daSBarry Smith .   sec_singular
15496e1639daSBarry Smith 
1550d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1551d9274352SBarry Smith S*/
155274637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1553d9274352SBarry Smith 
1554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1557d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
15595fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *);
15605fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace);
1561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1565014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1567014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
156874637425SBarry Smith 
1569014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1570014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1571014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1572014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
15733f1d51d7SBarry Smith 
1574014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1575014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1576014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1577c4f061fbSSatish Balay 
1578014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*);
1579b0a32e0cSBarry Smith 
1580014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
158104f1ad80SBarry Smith 
1582014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1583014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1584014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1585014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1586014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1587014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1588014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1589014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1590014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1591014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1592014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1593014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1594014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1595e884886fSBarry Smith 
15966370053bSBarry Smith /*S
15976370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
15986370053bSBarry Smith               Jacobian vector products
1599e884886fSBarry Smith 
16006370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
16016370053bSBarry Smith 
16026370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
16036370053bSBarry Smith 
16046370053bSBarry Smith     Level: developer
16056370053bSBarry Smith 
16066370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
16076370053bSBarry Smith S*/
1608e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1609e884886fSBarry Smith 
161076bdecfbSBarry Smith /*J
1611e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1612e884886fSBarry Smith 
1613e884886fSBarry Smith    Level: beginner
1614e884886fSBarry Smith 
1615e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
161676bdecfbSBarry Smith J*/
161719fd82e9SBarry Smith typedef const char* MatMFFDType;
1618e884886fSBarry Smith #define MATMFFD_DS  "ds"
1619e884886fSBarry Smith #define MATMFFD_WP  "wp"
1620e884886fSBarry Smith 
162119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1622bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));
1623e884886fSBarry Smith 
1624014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1625014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool );
1626e884886fSBarry Smith 
162761ab5df0SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType);
162861ab5df0SBarry Smith 
1629014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1630014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
16317dbadf16SMatthew Knepley 
163297969023SHong Zhang /*
163397969023SHong Zhang    PETSc interface to MUMPS
163497969023SHong Zhang */
163597969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1636014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1637bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*);
1638b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
1639bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*);
1640bc6112feSHong Zhang 
1641ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1642ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1643ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1644ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
164597969023SHong Zhang #endif
164623a5497aSJed Brown 
1647d954db57SHong Zhang /*
1648d615f992SSatish Balay    PETSc interface to Mkl_Pardiso
1649b500a6b7SJose David Bermeol */
1650b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO
1651d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1652b500a6b7SJose David Bermeol #endif
1653b500a6b7SJose David Bermeol 
1654b500a6b7SJose David Bermeol /*
1655d305a81bSVasiliy Kozyrev    PETSc interface to Mkl_CPardiso
1656d305a81bSVasiliy Kozyrev */
1657d305a81bSVasiliy Kozyrev #ifdef PETSC_HAVE_MKL_CPARDISO
1658d305a81bSVasiliy Kozyrev PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt);
1659d305a81bSVasiliy Kozyrev #endif
1660d305a81bSVasiliy Kozyrev 
1661d305a81bSVasiliy Kozyrev /*
1662d954db57SHong Zhang    PETSc interface to SUPERLU
1663d954db57SHong Zhang */
1664d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1665014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1666d954db57SHong Zhang #endif
1667d954db57SHong Zhang 
1668fb785984SHong Zhang /*
1669fb785984SHong Zhang    PETSc interface to SUPERLU_DIST
1670fb785984SHong Zhang */
1671fb785984SHong Zhang #ifdef PETSC_HAVE_SUPERLU_DIST
1672fb785984SHong Zhang PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*);
1673fb785984SHong Zhang #endif
1674fb785984SHong Zhang 
1675575704cbSPieter Ghysels /*
1676575704cbSPieter Ghysels    PETSc interface to STRUMPACK
1677575704cbSPieter Ghysels */
1678575704cbSPieter Ghysels #ifdef PETSC_HAVE_STRUMPACK
1679575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat,PetscReal);
1680575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat,PetscInt);
1681575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool);
1682575704cbSPieter Ghysels #endif
1683575704cbSPieter Ghysels 
1684575704cbSPieter Ghysels 
1685aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA
16863f9c0db1SPaul Mullowney /*E
1687e057df02SPaul Mullowney     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
16882692e278SPaul Mullowney     matrices.
1689e057df02SPaul Mullowney 
1690e057df02SPaul Mullowney     Not Collective
1691e057df02SPaul Mullowney 
16928468deeeSKarl Rupp +   MAT_CUSPARSE_CSR - Compressed Sparse Row
16932692e278SPaul Mullowney .   MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later).
16942692e278SPaul Mullowney -   MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).
1695e057df02SPaul Mullowney 
1696e057df02SPaul Mullowney     Level: intermediate
1697e057df02SPaul Mullowney 
1698af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1699e057df02SPaul Mullowney 
1700e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1701e057df02SPaul Mullowney E*/
1702e057df02SPaul Mullowney 
17033f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;
1704e057df02SPaul Mullowney 
1705e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1706e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
1707e057df02SPaul Mullowney 
17083f9c0db1SPaul Mullowney /*E
1709e057df02SPaul Mullowney     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
17102692e278SPaul Mullowney     matrices whose operation should use a particular storage format.
1711e057df02SPaul Mullowney 
1712e057df02SPaul Mullowney     Not Collective
1713e057df02SPaul Mullowney 
17148468deeeSKarl Rupp +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
17158468deeeSKarl Rupp .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
17168468deeeSKarl Rupp .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
17178468deeeSKarl Rupp -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices
1718e057df02SPaul Mullowney 
1719e057df02SPaul Mullowney     Level: intermediate
1720e057df02SPaul Mullowney 
1721e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1722e057df02SPaul Mullowney E*/
172336d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;
1724e057df02SPaul Mullowney 
172510e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
172610e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1727e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
17289ae82921SPaul Mullowney #endif
17299ae82921SPaul Mullowney 
173090c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
1731014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1732014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1733e057df02SPaul Mullowney 
17343f9c0db1SPaul Mullowney /*E
1735e057df02SPaul Mullowney     MatCUSPStorageFormat - indicates the storage format for CUSP (GPU)
173636d62e41SPaul Mullowney     matrices.
1737e057df02SPaul Mullowney 
1738e057df02SPaul Mullowney     Not Collective
1739e057df02SPaul Mullowney 
17408468deeeSKarl Rupp +   MAT_CUSP_CSR - Compressed Sparse Row
17418468deeeSKarl Rupp .   MAT_CUSP_DIA - Diagonal
17428468deeeSKarl Rupp -   MAT_CUSP_ELL - Ellpack
1743e057df02SPaul Mullowney 
1744e057df02SPaul Mullowney     Level: intermediate
1745e057df02SPaul Mullowney 
1746af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1747e057df02SPaul Mullowney 
1748e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
1749e057df02SPaul Mullowney E*/
17503f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat;
1751e057df02SPaul Mullowney 
1752e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1753e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[];
1754e057df02SPaul Mullowney 
17553f9c0db1SPaul Mullowney /*E
1756e057df02SPaul Mullowney     MatCUSPFormatOperation - indicates the operation of CUSP (GPU)
175736d62e41SPaul Mullowney     matrices whose operation should use a particular storage format.
1758e057df02SPaul Mullowney 
1759e057df02SPaul Mullowney     Not Collective
1760e057df02SPaul Mullowney 
17618468deeeSKarl Rupp +   MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
17628468deeeSKarl Rupp .   MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
17638468deeeSKarl Rupp .   MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
17648468deeeSKarl Rupp -   MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices
1765e057df02SPaul Mullowney 
1766e057df02SPaul Mullowney     Level: intermediate
1767e057df02SPaul Mullowney 
1768af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1769e057df02SPaul Mullowney 
1770e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat
1771e057df02SPaul Mullowney E*/
177236d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation;
1773e057df02SPaul Mullowney 
1774e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat);
1775e057df02SPaul Mullowney #endif
177690c53902SBarry Smith 
1777d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL)
1778d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1779d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1780d67ff14aSKarl Rupp #endif
1781d67ff14aSKarl Rupp 
178254efbe56SHong Zhang /*
178354efbe56SHong Zhang    PETSc interface to FFTW
178454efbe56SHong Zhang */
178554efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
1786014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1787014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
17882a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*);
178954efbe56SHong Zhang #endif
179054efbe56SHong Zhang 
1791014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1792014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1793014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1794014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1795014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1796014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
179719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1798014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1799014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1800d8588912SDave May 
18014325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
1802e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);
18034325cce7SMatthew G Knepley 
1804b541e6a4SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**);
180553dd7562SDmitry Karpeev 
1806c094ef40SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat);
1807c094ef40SMatthew G. Knepley 
1808539c167fSBarry Smith PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*);
1809539c167fSBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*);
1810539c167fSBarry Smith 
181123a5497aSJed Brown #endif
1812