xref: /petsc/include/petscmat.h (revision c1154cd5aacc273eaa4d2acbee03e5887dd0cfdb)
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"
84773941d6SBarry Smith 
8576bdecfbSBarry Smith /*J
86c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
879989ab13SBarry Smith 
88c2b89b5dSBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu_dist" the parallel SuperLU_DIST package etc
899989ab13SBarry Smith 
909989ab13SBarry Smith    Level: beginner
919989ab13SBarry Smith 
925c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
9376bdecfbSBarry Smith J*/
94c7393fdbSBarry Smith #define MatSolverPackage char*
952692d6eeSBarry Smith #define MATSOLVERSUPERLU      "superlu"
962692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist"
9708f5efcfSPieter Ghysels #define MATSOLVERSTRUMPACK    "strumpack"
982692d6eeSBarry Smith #define MATSOLVERUMFPACK      "umfpack"
9920db9a53SJed Brown #define MATSOLVERCHOLMOD      "cholmod"
100d89f5e7aSBarry Smith #define MATSOLVERKLU          "klu"
101ecddaf3cSBarry Smith #define MATSOLVERCLIQUE       "clique"
102d89f5e7aSBarry Smith #define MATSOLVERELEMENTAL    "elemental"
1032692d6eeSBarry Smith #define MATSOLVERESSL         "essl"
1042692d6eeSBarry Smith #define MATSOLVERLUSOL        "lusol"
1052692d6eeSBarry Smith #define MATSOLVERMUMPS        "mumps"
106d615f992SSatish Balay #define MATSOLVERMKL_PARDISO  "mkl_pardiso"
107d305a81bSVasiliy Kozyrev #define MATSOLVERMKL_CPARDISO "mkl_cpardiso"
1082692d6eeSBarry Smith #define MATSOLVERPASTIX       "pastix"
1092692d6eeSBarry Smith #define MATSOLVERMATLAB       "matlab"
1102692d6eeSBarry Smith #define MATSOLVERPETSC        "petsc"
1112692d6eeSBarry Smith #define MATSOLVERBAS          "bas"
1129ae82921SPaul Mullowney #define MATSOLVERCUSPARSE     "cusparse"
113c0cdd4a1SDahai Guo 
114b24902e0SBarry Smith /*E
115b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
116b24902e0SBarry Smith 
117b24902e0SBarry Smith     Level: beginner
118b24902e0SBarry Smith 
119af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
120b24902e0SBarry Smith 
121c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
122b24902e0SBarry Smith E*/
123599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
124014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[];
125e92e720dSBarry Smith 
126014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
127014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
128014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
129014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*);
13018713533SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister(const MatSolverPackage,const MatType,MatFactorType,PetscErrorCode(*)(Mat,MatFactorType,Mat*));
13142c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageGet(const MatSolverPackage,const MatType,MatFactorType,PetscBool*,PetscBool*,PetscErrorCode (**)(Mat,MatFactorType,Mat*));
1329989ab13SBarry Smith 
133c06d978dSMatthew Knepley /* Logging support */
1340700a824SBarry Smith #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
135014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID;
136335efc43SPeter Brune PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID;
137014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID;
138014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
139014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID;
140014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID;
141014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID;
142014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID;
143c06d978dSMatthew Knepley 
144ceb03754SKris Buschelman /*E
145ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
146d6eff37eSBarry Smith      or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate
147d6eff37eSBarry Smith      that the input matrix is to be replaced with the converted matrix.
148ceb03754SKris Buschelman 
149ceb03754SKris Buschelman     Level: beginner
150ceb03754SKris Buschelman 
151af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
152ceb03754SKris Buschelman 
1530c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
154ceb03754SKris Buschelman E*/
155511c6705SHong Zhang typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse;
156ceb03754SKris Buschelman 
1575494a064SHong Zhang /*E
1585494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
159829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1605494a064SHong Zhang 
1615494a064SHong Zhang     Level: beginner
1625494a064SHong Zhang 
163829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1645494a064SHong Zhang E*/
1655494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1665494a064SHong Zhang 
167607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void);
168c06d978dSMatthew Knepley 
169014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*);
170014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
17119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType);
172014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat);
173685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
174bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat));
175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]);
176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]);
177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]);
178014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]);
17984d44b13SHong Zhang PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat,PetscBool);
180f69a0ea3SMatthew Knepley 
181140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList;
182140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList;
183140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList;
18428988994SBarry Smith 
1853b224e63SBarry Smith /*E
186aa6c7ce3SBarry Smith     MatStructure - Indicates if two matrices have the same nonzero structure
1873b224e63SBarry Smith 
1883b224e63SBarry Smith     Level: beginner
1893b224e63SBarry Smith 
190af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1913b224e63SBarry Smith 
192aa6c7ce3SBarry Smith .seealso: MatCopy(), MatAXPY()
1933b224e63SBarry Smith E*/
194aa6c7ce3SBarry Smith typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure;
1953b224e63SBarry Smith 
196014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
197014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
198014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
199014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
200014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
201014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2028d7a6e47SBarry Smith 
203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
206d21a29f3SJed Brown 
207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
208014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
209d21a29f3SJed Brown 
210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
21238f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
2144cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]);
215dfb205c3SBarry Smith 
216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
218c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*);
219986c4d72SJose E. Roman PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Vec,Mat,Mat*);
220267ca84cSJose E. Roman PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat,Mat*,Mat*,Vec*,Mat*);
221c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*);
222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
224c0cdd4a1SDahai Guo 
225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2326d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
2346d7c1e57SBarry Smith 
23519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*);
236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
237dedccee8SHong Zhang 
238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
2398060fb66Sstefano_zampini PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat,Mat*);
240d0de2241SAndrew Spott PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*);
241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*);
242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS);
243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
2441472f72bSBarry Smith 
245978814f1SStefano Zampini #if defined(PETSC_HAVE_HYPRE)
246d975228cSstefano_zampini PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
247978814f1SStefano Zampini #endif
248978814f1SStefano Zampini 
249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]);
2501d6018f0SLisandro Dalcin 
251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
253e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*);
25421c89e3eSBarry Smith 
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
260713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
26199cafbc1SBarry Smith 
2628ed539a5SBarry Smith /* ------------------------------------------------------------*/
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
267014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
26873a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);
26984cb2905SBarry Smith 
2702ef4de8bSBarry Smith /*S
2712ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
27262ef0227SBarry Smith         column of a matrix as indexed on an associated grid. These are arguments to MatSetStencil() and MatSetBlockStencil()
27362ef0227SBarry Smith 
27462ef0227SBarry 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).
27562ef0227SBarry 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.
27662ef0227SBarry Smith 
27762ef0227SBarry Smith    For stencil access to vectors see DMDAVecGetArray(), DMDAVecGetArrayF90().
278be479b30SJed Brown 
279be479b30SJed Brown    Fortran usage is different, see MatSetValuesStencil() for details.
2802ef4de8bSBarry Smith 
2812ef4de8bSBarry Smith    Level: beginner
2822ef4de8bSBarry Smith 
2832ef4de8bSBarry Smith   Concepts: matrix; linear operator
2842ef4de8bSBarry Smith 
28562ef0227SBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90()
2862ef4de8bSBarry Smith S*/
287435da068SBarry Smith typedef struct {
288c1ac3661SBarry Smith   PetscInt k,j,i,c;
289435da068SBarry Smith } MatStencil;
2902ef4de8bSBarry Smith 
291014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
292014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
293014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
294435da068SBarry Smith 
295d91e6319SBarry Smith /*E
296d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
29762ef0227SBarry Smith      to continue to add or insert values to it
298d91e6319SBarry Smith 
299d91e6319SBarry Smith     Level: beginner
300d91e6319SBarry Smith 
301d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
302d91e6319SBarry Smith E*/
3036d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
304014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
305014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
306014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);
3074f9c727eSBarry Smith 
3081d73ed98SBarry Smith 
30930de9b25SBarry Smith 
310d91e6319SBarry Smith /*E
311d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
312d91e6319SBarry Smith 
313d91e6319SBarry Smith     Level: beginner
314d91e6319SBarry Smith 
315af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
316d91e6319SBarry Smith 
317382164b0SBarry Smith    Developer Notes: Entries that are negative need not be called collectively by all processes.
318382164b0SBarry Smith 
319d91e6319SBarry Smith .seealso: MatSetOption()
320d91e6319SBarry Smith E*/
321c5e4d11fSDmitry Karpeev typedef enum {MAT_OPTION_MIN = -3,
322c5e4d11fSDmitry Karpeev               MAT_UNUSED_NONZERO_LOCATION_ERR = -2,
32392d4d306SBarry Smith               MAT_ROW_ORIENTED = -1,
324382164b0SBarry Smith               MAT_SYMMETRIC = 1,
325382164b0SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC = 2,
326382164b0SBarry Smith               MAT_NEW_DIAGONALS = 3,
327382164b0SBarry Smith               MAT_IGNORE_OFF_PROC_ENTRIES = 4,
328382164b0SBarry Smith               MAT_USE_HASH_TABLE = 5,
329382164b0SBarry Smith               MAT_KEEP_NONZERO_PATTERN = 6,
330382164b0SBarry Smith               MAT_IGNORE_ZERO_ENTRIES = 7,
331382164b0SBarry Smith               MAT_USE_INODES = 8,
332382164b0SBarry Smith               MAT_HERMITIAN = 9,
333382164b0SBarry Smith               MAT_SYMMETRY_ETERNAL = 10,
334c5e4d11fSDmitry Karpeev               MAT_NEW_NONZERO_LOCATION_ERR = 11,
335382164b0SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR = 12,
336382164b0SBarry Smith               MAT_ERROR_LOWER_TRIANGULAR = 13,
337382164b0SBarry Smith               MAT_GETROW_UPPERTRIANGULAR = 14,
338382164b0SBarry Smith               MAT_SPD = 15,
33992d4d306SBarry Smith               MAT_NO_OFF_PROC_ZERO_ROWS = 16,
34092d4d306SBarry Smith               MAT_NO_OFF_PROC_ENTRIES = 17,
34192d4d306SBarry Smith               MAT_NEW_NONZERO_LOCATIONS = 18,
342c5e4d11fSDmitry Karpeev               MAT_NEW_NONZERO_ALLOCATION_ERR = 19,
343c5e4d11fSDmitry Karpeev               MAT_SUBSET_OFF_PROC_ENTRIES = 20,
344c10200c1SHong Zhang               MAT_SUBMAT_SINGLEIS = 21,
345c10200c1SHong Zhang               MAT_OPTION_MAX = 22} MatOption;
346e057df02SPaul Mullowney 
347014dd563SJed Brown PETSC_EXTERN const char *MatOptions[];
348014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool);
3497d68702bSBarry Smith PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*);
35019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);
35184cb2905SBarry Smith 
352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
3588c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
3598c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
36021e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*);
361bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3628c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
3638c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
364014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
365014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
366014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
367014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
36833d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat);
3691620fd73SBarry Smith 
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
382f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);
383f5edf698SHong Zhang 
384d91e6319SBarry Smith /*E
385d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
386d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
387d91e6319SBarry Smith 
388d91e6319SBarry Smith     Level: beginner
389d91e6319SBarry Smith 
390af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
391d91e6319SBarry Smith 
39270dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
39370dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
39470dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
39570dcbbb9SBarry Smith 
396d91e6319SBarry Smith .seealso: MatDuplicate()
397d91e6319SBarry Smith E*/
39870dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
3992e8a6d31SBarry Smith 
40019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
401014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);
40294a9d846SBarry Smith 
403cb5b572fSBarry Smith 
404014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
405014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
406014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
407014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
409014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
410014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
411014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
4137b80b807SBarry Smith 
4141a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4151a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
4161a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4171a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
418d4fbbf0eSBarry Smith 
419d91e6319SBarry Smith /*S
420d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
421d91e6319SBarry Smith 
422d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
423d91e6319SBarry Smith 
424d91e6319SBarry Smith    Level: intermediate
425d91e6319SBarry Smith 
426d91e6319SBarry Smith   Concepts: matrix^nonzero information
427d91e6319SBarry Smith 
428d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
429d91e6319SBarry Smith S*/
4304e220ebcSLois Curfman McInnes typedef struct {
431b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
432b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
433b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
434b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
435b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
436b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
437b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4384e220ebcSLois Curfman McInnes } MatInfo;
4394e220ebcSLois Curfman McInnes 
440d9274352SBarry Smith /*E
441d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
442d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
443d9274352SBarry Smith 
444d9274352SBarry Smith     Level: beginner
445d9274352SBarry Smith 
446af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
447d9274352SBarry Smith 
448d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
449d9274352SBarry Smith E*/
4507b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
452014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*);
461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);
463a52676ecSHong Zhang 
464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*);
465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*);
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*);
467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*);
468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*);
469a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
470a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
4717b80b807SBarry Smith 
472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*);
473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*);
474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
4817b80b807SBarry Smith 
482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
488566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);
4895ef9f2a5SBarry Smith 
490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
49153dd7562SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);
4987b80b807SBarry Smith 
499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5068efafbd8SBarry Smith 
507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
508aa1e27eaSFande Kong PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov);
509d2b2a242SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool);
5107b80b807SBarry Smith 
511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat);
51422440eb1SKris Buschelman 
5157bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
5167bab7c10SHong Zhang 
517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat);
520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat);
52322440eb1SKris Buschelman 
524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat);
527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat);
530bc011b1eSHong Zhang 
531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);
5331c741599SHong Zhang 
534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);
5367b80b807SBarry Smith 
537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
539a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*);
540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
546052efed2SBarry Smith 
547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
54990f02eecSBarry Smith 
550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
5532a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*);
55484cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);}
55553cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
5583a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
5599c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
560bd481603SBarry Smith 
561bd481603SBarry Smith /*MC
5621620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5631620fd73SBarry Smith 
5641620fd73SBarry Smith    Not collective
5651620fd73SBarry Smith 
566a9834a7dSSatish Balay    Synopsis:
567a9834a7dSSatish Balay      #include <petscmat.h>
568a9834a7dSSatish Balay      PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode)
569a9834a7dSSatish Balay 
5701620fd73SBarry Smith    Input Parameters:
5711620fd73SBarry Smith +  m - the matrix
5721620fd73SBarry Smith .  row - the row location of the entry
5731620fd73SBarry Smith .  col - the column location of the entry
5741620fd73SBarry Smith .  value - the value to insert
5751620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5761620fd73SBarry Smith 
5771620fd73SBarry Smith    Notes:
5781620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
5791620fd73SBarry Smith    values simultaneously if possible.
5801620fd73SBarry Smith 
5811620fd73SBarry Smith    Level: beginner
5821620fd73SBarry Smith 
5831620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
5841620fd73SBarry Smith M*/
5851620fd73SBarry 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);}
5861620fd73SBarry Smith 
5871620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
5881620fd73SBarry Smith 
5891620fd73SBarry 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);}
5901620fd73SBarry Smith 
5911620fd73SBarry Smith /*MC
5920d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
593bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
594bd481603SBarry Smith 
595bd481603SBarry Smith    Synopsis:
596aaa7dc30SBarry Smith    #include <petscmat.h>
597c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
598bd481603SBarry Smith 
599bd481603SBarry Smith    Collective on MPI_Comm
600bd481603SBarry Smith 
601bd481603SBarry Smith    Input Parameters:
602bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
603859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
604859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
605bd481603SBarry Smith 
606bd481603SBarry Smith    Output Parameters:
607bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
608bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
609bd481603SBarry Smith 
610bd481603SBarry Smith    Level: intermediate
611bd481603SBarry Smith 
612bd481603SBarry Smith    Notes:
613a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
614bd481603SBarry Smith 
6151d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
616bd481603SBarry Smith 
617bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
618bd481603SBarry Smith 
6191620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6201620fd73SBarry Smith 
621bd481603SBarry Smith   Concepts: preallocation^Matrix
622bd481603SBarry Smith 
623d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
624d6e23781SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock()
625bd481603SBarry Smith M*/
626c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6277c922b88SBarry Smith { \
628a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
6291795a4d1SJed Brown   _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \
6301795a4d1SJed Brown   __start = 0; __end = __start;                                         \
631c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
632a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6337c922b88SBarry Smith 
634bd481603SBarry Smith /*MC
635bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
636bd481603SBarry Smith        inserted using a local number of the rows and columns
637bd481603SBarry Smith 
638bd481603SBarry Smith    Synopsis:
639aaa7dc30SBarry Smith    #include <petscmat.h>
640c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
641bd481603SBarry Smith 
642bd481603SBarry Smith    Not Collective
643bd481603SBarry Smith 
644bd481603SBarry Smith    Input Parameters:
645784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
646bd481603SBarry Smith .  nrows - the number of rows indicated
6471d73ed98SBarry Smith .  rows - the indices of the rows
648784ac674SJed Brown .  cmap - the column mapping from local to global numbering
649bd481603SBarry Smith .  ncols - the number of columns in the matrix
650bd481603SBarry Smith .  cols - the columns indicated
651bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
652bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
653bd481603SBarry Smith 
654bd481603SBarry Smith    Level: intermediate
655bd481603SBarry Smith 
656bd481603SBarry Smith    Notes:
657a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
658bd481603SBarry Smith 
6591d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
660bd481603SBarry Smith 
661bd481603SBarry Smith   Concepts: preallocation^Matrix
662bd481603SBarry Smith 
663d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
664*c1154cd5SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups()
665bd481603SBarry Smith M*/
666784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
667c4f061fbSSatish Balay {\
668c1ac3661SBarry Smith   PetscInt __l;\
669784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
670784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
671c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
672ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
673c4f061fbSSatish Balay   }\
674c4f061fbSSatish Balay }
675c4f061fbSSatish Balay 
676bd481603SBarry Smith /*MC
677*c1154cd5SBarry Smith    MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be
678*c1154cd5SBarry Smith        inserted using a local number of the rows and columns. This version removes any duplicate columns in cols
679*c1154cd5SBarry Smith 
680*c1154cd5SBarry Smith    Synopsis:
681*c1154cd5SBarry Smith    #include <petscmat.h>
682*c1154cd5SBarry Smith    PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
683*c1154cd5SBarry Smith 
684*c1154cd5SBarry Smith    Not Collective
685*c1154cd5SBarry Smith 
686*c1154cd5SBarry Smith    Input Parameters:
687*c1154cd5SBarry Smith +  map - the row mapping from local numbering to global numbering
688*c1154cd5SBarry Smith .  nrows - the number of rows indicated
689*c1154cd5SBarry Smith .  rows - the indices of the rows (these values are mapped to the global values)
690*c1154cd5SBarry Smith .  cmap - the column mapping from local to global numbering
691*c1154cd5SBarry Smith .  ncols - the number of columns in the matrix   (this value will be changed if duplicate columns are found)
692*c1154cd5SBarry Smith .  cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed)
693*c1154cd5SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
694*c1154cd5SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
695*c1154cd5SBarry Smith 
696*c1154cd5SBarry Smith    Level: intermediate
697*c1154cd5SBarry Smith 
698*c1154cd5SBarry Smith    Notes:
699*c1154cd5SBarry Smith     See Users-Manual: ch_performance for more details.
700*c1154cd5SBarry Smith 
701*c1154cd5SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
702*c1154cd5SBarry Smith 
703*c1154cd5SBarry Smith   Concepts: preallocation^Matrix
704*c1154cd5SBarry Smith 
705*c1154cd5SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
706*c1154cd5SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
707*c1154cd5SBarry Smith M*/
708*c1154cd5SBarry Smith #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
709*c1154cd5SBarry Smith {\
710*c1154cd5SBarry Smith   PetscInt __l;\
711*c1154cd5SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
712*c1154cd5SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
713*c1154cd5SBarry Smith   _4_ierr = PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\
714*c1154cd5SBarry Smith   for (__l=0;__l<nrows;__l++) {\
715*c1154cd5SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
716*c1154cd5SBarry Smith   }\
717*c1154cd5SBarry Smith }
718*c1154cd5SBarry Smith 
719*c1154cd5SBarry Smith /*MC
720d6e23781SBarry Smith    MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
721bd481603SBarry Smith        inserted using a local number of the rows and columns
722bd481603SBarry Smith 
723bd481603SBarry Smith    Synopsis:
724aaa7dc30SBarry Smith    #include <petscmat.h>
725d6e23781SBarry Smith    PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
726d6e23781SBarry Smith 
727d6e23781SBarry Smith    Not Collective
728d6e23781SBarry Smith 
729d6e23781SBarry Smith    Input Parameters:
730d6e23781SBarry Smith +  map - the row mapping from local numbering to global numbering
731d6e23781SBarry Smith .  nrows - the number of rows indicated
732d6e23781SBarry Smith .  rows - the indices of the rows
733d6e23781SBarry Smith .  cmap - the column mapping from local to global numbering
734d6e23781SBarry Smith .  ncols - the number of columns in the matrix
735d6e23781SBarry Smith .  cols - the columns indicated
736d6e23781SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
737d6e23781SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
738d6e23781SBarry Smith 
739d6e23781SBarry Smith    Level: intermediate
740d6e23781SBarry Smith 
741d6e23781SBarry Smith    Notes:
742d6e23781SBarry Smith     See Users-Manual: ch_performance for more details.
743d6e23781SBarry Smith 
744d6e23781SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
745d6e23781SBarry Smith 
746d6e23781SBarry Smith   Concepts: preallocation^Matrix
747d6e23781SBarry Smith 
748d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
749d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
750d6e23781SBarry Smith M*/
751d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
752d6e23781SBarry Smith {\
753d6e23781SBarry Smith   PetscInt __l;\
754d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
755d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
756d6e23781SBarry Smith   for (__l=0;__l<nrows;__l++) {\
757d6e23781SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
758d6e23781SBarry Smith   }\
759d6e23781SBarry Smith }
760d6e23781SBarry Smith 
761d6e23781SBarry Smith /*MC
762d6e23781SBarry Smith    MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
763d6e23781SBarry Smith        inserted using a local number of the rows and columns
764d6e23781SBarry Smith 
765d6e23781SBarry Smith    Synopsis:
766d6e23781SBarry Smith    #include <petscmat.h>
767d6e23781SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
768bd481603SBarry Smith 
769bd481603SBarry Smith    Not Collective
770bd481603SBarry Smith 
771bd481603SBarry Smith    Input Parameters:
772bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
773bd481603SBarry Smith .  nrows - the number of rows indicated
7741d73ed98SBarry Smith .  rows - the indices of the rows
775bd481603SBarry Smith .  ncols - the number of columns in the matrix
776bd481603SBarry Smith .  cols - the columns indicated
777bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
778bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
779bd481603SBarry Smith 
780bd481603SBarry Smith    Level: intermediate
781bd481603SBarry Smith 
782bd481603SBarry Smith    Notes:
783a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
784bd481603SBarry Smith 
785bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
786bd481603SBarry Smith 
787bd481603SBarry Smith   Concepts: preallocation^Matrix
788bd481603SBarry Smith 
789d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(),
790d6e23781SBarry Smith           MatPreallocateInitialize(),  MatPreallocateSetLocal()
791bd481603SBarry Smith M*/
792d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\
793d3d32019SBarry Smith {\
794c1ac3661SBarry Smith   PetscInt __l;\
795d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
796d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
797d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
798d6e23781SBarry Smith     _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
799d3d32019SBarry Smith   }\
800d3d32019SBarry Smith }
801bd481603SBarry Smith /*MC
802bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
803bd481603SBarry Smith        inserted using a local number of the rows and columns
804bd481603SBarry Smith 
805bd481603SBarry Smith    Synopsis:
806aaa7dc30SBarry Smith    #include <petscmat.h>
807c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
808bd481603SBarry Smith 
809bd481603SBarry Smith    Not Collective
810bd481603SBarry Smith 
811bd481603SBarry Smith    Input Parameters:
81264075487SBarry Smith +  row - the row
813bd481603SBarry Smith .  ncols - the number of columns in the matrix
814a50f8bf6SHong Zhang -  cols - the columns indicated
815a50f8bf6SHong Zhang 
816a50f8bf6SHong Zhang    Output Parameters:
817a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
818bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
819bd481603SBarry Smith 
820bd481603SBarry Smith    Level: intermediate
821bd481603SBarry Smith 
822bd481603SBarry Smith    Notes:
823a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
824bd481603SBarry Smith 
825bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
826bd481603SBarry Smith 
8271620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8281620fd73SBarry Smith 
829bd481603SBarry Smith   Concepts: preallocation^Matrix
830bd481603SBarry Smith 
831d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
832d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSetLocal()
833bd481603SBarry Smith M*/
834c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
835c1ac3661SBarry Smith { PetscInt __i; \
836e32f2f54SBarry 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);\
837e32f2f54SBarry 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);\
8387c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
83964075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8407cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8417c922b88SBarry Smith   }\
8427c922b88SBarry Smith }
8437c922b88SBarry Smith 
844bd481603SBarry Smith /*MC
845d6e23781SBarry Smith    MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
846bd481603SBarry Smith        inserted using a local number of the rows and columns
847bd481603SBarry Smith 
848bd481603SBarry Smith    Synopsis:
849aaa7dc30SBarry Smith    #include <petscmat.h>
850d6e23781SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
851bd481603SBarry Smith 
852bd481603SBarry Smith    Not Collective
853bd481603SBarry Smith 
854bd481603SBarry Smith    Input Parameters:
855bd481603SBarry Smith +  nrows - the number of rows indicated
8561d73ed98SBarry Smith .  rows - the indices of the rows
857bd481603SBarry Smith .  ncols - the number of columns in the matrix
858bd481603SBarry Smith .  cols - the columns indicated
859bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
860bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
861bd481603SBarry Smith 
862bd481603SBarry Smith    Level: intermediate
863bd481603SBarry Smith 
864bd481603SBarry Smith    Notes:
865a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
866bd481603SBarry Smith 
867bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
868bd481603SBarry Smith 
8691620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8701620fd73SBarry Smith 
871bd481603SBarry Smith   Concepts: preallocation^Matrix
872bd481603SBarry Smith 
873d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(),  MatPreallocateInitialize(),
874d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
875bd481603SBarry Smith M*/
876d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\
877c1ac3661SBarry Smith { PetscInt __i; \
878d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
879d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
880d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
881d3d32019SBarry Smith   }\
882d3d32019SBarry Smith }
883d3d32019SBarry Smith 
884bd481603SBarry Smith /*MC
88516371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
88616371a99SBarry Smith 
88716371a99SBarry Smith    Synopsis:
888aaa7dc30SBarry Smith    #include <petscmat.h>
88916371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
89016371a99SBarry Smith 
89116371a99SBarry Smith    Not Collective
89216371a99SBarry Smith 
89316371a99SBarry Smith    Input Parameters:
89416371a99SBarry Smith .  A - matrix
89516371a99SBarry Smith .  row - row where values exist (must be local to this process)
89616371a99SBarry Smith .  ncols - number of columns
89716371a99SBarry Smith .  cols - columns with nonzeros
89816371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
89916371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
90016371a99SBarry Smith 
90116371a99SBarry Smith    Level: intermediate
90216371a99SBarry Smith 
90316371a99SBarry Smith    Notes:
904a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
90516371a99SBarry Smith 
90616371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
90716371a99SBarry Smith 
90816371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
90916371a99SBarry Smith 
91016371a99SBarry Smith   Concepts: preallocation^Matrix
91116371a99SBarry Smith 
912d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
913d6e23781SBarry Smith           MatPreallocateSymmetricSetLocalBlock()
91416371a99SBarry Smith M*/
9150298fd71SBarry 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);}
91616371a99SBarry Smith 
91716371a99SBarry Smith 
91816371a99SBarry Smith /*MC
9190d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
920bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
921bd481603SBarry Smith 
922bd481603SBarry Smith    Synopsis:
923aaa7dc30SBarry Smith    #include <petscmat.h>
924c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
925bd481603SBarry Smith 
926bd481603SBarry Smith    Collective on MPI_Comm
927bd481603SBarry Smith 
928bd481603SBarry Smith    Input Parameters:
92916371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
930bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
931bd481603SBarry Smith 
932bd481603SBarry Smith    Level: intermediate
933bd481603SBarry Smith 
934bd481603SBarry Smith    Notes:
935a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
936bd481603SBarry Smith 
937bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
938bd481603SBarry Smith 
9391620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9401620fd73SBarry Smith 
941bd481603SBarry Smith   Concepts: preallocation^Matrix
942bd481603SBarry Smith 
943d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
944d6e23781SBarry Smith           MatPreallocateSymmetricSetLocalBlock()
945bd481603SBarry Smith M*/
946a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9477c922b88SBarry Smith 
9487b80b807SBarry Smith /* Routines unique to particular data structures */
949014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
950698d4c6aSKris Buschelman 
951014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
952014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
953698d4c6aSKris Buschelman 
954014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
955014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
956014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
957014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
958014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
959014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
9607b80b807SBarry Smith 
961a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
962a96a251dSBarry Smith 
963014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
964014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
965014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
966273d9f13SBarry Smith 
967014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
968014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
969014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
970014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
971014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
972014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
973014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
974014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
975014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
976014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
9779230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
9789230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
979014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
980273d9f13SBarry Smith 
9812e1947a5SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
982014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt);
983014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
9841b807ce4Svictorle 
985014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
986014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
9872e8a6d31SBarry Smith 
988014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
9893a7fca6bSBarry Smith 
990014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
9917b80b807SBarry Smith /*
9927b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
99394b7f48cSBarry Smith   done through the KSP and PC interfaces.
9947b80b807SBarry Smith */
9957b80b807SBarry Smith 
99676bdecfbSBarry Smith /*J
9978f6c3df8SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering
998d9274352SBarry Smith 
999d9274352SBarry Smith    Level: beginner
1000d9274352SBarry Smith 
1001d9274352SBarry Smith .seealso: MatGetOrdering()
100276bdecfbSBarry Smith J*/
100319fd82e9SBarry Smith typedef const char* MatOrderingType;
10042692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
10052692d6eeSBarry Smith #define MATORDERINGND          "nd"
10062692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
10072692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
10082692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
10092692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
1010510184a7SJed Brown #define MATORDERINGWBM         "wbm"
1011fc1bef75SJed Brown #define MATORDERINGSPECTRAL    "spectral"
1012312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
1013b12f92e5SBarry Smith 
101419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
1015140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
1016bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
1017140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList;
1018d4fbbf0eSBarry Smith 
1019014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1020fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);
1021a2ce50c7SBarry Smith 
1022d91e6319SBarry Smith /*S
1023d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
1024d90ac83dSHong Zhang 
1025d90ac83dSHong Zhang    Level: beginner
1026d90ac83dSHong Zhang 
1027d90ac83dSHong Zhang S*/
1028d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
10296a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[];
10305e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];
1031d90ac83dSHong Zhang 
1032d90ac83dSHong Zhang /*S
10339aa7eafdSHong Zhang     MatFactorError - indicates what type of error in matrix factor
10349aa7eafdSHong Zhang 
10359aa7eafdSHong Zhang     Level: beginner
10369aa7eafdSHong Zhang 
103700e125f8SBarry Smith     Developer Notes: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
103800e125f8SBarry Smith 
103900e125f8SBarry Smith .seealso: MatGetFactor()
10409aa7eafdSHong Zhang S*/
10415cd7cf9dSHong Zhang typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError;
10429aa7eafdSHong Zhang 
104300e125f8SBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*);
1044b8b68cfdSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat);
10457b6c816cSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*);
104600e125f8SBarry Smith 
10479aa7eafdSHong Zhang /*S
1048422a814eSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization
1049b00f7748SHong Zhang 
105061cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
105161cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1052b00f7748SHong Zhang 
105315e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1054b00f7748SHong Zhang 
105561cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
105661cad860SBarry Smith 
1057b00f7748SHong Zhang    Level: developer
1058b00f7748SHong Zhang 
1059d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1060d7d82daaSBarry Smith           MatFactorInfoInitialize()
1061b00f7748SHong Zhang 
1062b00f7748SHong Zhang S*/
1063b00f7748SHong Zhang typedef struct {
106415e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
106585317021SBarry Smith   PetscReal     usedt;
106615e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1067b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
106815e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
106967e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1070348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1071bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1072bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
107315e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1074f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1075f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
107615e8a5b3SHong Zhang } MatFactorInfo;
1077ffa6d0a5SLois Curfman McInnes 
1078014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
10799a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
10809a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
10819a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
10829a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
10839a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
10849a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
10859a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
10869a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
10879a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
10889a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1089014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1090014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1091014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1092014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1093014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1094014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1095014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1096014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
10978e7ba810SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS);
10985a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*);
10995a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*);
11005a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat);
11015a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*);
11025a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec);
11035a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec);
1104014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
1105bb5a7306SBarry Smith 
1106d91e6319SBarry Smith /*E
1107d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1108bb1eb677SSatish Balay 
1109d91e6319SBarry Smith     Level: beginner
1110d91e6319SBarry Smith 
1111d9274352SBarry Smith    May be bitwise ORd together
1112d9274352SBarry Smith 
1113af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1114d91e6319SBarry Smith 
11154e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11164e7234bfSBarry Smith 
111741f059aeSBarry Smith .seealso: MatSOR()
1118d91e6319SBarry Smith E*/
1119ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1120ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1121ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
112284cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1123014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11248ed539a5SBarry Smith 
1125d4fbbf0eSBarry Smith /*
1126639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1127639f9d9dSBarry Smith */
1128b12f92e5SBarry Smith 
11297086a01eSPeter Brune /*S
11307086a01eSPeter Brune      MatColoring - Object for managing the coloring of matrices.
11317086a01eSPeter Brune 
11327086a01eSPeter Brune    Level: beginner
11337086a01eSPeter Brune 
11347086a01eSPeter Brune   Concepts: matrix, coloring
11357086a01eSPeter Brune 
11367086a01eSPeter Brune .seealso:  MatFDColoringCreate() ISColoring MatFDColoring
11377086a01eSPeter Brune S*/
11387086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring;
113976bdecfbSBarry Smith /*J
11408f6c3df8SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring
1141d9274352SBarry Smith 
1142d9274352SBarry Smith    Level: beginner
1143d9274352SBarry Smith 
11447086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring
114576bdecfbSBarry Smith J*/
11467086a01eSPeter Brune 
114719fd82e9SBarry Smith typedef const  char*           MatColoringType;
11485567d71aSPeter Brune #define MATCOLORINGJP      "jp"
1149b9acb617SPeter Brune #define MATCOLORINGPOWER   "power"
11502692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
11512692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
11522692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
11532692d6eeSBarry Smith #define MATCOLORINGID      "id"
11548121bdceSPeter Brune #define MATCOLORINGGREEDY  "greedy"
1155b12f92e5SBarry Smith 
11568ac6da0aSPeter Brune /*E
11578ac6da0aSPeter Brune    MatColoringWeightType - Type of weight scheme
11588ac6da0aSPeter Brune 
11598ac6da0aSPeter Brune     Not Collective
11608ac6da0aSPeter Brune 
11618ac6da0aSPeter Brune +   MAT_COLORING_RANDOM  - Random weights
11628ac6da0aSPeter Brune .   MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
11638ac6da0aSPeter Brune -   MAT_COLORING_LF      - Last-first weighting.
11648ac6da0aSPeter Brune 
11658ac6da0aSPeter Brune     Level: intermediate
11668ac6da0aSPeter Brune 
1167af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
11688ac6da0aSPeter Brune 
11698ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
11708ac6da0aSPeter Brune E*/
1171301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;
11728ac6da0aSPeter Brune 
1173335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1174d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*);
1175335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1176335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1177335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1178335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1179335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1180335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1181335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1182335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1183335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1184335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1185014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
11868ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType);
1187c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*);
11888ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm);
1189cb394dc5SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringTestValid(MatColoring,ISColoring);
1190639f9d9dSBarry Smith 
1191d9274352SBarry Smith /*S
1192d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1193d9274352SBarry Smith         and coloring
1194639f9d9dSBarry Smith 
1195d9274352SBarry Smith    Level: beginner
1196d9274352SBarry Smith 
1197d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1198d9274352SBarry Smith 
1199d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1200d9274352SBarry Smith S*/
1201e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1202639f9d9dSBarry Smith 
1203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1206014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1208014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1210d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *);
1211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1213f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1214f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1215f86b9fbaSHong Zhang 
1216b1683b59SHong Zhang 
1217b1683b59SHong Zhang /*S
1218b9af6bddSHong Zhang      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1219b1683b59SHong Zhang 
1220b1683b59SHong Zhang    Level: beginner
1221b1683b59SHong Zhang 
1222b1683b59SHong Zhang   Concepts: coloring, sparse matrix product
1223b1683b59SHong Zhang 
1224b9af6bddSHong Zhang .seealso:  MatTransposeColoringCreate()
1225b1683b59SHong Zhang S*/
1226b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1227b1683b59SHong Zhang 
1228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1232b1683b59SHong Zhang 
1233639f9d9dSBarry Smith /*
12340752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12353eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12360752156aSBarry Smith */
1237ca161407SBarry Smith 
1238d9274352SBarry Smith /*S
1239d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1240d9274352SBarry Smith 
1241d9274352SBarry Smith    Level: beginner
1242d9274352SBarry Smith 
1243d9274352SBarry Smith   Concepts: partitioning
1244d9274352SBarry Smith 
1245743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1246d9274352SBarry Smith S*/
124791e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1248d9274352SBarry Smith 
124976bdecfbSBarry Smith /*J
12508f6c3df8SBarry Smith     MatPartitioningType - String with the name of a PETSc matrix partitioning
1251d9274352SBarry Smith 
1252d9274352SBarry Smith    Level: beginner
12530d04baf8SBarry Smith dm
1254b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
125576bdecfbSBarry Smith J*/
125619fd82e9SBarry Smith typedef const char* MatPartitioningType;
12572692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
1258aa1e27eaSFande Kong #define MATPARTITIONINGAVERAGE   "average"
12592692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
12602692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
12612692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
12622692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
126350d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
126488d2ac2bSFande Kong #define MATPARTITIONINGHIERARCH  "hierarch"
126571306c60Sjroman 
1266ca161407SBarry Smith 
1267014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
126819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1271014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1272014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1273014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1274014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
12752aabb6bbSBarry Smith 
1276bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
127730de9b25SBarry Smith 
1278f1144a30SSatish Balay 
12792bad1931SBarry Smith 
1280014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1281014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
128219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1283ca161407SBarry Smith 
128427538973SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1285014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1286014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
12870752156aSBarry Smith 
1288b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
12896a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1290b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
12916a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[];
1292b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
12936a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[];
1294b6956c12SJose Roman 
1295014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1296014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1297014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1298014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1299014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1300014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1301014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1302014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1303014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1304014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1305014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
130671306c60Sjroman 
130771306c60Sjroman #define MP_PARTY_OPT "opt"
130871306c60Sjroman #define MP_PARTY_LIN "lin"
130971306c60Sjroman #define MP_PARTY_SCA "sca"
131071306c60Sjroman #define MP_PARTY_RAN "ran"
131171306c60Sjroman #define MP_PARTY_GBF "gbf"
131271306c60Sjroman #define MP_PARTY_GCF "gcf"
131371306c60Sjroman #define MP_PARTY_BUB "bub"
131471306c60Sjroman #define MP_PARTY_DEF "def"
1315014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
131671306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
131771306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
131871306c60Sjroman #define MP_PARTY_NONE "no"
1319014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1320014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1321014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1322014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
132371306c60Sjroman 
132450d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
13256a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1326e0f1cffaSJose Roman 
1327014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1328014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1329014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1330014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
133171306c60Sjroman 
1332b43b03e9SMark F. Adams /*
13336294fa2bSFande Kong  * hierarchical partitioning
13346294fa2bSFande Kong  */
13354e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*);
13364e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*);
13374e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt);
13384e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);
13396294fa2bSFande Kong 
1340014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1341014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1342591294e4SBarry Smith 
13430752156aSBarry Smith /*
1344af0996ceSBarry Smith     If you add entries here you must also add them to petsc/finclude/petscmat.h
1345d4fbbf0eSBarry Smith */
13461c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13471c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13481c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13491c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13501c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13517c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13527c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13531c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13541c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13557c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13567c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13571c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13581c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1359a714c06dSBarry Smith                MATOP_SOR=13,
13601c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13611c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13621c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13631c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13641c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13651c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13661c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13671c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1368d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1369d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1370d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1371d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1372d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1373d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1374d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1375d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1376d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1377d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1378a5b7ff6bSBarry Smith                MATOP_GET_DIAGONAL_BLOCK=32,
13799d39f69dSJed Brown                /* MATOP_PLACEHOLDER_33=33, */
1380d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1381d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1382d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1383d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1384d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1385d519adbfSMatthew Knepley                MATOP_AXPY=39,
1386d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1387d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1388d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1389d519adbfSMatthew Knepley                MATOP_COPY=43,
1390d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1391d519adbfSMatthew Knepley                MATOP_SCALE=45,
1392d519adbfSMatthew Knepley                MATOP_SHIFT=46,
139335153367SBarry Smith                MATOP_DIAGONAL_SET=47,
13949d39f69dSJed Brown                MATOP_ZERO_ROWS_COLUMNS=48,
13959d39f69dSJed Brown                MATOP_SET_RANDOM=49,
1396d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1397d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1398d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1399d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1400d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1401d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1402d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1403d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1404d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1405d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1406d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1407d519adbfSMatthew Knepley                MATOP_VIEW=61,
1408d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
14097bab7c10SHong Zhang                MATOP_MATMAT_MULT=63,
14107bab7c10SHong Zhang                MATOP_MATMAT_MULT_SYMBOLIC=64,
14117bab7c10SHong Zhang                MATOP_MATMAT_MULT_NUMERIC=65,
1412d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1413d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1414d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1415d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1416d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1417d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1418d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
14199d39f69dSJed Brown                /* MATOP_PLACEHOLDER_73=73, */
1420d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1421d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1422d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1423bda6c4cbSBarry Smith                MATOP_MULT_CONSTRAINED=77,
1424bda6c4cbSBarry Smith                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
14259d39f69dSJed Brown                MATOP_FIND_ZERO_DIAGONALS=79,
1426d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1427d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1428d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1429d519adbfSMatthew Knepley                MATOP_LOAD=83,
1430d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1431d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1432d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1433820469bcSHong Zhang                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1434d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1435d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1436d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1437d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1438d519adbfSMatthew Knepley                MATOP_PTAP=92,
1439d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1440d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1441bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT=95,
1442bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1443bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
14449d39f69dSJed Brown                /* MATOP_PLACEHOLDER_98=98, */
14459d39f69dSJed Brown                /* MATOP_PLACEHOLDER_99=99, */
14469d39f69dSJed Brown                /* MATOP_PLACEHOLDER_100=100, */
14479d39f69dSJed Brown                /* MATOP_PLACEHOLDER_101=101, */
1448d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
14499d39f69dSJed Brown                /* MATOP_PLACEHOLDER_103=103, */
1450d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1451d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1452bda6c4cbSBarry Smith                MATOP_IMAGINARY_PART=106,
1453bda6c4cbSBarry Smith                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1454bda6c4cbSBarry Smith                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1455bda6c4cbSBarry Smith                MATOP_MAT_SOLVE=109,
14560e8d04f7SBarry Smith                MATOP_GET_REDUNDANT_MATRIX=110,
1457d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
14580e8d04f7SBarry Smith                MATOP_GET_COLUMN_VECTOR=112,
1459d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
14600e8d04f7SBarry Smith                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
146189c6957cSBarry Smith                MATOP_CREATE=115,
146289c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
14630e8d04f7SBarry Smith                MATOP_GET_LOCAL_SUB_MATRIX=117,
14640e8d04f7SBarry Smith                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1465eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
14660e8d04f7SBarry Smith                MATOP_HERMITIAN_TRANSPOSE=120,
14670e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
14680e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
14690e8d04f7SBarry Smith                MATOP_GET_MULTI_PROC_BLOCK=123,
14709d39f69dSJed Brown                MATOP_FIND_NONZERO_ROWS=124,
14710e8d04f7SBarry Smith                MATOP_GET_COLUMN_NORMS=125,
14729d39f69dSJed Brown                MATOP_INVERT_BLOCK_DIAGONAL=126,
14739d39f69dSJed Brown                /* MATOP_PLACEHOLDER_127=127, */
14740e8d04f7SBarry Smith                MATOP_GET_SUB_MATRICES_PARALLE=128,
14752b8ad9a3SHong Zhang                MATOP_SET_VALUES_BATCH=129,
14760e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT=130,
14770e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
1478e9b602ebSSatish Balay                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
14790e8d04f7SBarry Smith                MATOP_TRANSPOSE_COLORING_CREAT=133,
14800e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_SPT=134,
14810e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_DEN=135,
14820e8d04f7SBarry Smith                MATOP_RART=136,
14830e8d04f7SBarry Smith                MATOP_RART_SYMBOLIC=137,
14840e8d04f7SBarry Smith                MATOP_RART_NUMERIC=138,
1485e09a3074SHong Zhang                MATOP_SET_BLOCK_SIZES=139,
1486f9426fe0SMark Adams                MATOP_AYPX=140,
14871919a2e2SJed Brown                MATOP_RESIDUAL=141,
14889c8f2541SHong Zhang                MATOP_FDCOLORING_SETUP=142,
14899c8f2541SHong Zhang                MATOP_MPICONCATENATESEQ=144
1490fae171e0SBarry Smith              } MatOperation;
1491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *);
1492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1495112a2221SBarry Smith 
149690ace30eSBarry Smith /*
149790ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
149890ace30eSBarry Smith    stored in a universal format. By changing the format with
14996a9046bcSBarry Smith    PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
150090ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
150190ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1502f4403165SShri Abhyankar    read into matrices of the same type.
150390ace30eSBarry Smith */
150490ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
150590ace30eSBarry Smith 
1506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
1509b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*);
15101f4e1ec7SBarry Smith 
1511d9274352SBarry Smith /*S
1512d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1513d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1514d9274352SBarry Smith 
1515f7a9e4ceSBarry Smith    Level: advanced
1516d9274352SBarry Smith 
1517d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1518d9274352SBarry Smith 
15196e1639daSBarry Smith   Users manual sections:
15206e1639daSBarry Smith .   sec_singular
15216e1639daSBarry Smith 
1522d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1523d9274352SBarry Smith S*/
152474637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1525d9274352SBarry Smith 
1526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1529d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
15315fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *);
15325fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace);
1533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
154074637425SBarry Smith 
1541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
15453f1d51d7SBarry Smith 
1546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1549c4f061fbSSatish Balay 
1550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*);
1551b0a32e0cSBarry Smith 
1552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
155304f1ad80SBarry Smith 
1554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1560014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1565014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1567e884886fSBarry Smith 
15686370053bSBarry Smith /*S
15696370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
15706370053bSBarry Smith               Jacobian vector products
1571e884886fSBarry Smith 
15726370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
15736370053bSBarry Smith 
15746370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
15756370053bSBarry Smith 
15766370053bSBarry Smith     Level: developer
15776370053bSBarry Smith 
15786370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
15796370053bSBarry Smith S*/
1580e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1581e884886fSBarry Smith 
158276bdecfbSBarry Smith /*J
1583e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1584e884886fSBarry Smith 
1585e884886fSBarry Smith    Level: beginner
1586e884886fSBarry Smith 
1587e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
158876bdecfbSBarry Smith J*/
158919fd82e9SBarry Smith typedef const char* MatMFFDType;
1590e884886fSBarry Smith #define MATMFFD_DS  "ds"
1591e884886fSBarry Smith #define MATMFFD_WP  "wp"
1592e884886fSBarry Smith 
159319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1594bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));
1595e884886fSBarry Smith 
1596014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1597014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool );
1598e884886fSBarry Smith 
159961ab5df0SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType);
160061ab5df0SBarry Smith 
1601014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1602014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
16037dbadf16SMatthew Knepley 
160497969023SHong Zhang /*
160597969023SHong Zhang    PETSc interface to MUMPS
160697969023SHong Zhang */
160797969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1608014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1609bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*);
1610b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
1611bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*);
1612bc6112feSHong Zhang 
1613ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1614ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1615ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1616ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
161797969023SHong Zhang #endif
161823a5497aSJed Brown 
1619d954db57SHong Zhang /*
1620d615f992SSatish Balay    PETSc interface to Mkl_Pardiso
1621b500a6b7SJose David Bermeol */
1622b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO
1623d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1624b500a6b7SJose David Bermeol #endif
1625b500a6b7SJose David Bermeol 
1626b500a6b7SJose David Bermeol /*
1627d305a81bSVasiliy Kozyrev    PETSc interface to Mkl_CPardiso
1628d305a81bSVasiliy Kozyrev */
1629d305a81bSVasiliy Kozyrev #ifdef PETSC_HAVE_MKL_CPARDISO
1630d305a81bSVasiliy Kozyrev PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt);
1631d305a81bSVasiliy Kozyrev #endif
1632d305a81bSVasiliy Kozyrev 
1633d305a81bSVasiliy Kozyrev /*
1634d954db57SHong Zhang    PETSc interface to SUPERLU
1635d954db57SHong Zhang */
1636d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1637014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1638d954db57SHong Zhang #endif
1639d954db57SHong Zhang 
1640fb785984SHong Zhang /*
1641fb785984SHong Zhang    PETSc interface to SUPERLU_DIST
1642fb785984SHong Zhang */
1643fb785984SHong Zhang #ifdef PETSC_HAVE_SUPERLU_DIST
1644fb785984SHong Zhang PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*);
1645fb785984SHong Zhang #endif
1646fb785984SHong Zhang 
1647575704cbSPieter Ghysels /*
1648575704cbSPieter Ghysels    PETSc interface to STRUMPACK
1649575704cbSPieter Ghysels */
1650575704cbSPieter Ghysels #ifdef PETSC_HAVE_STRUMPACK
1651575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat,PetscReal);
1652575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat,PetscInt);
1653575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool);
1654575704cbSPieter Ghysels #endif
1655575704cbSPieter Ghysels 
1656575704cbSPieter Ghysels 
1657aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA
16583f9c0db1SPaul Mullowney /*E
1659e057df02SPaul Mullowney     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
16602692e278SPaul Mullowney     matrices.
1661e057df02SPaul Mullowney 
1662e057df02SPaul Mullowney     Not Collective
1663e057df02SPaul Mullowney 
16648468deeeSKarl Rupp +   MAT_CUSPARSE_CSR - Compressed Sparse Row
16652692e278SPaul Mullowney .   MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later).
16662692e278SPaul Mullowney -   MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).
1667e057df02SPaul Mullowney 
1668e057df02SPaul Mullowney     Level: intermediate
1669e057df02SPaul Mullowney 
1670af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1671e057df02SPaul Mullowney 
1672e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1673e057df02SPaul Mullowney E*/
1674e057df02SPaul Mullowney 
16753f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;
1676e057df02SPaul Mullowney 
1677e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1678e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
1679e057df02SPaul Mullowney 
16803f9c0db1SPaul Mullowney /*E
1681e057df02SPaul Mullowney     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
16822692e278SPaul Mullowney     matrices whose operation should use a particular storage format.
1683e057df02SPaul Mullowney 
1684e057df02SPaul Mullowney     Not Collective
1685e057df02SPaul Mullowney 
16868468deeeSKarl Rupp +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
16878468deeeSKarl Rupp .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
16888468deeeSKarl Rupp .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
16898468deeeSKarl Rupp -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices
1690e057df02SPaul Mullowney 
1691e057df02SPaul Mullowney     Level: intermediate
1692e057df02SPaul Mullowney 
1693e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1694e057df02SPaul Mullowney E*/
169536d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;
1696e057df02SPaul Mullowney 
169710e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
169810e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1699e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
17009ae82921SPaul Mullowney #endif
17019ae82921SPaul Mullowney 
170290c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
1703014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1704014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1705e057df02SPaul Mullowney 
17063f9c0db1SPaul Mullowney /*E
1707e057df02SPaul Mullowney     MatCUSPStorageFormat - indicates the storage format for CUSP (GPU)
170836d62e41SPaul Mullowney     matrices.
1709e057df02SPaul Mullowney 
1710e057df02SPaul Mullowney     Not Collective
1711e057df02SPaul Mullowney 
17128468deeeSKarl Rupp +   MAT_CUSP_CSR - Compressed Sparse Row
17138468deeeSKarl Rupp .   MAT_CUSP_DIA - Diagonal
17148468deeeSKarl Rupp -   MAT_CUSP_ELL - Ellpack
1715e057df02SPaul Mullowney 
1716e057df02SPaul Mullowney     Level: intermediate
1717e057df02SPaul Mullowney 
1718af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1719e057df02SPaul Mullowney 
1720e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
1721e057df02SPaul Mullowney E*/
17223f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat;
1723e057df02SPaul Mullowney 
1724e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1725e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[];
1726e057df02SPaul Mullowney 
17273f9c0db1SPaul Mullowney /*E
1728e057df02SPaul Mullowney     MatCUSPFormatOperation - indicates the operation of CUSP (GPU)
172936d62e41SPaul Mullowney     matrices whose operation should use a particular storage format.
1730e057df02SPaul Mullowney 
1731e057df02SPaul Mullowney     Not Collective
1732e057df02SPaul Mullowney 
17338468deeeSKarl Rupp +   MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
17348468deeeSKarl Rupp .   MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
17358468deeeSKarl Rupp .   MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
17368468deeeSKarl Rupp -   MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices
1737e057df02SPaul Mullowney 
1738e057df02SPaul Mullowney     Level: intermediate
1739e057df02SPaul Mullowney 
1740af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1741e057df02SPaul Mullowney 
1742e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat
1743e057df02SPaul Mullowney E*/
174436d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation;
1745e057df02SPaul Mullowney 
1746e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat);
1747e057df02SPaul Mullowney #endif
174890c53902SBarry Smith 
1749d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL)
1750d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1751d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1752d67ff14aSKarl Rupp #endif
1753d67ff14aSKarl Rupp 
175454efbe56SHong Zhang /*
175554efbe56SHong Zhang    PETSc interface to FFTW
175654efbe56SHong Zhang */
175754efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
1758014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1759014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
17602a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*);
176154efbe56SHong Zhang #endif
176254efbe56SHong Zhang 
1763014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1764014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1765014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1766014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1767014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1768014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
176919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1770014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1771014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1772d8588912SDave May 
17734325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
1774e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);
17754325cce7SMatthew G Knepley 
1776b541e6a4SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**);
177753dd7562SDmitry Karpeev 
1778c094ef40SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat);
1779c094ef40SMatthew G. Knepley 
1780539c167fSBarry Smith PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*);
1781539c167fSBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*);
1782539c167fSBarry Smith 
178323a5497aSJed Brown #endif
1784