xref: /petsc/include/petscmat.h (revision d975228c4e81a8c2e35afc29ff2de123c246fb54)
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*);
219014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*);
220c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*);
221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
223c0cdd4a1SDahai Guo 
224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2316d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
2336d7c1e57SBarry Smith 
23419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*);
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
236dedccee8SHong Zhang 
237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
238d0de2241SAndrew Spott PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*);
239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*);
240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS);
241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
2421472f72bSBarry Smith 
243978814f1SStefano Zampini #if defined(PETSC_HAVE_HYPRE)
244225daaf8SStefano Zampini PETSC_EXTERN PetscErrorCode MatCreateFromParCSR(void*,MatType,PetscCopyMode,Mat*);
245*d975228cSstefano_zampini PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
246978814f1SStefano Zampini #endif
247978814f1SStefano Zampini 
248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]);
2491d6018f0SLisandro Dalcin 
250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
252e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*);
25321c89e3eSBarry Smith 
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
259713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
26099cafbc1SBarry Smith 
2618ed539a5SBarry Smith /* ------------------------------------------------------------*/
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
26773a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);
26884cb2905SBarry Smith 
2692ef4de8bSBarry Smith /*S
2702ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
27162ef0227SBarry Smith         column of a matrix as indexed on an associated grid. These are arguments to MatSetStencil() and MatSetBlockStencil()
27262ef0227SBarry Smith 
27362ef0227SBarry 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).
27462ef0227SBarry 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.
27562ef0227SBarry Smith 
27662ef0227SBarry Smith    For stencil access to vectors see DMDAVecGetArray(), DMDAVecGetArrayF90().
277be479b30SJed Brown 
278be479b30SJed Brown    Fortran usage is different, see MatSetValuesStencil() for details.
2792ef4de8bSBarry Smith 
2802ef4de8bSBarry Smith    Level: beginner
2812ef4de8bSBarry Smith 
2822ef4de8bSBarry Smith   Concepts: matrix; linear operator
2832ef4de8bSBarry Smith 
28462ef0227SBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90()
2852ef4de8bSBarry Smith S*/
286435da068SBarry Smith typedef struct {
287c1ac3661SBarry Smith   PetscInt k,j,i,c;
288435da068SBarry Smith } MatStencil;
2892ef4de8bSBarry Smith 
290014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
291014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
292014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
293435da068SBarry Smith 
294d91e6319SBarry Smith /*E
295d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
29662ef0227SBarry Smith      to continue to add or insert values to it
297d91e6319SBarry Smith 
298d91e6319SBarry Smith     Level: beginner
299d91e6319SBarry Smith 
300d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
301d91e6319SBarry Smith E*/
3026d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
303014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
304014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
305014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);
3064f9c727eSBarry Smith 
3071d73ed98SBarry Smith 
30830de9b25SBarry Smith 
309d91e6319SBarry Smith /*E
310d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
311d91e6319SBarry Smith 
312d91e6319SBarry Smith     Level: beginner
313d91e6319SBarry Smith 
314af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
315d91e6319SBarry Smith 
316382164b0SBarry Smith    Developer Notes: Entries that are negative need not be called collectively by all processes.
317382164b0SBarry Smith 
318d91e6319SBarry Smith .seealso: MatSetOption()
319d91e6319SBarry Smith E*/
320c5e4d11fSDmitry Karpeev typedef enum {MAT_OPTION_MIN = -3,
321c5e4d11fSDmitry Karpeev               MAT_UNUSED_NONZERO_LOCATION_ERR = -2,
32292d4d306SBarry Smith               MAT_ROW_ORIENTED = -1,
323382164b0SBarry Smith               MAT_SYMMETRIC = 1,
324382164b0SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC = 2,
325382164b0SBarry Smith               MAT_NEW_DIAGONALS = 3,
326382164b0SBarry Smith               MAT_IGNORE_OFF_PROC_ENTRIES = 4,
327382164b0SBarry Smith               MAT_USE_HASH_TABLE = 5,
328382164b0SBarry Smith               MAT_KEEP_NONZERO_PATTERN = 6,
329382164b0SBarry Smith               MAT_IGNORE_ZERO_ENTRIES = 7,
330382164b0SBarry Smith               MAT_USE_INODES = 8,
331382164b0SBarry Smith               MAT_HERMITIAN = 9,
332382164b0SBarry Smith               MAT_SYMMETRY_ETERNAL = 10,
333c5e4d11fSDmitry Karpeev               MAT_NEW_NONZERO_LOCATION_ERR = 11,
334382164b0SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR = 12,
335382164b0SBarry Smith               MAT_ERROR_LOWER_TRIANGULAR = 13,
336382164b0SBarry Smith               MAT_GETROW_UPPERTRIANGULAR = 14,
337382164b0SBarry Smith               MAT_SPD = 15,
33892d4d306SBarry Smith               MAT_NO_OFF_PROC_ZERO_ROWS = 16,
33992d4d306SBarry Smith               MAT_NO_OFF_PROC_ENTRIES = 17,
34092d4d306SBarry Smith               MAT_NEW_NONZERO_LOCATIONS = 18,
341c5e4d11fSDmitry Karpeev               MAT_NEW_NONZERO_ALLOCATION_ERR = 19,
342c5e4d11fSDmitry Karpeev               MAT_SUBSET_OFF_PROC_ENTRIES = 20,
343c5e4d11fSDmitry Karpeev               MAT_OPTION_MAX = 21} MatOption;
344e057df02SPaul Mullowney 
345014dd563SJed Brown PETSC_EXTERN const char *MatOptions[];
346014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool);
3477d68702bSBarry Smith PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*);
34819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);
34984cb2905SBarry Smith 
350014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
351014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
3568c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
3578c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
35821e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*);
359bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3608c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
3618c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
363014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
364014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
365014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
36633d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat);
3671620fd73SBarry Smith 
368014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
369014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
380f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);
381f5edf698SHong Zhang 
382d91e6319SBarry Smith /*E
383d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
384d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
385d91e6319SBarry Smith 
386d91e6319SBarry Smith     Level: beginner
387d91e6319SBarry Smith 
388af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
389d91e6319SBarry Smith 
39070dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
39170dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
39270dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
39370dcbbb9SBarry Smith 
394d91e6319SBarry Smith .seealso: MatDuplicate()
395d91e6319SBarry Smith E*/
39670dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
3972e8a6d31SBarry Smith 
39819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
399014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);
40094a9d846SBarry Smith 
401cb5b572fSBarry Smith 
402014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
403014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
404014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
405014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
406014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
407014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
409014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
410014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
4117b80b807SBarry Smith 
4121a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4131a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
4141a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4151a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
416d4fbbf0eSBarry Smith 
417d91e6319SBarry Smith /*S
418d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
419d91e6319SBarry Smith 
420d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
421d91e6319SBarry Smith 
422d91e6319SBarry Smith    Level: intermediate
423d91e6319SBarry Smith 
424d91e6319SBarry Smith   Concepts: matrix^nonzero information
425d91e6319SBarry Smith 
426d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
427d91e6319SBarry Smith S*/
4284e220ebcSLois Curfman McInnes typedef struct {
429b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
430b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
431b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
432b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
433b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
434b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
435b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4364e220ebcSLois Curfman McInnes } MatInfo;
4374e220ebcSLois Curfman McInnes 
438d9274352SBarry Smith /*E
439d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
440d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
441d9274352SBarry Smith 
442d9274352SBarry Smith     Level: beginner
443d9274352SBarry Smith 
444af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
445d9274352SBarry Smith 
446d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
447d9274352SBarry Smith E*/
4487b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
449014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
450014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
452014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*);
459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);
461a52676ecSHong Zhang 
462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*);
463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*);
464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*);
465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*);
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*);
467a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
468a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
4697b80b807SBarry Smith 
470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*);
471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*);
472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
4797b80b807SBarry Smith 
480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
486566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);
4875ef9f2a5SBarry Smith 
488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
48953dd7562SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);
4967b80b807SBarry Smith 
497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5048efafbd8SBarry Smith 
505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
506aa1e27eaSFande Kong PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov);
507d2b2a242SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool);
5087b80b807SBarry Smith 
509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat);
51222440eb1SKris Buschelman 
5137bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
5147bab7c10SHong Zhang 
515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat);
518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat);
52122440eb1SKris Buschelman 
522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat);
525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat);
528bc011b1eSHong Zhang 
529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);
5311c741599SHong Zhang 
532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);
5347b80b807SBarry Smith 
535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
537a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*);
538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
544052efed2SBarry Smith 
545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
54790f02eecSBarry Smith 
548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
5512a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*);
55284cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);}
55353cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
5563a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
5579c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
558bd481603SBarry Smith 
559bd481603SBarry Smith /*MC
5601620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5611620fd73SBarry Smith 
5621620fd73SBarry Smith    Not collective
5631620fd73SBarry Smith 
564a9834a7dSSatish Balay    Synopsis:
565a9834a7dSSatish Balay      #include <petscmat.h>
566a9834a7dSSatish Balay      PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode)
567a9834a7dSSatish Balay 
5681620fd73SBarry Smith    Input Parameters:
5691620fd73SBarry Smith +  m - the matrix
5701620fd73SBarry Smith .  row - the row location of the entry
5711620fd73SBarry Smith .  col - the column location of the entry
5721620fd73SBarry Smith .  value - the value to insert
5731620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5741620fd73SBarry Smith 
5751620fd73SBarry Smith    Notes:
5761620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
5771620fd73SBarry Smith    values simultaneously if possible.
5781620fd73SBarry Smith 
5791620fd73SBarry Smith    Level: beginner
5801620fd73SBarry Smith 
5811620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
5821620fd73SBarry Smith M*/
5831620fd73SBarry 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);}
5841620fd73SBarry Smith 
5851620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
5861620fd73SBarry Smith 
5871620fd73SBarry 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);}
5881620fd73SBarry Smith 
5891620fd73SBarry Smith /*MC
5900d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
591bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
592bd481603SBarry Smith 
593bd481603SBarry Smith    Synopsis:
594aaa7dc30SBarry Smith    #include <petscmat.h>
595c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
596bd481603SBarry Smith 
597bd481603SBarry Smith    Collective on MPI_Comm
598bd481603SBarry Smith 
599bd481603SBarry Smith    Input Parameters:
600bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
601859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
602859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
603bd481603SBarry Smith 
604bd481603SBarry Smith    Output Parameters:
605bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
606bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
607bd481603SBarry Smith 
608bd481603SBarry Smith    Level: intermediate
609bd481603SBarry Smith 
610bd481603SBarry Smith    Notes:
611a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
612bd481603SBarry Smith 
6131d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
614bd481603SBarry Smith 
615bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
616bd481603SBarry Smith 
6171620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6181620fd73SBarry Smith 
619bd481603SBarry Smith   Concepts: preallocation^Matrix
620bd481603SBarry Smith 
621d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
622d6e23781SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock()
623bd481603SBarry Smith M*/
624c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6257c922b88SBarry Smith { \
626a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
6271795a4d1SJed Brown   _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \
6281795a4d1SJed Brown   __start = 0; __end = __start;                                         \
629c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
630a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6317c922b88SBarry Smith 
632bd481603SBarry Smith /*MC
633bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
634bd481603SBarry Smith        inserted using a local number of the rows and columns
635bd481603SBarry Smith 
636bd481603SBarry Smith    Synopsis:
637aaa7dc30SBarry Smith    #include <petscmat.h>
638c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
639bd481603SBarry Smith 
640bd481603SBarry Smith    Not Collective
641bd481603SBarry Smith 
642bd481603SBarry Smith    Input Parameters:
643784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
644bd481603SBarry Smith .  nrows - the number of rows indicated
6451d73ed98SBarry Smith .  rows - the indices of the rows
646784ac674SJed Brown .  cmap - the column mapping from local to global numbering
647bd481603SBarry Smith .  ncols - the number of columns in the matrix
648bd481603SBarry Smith .  cols - the columns indicated
649bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
650bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
651bd481603SBarry Smith 
652bd481603SBarry Smith    Level: intermediate
653bd481603SBarry Smith 
654bd481603SBarry Smith    Notes:
655a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
656bd481603SBarry Smith 
6571d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
658bd481603SBarry Smith 
659bd481603SBarry Smith   Concepts: preallocation^Matrix
660bd481603SBarry Smith 
661d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
662d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
663bd481603SBarry Smith M*/
664784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
665c4f061fbSSatish Balay {\
666c1ac3661SBarry Smith   PetscInt __l;\
667784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
668784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
669c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
670ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
671c4f061fbSSatish Balay   }\
672c4f061fbSSatish Balay }
673c4f061fbSSatish Balay 
674bd481603SBarry Smith /*MC
675d6e23781SBarry Smith    MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
676bd481603SBarry Smith        inserted using a local number of the rows and columns
677bd481603SBarry Smith 
678bd481603SBarry Smith    Synopsis:
679aaa7dc30SBarry Smith    #include <petscmat.h>
680d6e23781SBarry Smith    PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
681d6e23781SBarry Smith 
682d6e23781SBarry Smith    Not Collective
683d6e23781SBarry Smith 
684d6e23781SBarry Smith    Input Parameters:
685d6e23781SBarry Smith +  map - the row mapping from local numbering to global numbering
686d6e23781SBarry Smith .  nrows - the number of rows indicated
687d6e23781SBarry Smith .  rows - the indices of the rows
688d6e23781SBarry Smith .  cmap - the column mapping from local to global numbering
689d6e23781SBarry Smith .  ncols - the number of columns in the matrix
690d6e23781SBarry Smith .  cols - the columns indicated
691d6e23781SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
692d6e23781SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
693d6e23781SBarry Smith 
694d6e23781SBarry Smith    Level: intermediate
695d6e23781SBarry Smith 
696d6e23781SBarry Smith    Notes:
697d6e23781SBarry Smith     See Users-Manual: ch_performance for more details.
698d6e23781SBarry Smith 
699d6e23781SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
700d6e23781SBarry Smith 
701d6e23781SBarry Smith   Concepts: preallocation^Matrix
702d6e23781SBarry Smith 
703d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
704d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
705d6e23781SBarry Smith M*/
706d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
707d6e23781SBarry Smith {\
708d6e23781SBarry Smith   PetscInt __l;\
709d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
710d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
711d6e23781SBarry Smith   for (__l=0;__l<nrows;__l++) {\
712d6e23781SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
713d6e23781SBarry Smith   }\
714d6e23781SBarry Smith }
715d6e23781SBarry Smith 
716d6e23781SBarry Smith /*MC
717d6e23781SBarry Smith    MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
718d6e23781SBarry Smith        inserted using a local number of the rows and columns
719d6e23781SBarry Smith 
720d6e23781SBarry Smith    Synopsis:
721d6e23781SBarry Smith    #include <petscmat.h>
722d6e23781SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
723bd481603SBarry Smith 
724bd481603SBarry Smith    Not Collective
725bd481603SBarry Smith 
726bd481603SBarry Smith    Input Parameters:
727bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
728bd481603SBarry Smith .  nrows - the number of rows indicated
7291d73ed98SBarry Smith .  rows - the indices of the rows
730bd481603SBarry Smith .  ncols - the number of columns in the matrix
731bd481603SBarry Smith .  cols - the columns indicated
732bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
733bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
734bd481603SBarry Smith 
735bd481603SBarry Smith    Level: intermediate
736bd481603SBarry Smith 
737bd481603SBarry Smith    Notes:
738a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
739bd481603SBarry Smith 
740bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
741bd481603SBarry Smith 
742bd481603SBarry Smith   Concepts: preallocation^Matrix
743bd481603SBarry Smith 
744d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(),
745d6e23781SBarry Smith           MatPreallocateInitialize(),  MatPreallocateSetLocal()
746bd481603SBarry Smith M*/
747d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\
748d3d32019SBarry Smith {\
749c1ac3661SBarry Smith   PetscInt __l;\
750d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
751d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
752d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
753d6e23781SBarry Smith     _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
754d3d32019SBarry Smith   }\
755d3d32019SBarry Smith }
756bd481603SBarry Smith /*MC
757bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
758bd481603SBarry Smith        inserted using a local number of the rows and columns
759bd481603SBarry Smith 
760bd481603SBarry Smith    Synopsis:
761aaa7dc30SBarry Smith    #include <petscmat.h>
762c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
763bd481603SBarry Smith 
764bd481603SBarry Smith    Not Collective
765bd481603SBarry Smith 
766bd481603SBarry Smith    Input Parameters:
76764075487SBarry Smith +  row - the row
768bd481603SBarry Smith .  ncols - the number of columns in the matrix
769a50f8bf6SHong Zhang -  cols - the columns indicated
770a50f8bf6SHong Zhang 
771a50f8bf6SHong Zhang    Output Parameters:
772a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
773bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
774bd481603SBarry Smith 
775bd481603SBarry Smith    Level: intermediate
776bd481603SBarry Smith 
777bd481603SBarry Smith    Notes:
778a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
779bd481603SBarry Smith 
780bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
781bd481603SBarry Smith 
7821620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7831620fd73SBarry Smith 
784bd481603SBarry Smith   Concepts: preallocation^Matrix
785bd481603SBarry Smith 
786d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
787d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSetLocal()
788bd481603SBarry Smith M*/
789c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
790c1ac3661SBarry Smith { PetscInt __i; \
791e32f2f54SBarry 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);\
792e32f2f54SBarry 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);\
7937c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
79464075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
7957cc688d7SBarry Smith     else dnz[row - __rstart]++;\
7967c922b88SBarry Smith   }\
7977c922b88SBarry Smith }
7987c922b88SBarry Smith 
799bd481603SBarry Smith /*MC
800d6e23781SBarry Smith    MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
801bd481603SBarry Smith        inserted using a local number of the rows and columns
802bd481603SBarry Smith 
803bd481603SBarry Smith    Synopsis:
804aaa7dc30SBarry Smith    #include <petscmat.h>
805d6e23781SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
806bd481603SBarry Smith 
807bd481603SBarry Smith    Not Collective
808bd481603SBarry Smith 
809bd481603SBarry Smith    Input Parameters:
810bd481603SBarry Smith +  nrows - the number of rows indicated
8111d73ed98SBarry Smith .  rows - the indices of the rows
812bd481603SBarry Smith .  ncols - the number of columns in the matrix
813bd481603SBarry Smith .  cols - the columns indicated
814bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
815bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
816bd481603SBarry Smith 
817bd481603SBarry Smith    Level: intermediate
818bd481603SBarry Smith 
819bd481603SBarry Smith    Notes:
820a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
821bd481603SBarry Smith 
822bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
823bd481603SBarry Smith 
8241620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8251620fd73SBarry Smith 
826bd481603SBarry Smith   Concepts: preallocation^Matrix
827bd481603SBarry Smith 
828d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(),  MatPreallocateInitialize(),
829d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
830bd481603SBarry Smith M*/
831d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\
832c1ac3661SBarry Smith { PetscInt __i; \
833d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
834d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
835d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
836d3d32019SBarry Smith   }\
837d3d32019SBarry Smith }
838d3d32019SBarry Smith 
839bd481603SBarry Smith /*MC
84016371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
84116371a99SBarry Smith 
84216371a99SBarry Smith    Synopsis:
843aaa7dc30SBarry Smith    #include <petscmat.h>
84416371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
84516371a99SBarry Smith 
84616371a99SBarry Smith    Not Collective
84716371a99SBarry Smith 
84816371a99SBarry Smith    Input Parameters:
84916371a99SBarry Smith .  A - matrix
85016371a99SBarry Smith .  row - row where values exist (must be local to this process)
85116371a99SBarry Smith .  ncols - number of columns
85216371a99SBarry Smith .  cols - columns with nonzeros
85316371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
85416371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
85516371a99SBarry Smith 
85616371a99SBarry Smith    Level: intermediate
85716371a99SBarry Smith 
85816371a99SBarry Smith    Notes:
859a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
86016371a99SBarry Smith 
86116371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
86216371a99SBarry Smith 
86316371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
86416371a99SBarry Smith 
86516371a99SBarry Smith   Concepts: preallocation^Matrix
86616371a99SBarry Smith 
867d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
868d6e23781SBarry Smith           MatPreallocateSymmetricSetLocalBlock()
86916371a99SBarry Smith M*/
8700298fd71SBarry 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);}
87116371a99SBarry Smith 
87216371a99SBarry Smith 
87316371a99SBarry Smith /*MC
8740d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
875bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
876bd481603SBarry Smith 
877bd481603SBarry Smith    Synopsis:
878aaa7dc30SBarry Smith    #include <petscmat.h>
879c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
880bd481603SBarry Smith 
881bd481603SBarry Smith    Collective on MPI_Comm
882bd481603SBarry Smith 
883bd481603SBarry Smith    Input Parameters:
88416371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
885bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
886bd481603SBarry Smith 
887bd481603SBarry Smith    Level: intermediate
888bd481603SBarry Smith 
889bd481603SBarry Smith    Notes:
890a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
891bd481603SBarry Smith 
892bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
893bd481603SBarry Smith 
8941620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
8951620fd73SBarry Smith 
896bd481603SBarry Smith   Concepts: preallocation^Matrix
897bd481603SBarry Smith 
898d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
899d6e23781SBarry Smith           MatPreallocateSymmetricSetLocalBlock()
900bd481603SBarry Smith M*/
901a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9027c922b88SBarry Smith 
9037b80b807SBarry Smith /* Routines unique to particular data structures */
904014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
905698d4c6aSKris Buschelman 
906014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
907014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
908698d4c6aSKris Buschelman 
909014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
910014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
911014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
912014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
913014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
914014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
9157b80b807SBarry Smith 
916a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
917a96a251dSBarry Smith 
918014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
919014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
920014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
921273d9f13SBarry Smith 
922014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
923014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
924014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
925014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
926014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
927014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
928014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
929014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
930014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
931014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
9329230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
9339230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
934014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
935273d9f13SBarry Smith 
9362e1947a5SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
937014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt);
938014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
9391b807ce4Svictorle 
940014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
941014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
9422e8a6d31SBarry Smith 
943014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
9443a7fca6bSBarry Smith 
945014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
9467b80b807SBarry Smith /*
9477b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
94894b7f48cSBarry Smith   done through the KSP and PC interfaces.
9497b80b807SBarry Smith */
9507b80b807SBarry Smith 
95176bdecfbSBarry Smith /*J
9528f6c3df8SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering
953d9274352SBarry Smith 
954d9274352SBarry Smith    Level: beginner
955d9274352SBarry Smith 
956d9274352SBarry Smith .seealso: MatGetOrdering()
95776bdecfbSBarry Smith J*/
95819fd82e9SBarry Smith typedef const char* MatOrderingType;
9592692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
9602692d6eeSBarry Smith #define MATORDERINGND          "nd"
9612692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
9622692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
9632692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
9642692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
965510184a7SJed Brown #define MATORDERINGWBM         "wbm"
966fc1bef75SJed Brown #define MATORDERINGSPECTRAL    "spectral"
967312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
968b12f92e5SBarry Smith 
96919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
970140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
971bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
972140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList;
973d4fbbf0eSBarry Smith 
974014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
975fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);
976a2ce50c7SBarry Smith 
977d91e6319SBarry Smith /*S
978d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
979d90ac83dSHong Zhang 
980d90ac83dSHong Zhang    Level: beginner
981d90ac83dSHong Zhang 
982d90ac83dSHong Zhang S*/
983d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
9846a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[];
9855e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];
986d90ac83dSHong Zhang 
987d90ac83dSHong Zhang /*S
9889aa7eafdSHong Zhang     MatFactorError - indicates what type of error in matrix factor
9899aa7eafdSHong Zhang 
9909aa7eafdSHong Zhang     Level: beginner
9919aa7eafdSHong Zhang 
99200e125f8SBarry Smith     Developer Notes: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
99300e125f8SBarry Smith 
99400e125f8SBarry Smith .seealso: MatGetFactor()
9959aa7eafdSHong Zhang S*/
9965cd7cf9dSHong Zhang typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError;
9979aa7eafdSHong Zhang 
99800e125f8SBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*);
999b8b68cfdSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat);
10007b6c816cSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*);
100100e125f8SBarry Smith 
10029aa7eafdSHong Zhang /*S
1003422a814eSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization
1004b00f7748SHong Zhang 
100561cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
100661cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1007b00f7748SHong Zhang 
100815e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1009b00f7748SHong Zhang 
101061cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
101161cad860SBarry Smith 
1012b00f7748SHong Zhang    Level: developer
1013b00f7748SHong Zhang 
1014d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1015d7d82daaSBarry Smith           MatFactorInfoInitialize()
1016b00f7748SHong Zhang 
1017b00f7748SHong Zhang S*/
1018b00f7748SHong Zhang typedef struct {
101915e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
102085317021SBarry Smith   PetscReal     usedt;
102115e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1022b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
102315e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
102467e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1025348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1026bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1027bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
102815e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1029f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1030f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
103115e8a5b3SHong Zhang } MatFactorInfo;
1032ffa6d0a5SLois Curfman McInnes 
1033014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
10349a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
10359a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
10369a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
10379a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
10389a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
10399a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
10409a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
10419a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
10429a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
10439a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1044014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1045014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1046014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1047014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1048014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1049014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1050014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1051014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
10528e7ba810SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS);
10535a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*);
10545a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*);
10555a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat);
10565a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*);
10575a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec);
10585a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec);
1059014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
1060bb5a7306SBarry Smith 
1061d91e6319SBarry Smith /*E
1062d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1063bb1eb677SSatish Balay 
1064d91e6319SBarry Smith     Level: beginner
1065d91e6319SBarry Smith 
1066d9274352SBarry Smith    May be bitwise ORd together
1067d9274352SBarry Smith 
1068af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1069d91e6319SBarry Smith 
10704e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10714e7234bfSBarry Smith 
107241f059aeSBarry Smith .seealso: MatSOR()
1073d91e6319SBarry Smith E*/
1074ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1075ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1076ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
107784cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1078014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10798ed539a5SBarry Smith 
1080d4fbbf0eSBarry Smith /*
1081639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1082639f9d9dSBarry Smith */
1083b12f92e5SBarry Smith 
10847086a01eSPeter Brune /*S
10857086a01eSPeter Brune      MatColoring - Object for managing the coloring of matrices.
10867086a01eSPeter Brune 
10877086a01eSPeter Brune    Level: beginner
10887086a01eSPeter Brune 
10897086a01eSPeter Brune   Concepts: matrix, coloring
10907086a01eSPeter Brune 
10917086a01eSPeter Brune .seealso:  MatFDColoringCreate() ISColoring MatFDColoring
10927086a01eSPeter Brune S*/
10937086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring;
109476bdecfbSBarry Smith /*J
10958f6c3df8SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring
1096d9274352SBarry Smith 
1097d9274352SBarry Smith    Level: beginner
1098d9274352SBarry Smith 
10997086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring
110076bdecfbSBarry Smith J*/
11017086a01eSPeter Brune 
110219fd82e9SBarry Smith typedef const  char*           MatColoringType;
11035567d71aSPeter Brune #define MATCOLORINGJP      "jp"
1104b9acb617SPeter Brune #define MATCOLORINGPOWER   "power"
11052692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
11062692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
11072692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
11082692d6eeSBarry Smith #define MATCOLORINGID      "id"
11098121bdceSPeter Brune #define MATCOLORINGGREEDY  "greedy"
1110b12f92e5SBarry Smith 
11118ac6da0aSPeter Brune /*E
11128ac6da0aSPeter Brune    MatColoringWeightType - Type of weight scheme
11138ac6da0aSPeter Brune 
11148ac6da0aSPeter Brune     Not Collective
11158ac6da0aSPeter Brune 
11168ac6da0aSPeter Brune +   MAT_COLORING_RANDOM  - Random weights
11178ac6da0aSPeter Brune .   MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
11188ac6da0aSPeter Brune -   MAT_COLORING_LF      - Last-first weighting.
11198ac6da0aSPeter Brune 
11208ac6da0aSPeter Brune     Level: intermediate
11218ac6da0aSPeter Brune 
1122af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
11238ac6da0aSPeter Brune 
11248ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
11258ac6da0aSPeter Brune E*/
1126301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;
11278ac6da0aSPeter Brune 
1128335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1129d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*);
1130335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1131335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1132335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1133335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1134335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1135335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1136335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1137335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1138335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1139335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1140014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
11418ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType);
1142c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*);
11438ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm);
1144cb394dc5SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringTestValid(MatColoring,ISColoring);
1145639f9d9dSBarry Smith 
1146d9274352SBarry Smith /*S
1147d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1148d9274352SBarry Smith         and coloring
1149639f9d9dSBarry Smith 
1150d9274352SBarry Smith    Level: beginner
1151d9274352SBarry Smith 
1152d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1153d9274352SBarry Smith 
1154d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1155d9274352SBarry Smith S*/
1156e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1157639f9d9dSBarry Smith 
1158014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1159014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1160014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1161014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1162014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1163014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1164014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1165d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *);
1166014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1167014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1168f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1169f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1170f86b9fbaSHong Zhang 
1171b1683b59SHong Zhang 
1172b1683b59SHong Zhang /*S
1173b9af6bddSHong Zhang      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1174b1683b59SHong Zhang 
1175b1683b59SHong Zhang    Level: beginner
1176b1683b59SHong Zhang 
1177b1683b59SHong Zhang   Concepts: coloring, sparse matrix product
1178b1683b59SHong Zhang 
1179b9af6bddSHong Zhang .seealso:  MatTransposeColoringCreate()
1180b1683b59SHong Zhang S*/
1181b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1182b1683b59SHong Zhang 
1183014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1184014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1185014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1186014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1187b1683b59SHong Zhang 
1188639f9d9dSBarry Smith /*
11890752156aSBarry Smith     These routines are for partitioning matrices: currently used only
11903eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
11910752156aSBarry Smith */
1192ca161407SBarry Smith 
1193d9274352SBarry Smith /*S
1194d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1195d9274352SBarry Smith 
1196d9274352SBarry Smith    Level: beginner
1197d9274352SBarry Smith 
1198d9274352SBarry Smith   Concepts: partitioning
1199d9274352SBarry Smith 
1200743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1201d9274352SBarry Smith S*/
120291e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1203d9274352SBarry Smith 
120476bdecfbSBarry Smith /*J
12058f6c3df8SBarry Smith     MatPartitioningType - String with the name of a PETSc matrix partitioning
1206d9274352SBarry Smith 
1207d9274352SBarry Smith    Level: beginner
12080d04baf8SBarry Smith dm
1209b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
121076bdecfbSBarry Smith J*/
121119fd82e9SBarry Smith typedef const char* MatPartitioningType;
12122692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
1213aa1e27eaSFande Kong #define MATPARTITIONINGAVERAGE   "average"
12142692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
12152692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
12162692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
12172692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
121850d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
121988d2ac2bSFande Kong #define MATPARTITIONINGHIERARCH  "hierarch"
122071306c60Sjroman 
1221ca161407SBarry Smith 
1222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
122319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
12302aabb6bbSBarry Smith 
1231bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
123230de9b25SBarry Smith 
1233f1144a30SSatish Balay 
12342bad1931SBarry Smith 
1235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
123719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1238ca161407SBarry Smith 
123927538973SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
12420752156aSBarry Smith 
1243b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
12446a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1245b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
12466a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[];
1247b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
12486a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[];
1249b6956c12SJose Roman 
1250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
126171306c60Sjroman 
126271306c60Sjroman #define MP_PARTY_OPT "opt"
126371306c60Sjroman #define MP_PARTY_LIN "lin"
126471306c60Sjroman #define MP_PARTY_SCA "sca"
126571306c60Sjroman #define MP_PARTY_RAN "ran"
126671306c60Sjroman #define MP_PARTY_GBF "gbf"
126771306c60Sjroman #define MP_PARTY_GCF "gcf"
126871306c60Sjroman #define MP_PARTY_BUB "bub"
126971306c60Sjroman #define MP_PARTY_DEF "def"
1270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
127171306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
127271306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
127371306c60Sjroman #define MP_PARTY_NONE "no"
1274014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1275014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1276014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1277014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
127871306c60Sjroman 
127950d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
12806a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1281e0f1cffaSJose Roman 
1282014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1283014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1284014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1285014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
128671306c60Sjroman 
1287b43b03e9SMark F. Adams /*
12886294fa2bSFande Kong  * hierarchical partitioning
12896294fa2bSFande Kong  */
12904e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*);
12914e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*);
12924e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt);
12934e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);
12946294fa2bSFande Kong 
1295014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1296014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1297591294e4SBarry Smith 
12980752156aSBarry Smith /*
1299af0996ceSBarry Smith     If you add entries here you must also add them to petsc/finclude/petscmat.h
1300d4fbbf0eSBarry Smith */
13011c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13021c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13031c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13041c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13051c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13067c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13077c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13081c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13091c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13107c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13117c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13121c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13131c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1314a714c06dSBarry Smith                MATOP_SOR=13,
13151c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13161c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13171c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13181c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13191c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13201c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13211c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13221c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1323d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1324d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1325d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1326d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1327d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1328d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1329d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1330d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1331d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1332d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1333a5b7ff6bSBarry Smith                MATOP_GET_DIAGONAL_BLOCK=32,
13349d39f69dSJed Brown                /* MATOP_PLACEHOLDER_33=33, */
1335d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1336d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1337d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1338d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1339d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1340d519adbfSMatthew Knepley                MATOP_AXPY=39,
1341d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1342d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1343d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1344d519adbfSMatthew Knepley                MATOP_COPY=43,
1345d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1346d519adbfSMatthew Knepley                MATOP_SCALE=45,
1347d519adbfSMatthew Knepley                MATOP_SHIFT=46,
134835153367SBarry Smith                MATOP_DIAGONAL_SET=47,
13499d39f69dSJed Brown                MATOP_ZERO_ROWS_COLUMNS=48,
13509d39f69dSJed Brown                MATOP_SET_RANDOM=49,
1351d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1352d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1353d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1354d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1355d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1356d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1357d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1358d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1359d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1360d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1361d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1362d519adbfSMatthew Knepley                MATOP_VIEW=61,
1363d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
13647bab7c10SHong Zhang                MATOP_MATMAT_MULT=63,
13657bab7c10SHong Zhang                MATOP_MATMAT_MULT_SYMBOLIC=64,
13667bab7c10SHong Zhang                MATOP_MATMAT_MULT_NUMERIC=65,
1367d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1368d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1369d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1370d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1371d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1372d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1373d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
13749d39f69dSJed Brown                /* MATOP_PLACEHOLDER_73=73, */
1375d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1376d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1377d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1378bda6c4cbSBarry Smith                MATOP_MULT_CONSTRAINED=77,
1379bda6c4cbSBarry Smith                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
13809d39f69dSJed Brown                MATOP_FIND_ZERO_DIAGONALS=79,
1381d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1382d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1383d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1384d519adbfSMatthew Knepley                MATOP_LOAD=83,
1385d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1386d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1387d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1388820469bcSHong Zhang                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1389d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1390d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1391d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1392d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1393d519adbfSMatthew Knepley                MATOP_PTAP=92,
1394d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1395d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1396bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT=95,
1397bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1398bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
13999d39f69dSJed Brown                /* MATOP_PLACEHOLDER_98=98, */
14009d39f69dSJed Brown                /* MATOP_PLACEHOLDER_99=99, */
14019d39f69dSJed Brown                /* MATOP_PLACEHOLDER_100=100, */
14029d39f69dSJed Brown                /* MATOP_PLACEHOLDER_101=101, */
1403d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
14049d39f69dSJed Brown                /* MATOP_PLACEHOLDER_103=103, */
1405d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1406d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1407bda6c4cbSBarry Smith                MATOP_IMAGINARY_PART=106,
1408bda6c4cbSBarry Smith                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1409bda6c4cbSBarry Smith                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1410bda6c4cbSBarry Smith                MATOP_MAT_SOLVE=109,
14110e8d04f7SBarry Smith                MATOP_GET_REDUNDANT_MATRIX=110,
1412d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
14130e8d04f7SBarry Smith                MATOP_GET_COLUMN_VECTOR=112,
1414d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
14150e8d04f7SBarry Smith                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
141689c6957cSBarry Smith                MATOP_CREATE=115,
141789c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
14180e8d04f7SBarry Smith                MATOP_GET_LOCAL_SUB_MATRIX=117,
14190e8d04f7SBarry Smith                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1420eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
14210e8d04f7SBarry Smith                MATOP_HERMITIAN_TRANSPOSE=120,
14220e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
14230e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
14240e8d04f7SBarry Smith                MATOP_GET_MULTI_PROC_BLOCK=123,
14259d39f69dSJed Brown                MATOP_FIND_NONZERO_ROWS=124,
14260e8d04f7SBarry Smith                MATOP_GET_COLUMN_NORMS=125,
14279d39f69dSJed Brown                MATOP_INVERT_BLOCK_DIAGONAL=126,
14289d39f69dSJed Brown                /* MATOP_PLACEHOLDER_127=127, */
14290e8d04f7SBarry Smith                MATOP_GET_SUB_MATRICES_PARALLE=128,
14302b8ad9a3SHong Zhang                MATOP_SET_VALUES_BATCH=129,
14310e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT=130,
14320e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
1433e9b602ebSSatish Balay                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
14340e8d04f7SBarry Smith                MATOP_TRANSPOSE_COLORING_CREAT=133,
14350e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_SPT=134,
14360e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_DEN=135,
14370e8d04f7SBarry Smith                MATOP_RART=136,
14380e8d04f7SBarry Smith                MATOP_RART_SYMBOLIC=137,
14390e8d04f7SBarry Smith                MATOP_RART_NUMERIC=138,
1440e09a3074SHong Zhang                MATOP_SET_BLOCK_SIZES=139,
1441f9426fe0SMark Adams                MATOP_AYPX=140,
14421919a2e2SJed Brown                MATOP_RESIDUAL=141,
14439c8f2541SHong Zhang                MATOP_FDCOLORING_SETUP=142,
14449c8f2541SHong Zhang                MATOP_MPICONCATENATESEQ=144
1445fae171e0SBarry Smith              } MatOperation;
1446014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *);
1447014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1448014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1449014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1450112a2221SBarry Smith 
145190ace30eSBarry Smith /*
145290ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
145390ace30eSBarry Smith    stored in a universal format. By changing the format with
14546a9046bcSBarry Smith    PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
145590ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
145690ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1457f4403165SShri Abhyankar    read into matrices of the same type.
145890ace30eSBarry Smith */
145990ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
146090ace30eSBarry Smith 
1461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
1464b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*);
14651f4e1ec7SBarry Smith 
1466d9274352SBarry Smith /*S
1467d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1468d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1469d9274352SBarry Smith 
1470f7a9e4ceSBarry Smith    Level: advanced
1471d9274352SBarry Smith 
1472d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1473d9274352SBarry Smith 
14746e1639daSBarry Smith   Users manual sections:
14756e1639daSBarry Smith .   sec_singular
14766e1639daSBarry Smith 
1477d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1478d9274352SBarry Smith S*/
147974637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1480d9274352SBarry Smith 
1481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1484d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
14865fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *);
14875fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace);
1488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
149574637425SBarry Smith 
1496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
15003f1d51d7SBarry Smith 
1501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1504c4f061fbSSatish Balay 
1505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*);
1506b0a32e0cSBarry Smith 
1507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
150804f1ad80SBarry Smith 
1509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1522e884886fSBarry Smith 
15236370053bSBarry Smith /*S
15246370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
15256370053bSBarry Smith               Jacobian vector products
1526e884886fSBarry Smith 
15276370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
15286370053bSBarry Smith 
15296370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
15306370053bSBarry Smith 
15316370053bSBarry Smith     Level: developer
15326370053bSBarry Smith 
15336370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
15346370053bSBarry Smith S*/
1535e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1536e884886fSBarry Smith 
153776bdecfbSBarry Smith /*J
1538e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1539e884886fSBarry Smith 
1540e884886fSBarry Smith    Level: beginner
1541e884886fSBarry Smith 
1542e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
154376bdecfbSBarry Smith J*/
154419fd82e9SBarry Smith typedef const char* MatMFFDType;
1545e884886fSBarry Smith #define MATMFFD_DS  "ds"
1546e884886fSBarry Smith #define MATMFFD_WP  "wp"
1547e884886fSBarry Smith 
154819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1549bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));
1550e884886fSBarry Smith 
1551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool );
1553e884886fSBarry Smith 
155461ab5df0SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType);
155561ab5df0SBarry Smith 
1556014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1557014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
15587dbadf16SMatthew Knepley 
155997969023SHong Zhang /*
156097969023SHong Zhang    PETSc interface to MUMPS
156197969023SHong Zhang */
156297969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1564bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*);
1565b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
1566bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*);
1567bc6112feSHong Zhang 
1568ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1569ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1570ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1571ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
157297969023SHong Zhang #endif
157323a5497aSJed Brown 
1574d954db57SHong Zhang /*
1575d615f992SSatish Balay    PETSc interface to Mkl_Pardiso
1576b500a6b7SJose David Bermeol */
1577b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO
1578d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1579b500a6b7SJose David Bermeol #endif
1580b500a6b7SJose David Bermeol 
1581b500a6b7SJose David Bermeol /*
1582d305a81bSVasiliy Kozyrev    PETSc interface to Mkl_CPardiso
1583d305a81bSVasiliy Kozyrev */
1584d305a81bSVasiliy Kozyrev #ifdef PETSC_HAVE_MKL_CPARDISO
1585d305a81bSVasiliy Kozyrev PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt);
1586d305a81bSVasiliy Kozyrev #endif
1587d305a81bSVasiliy Kozyrev 
1588d305a81bSVasiliy Kozyrev /*
1589d954db57SHong Zhang    PETSc interface to SUPERLU
1590d954db57SHong Zhang */
1591d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1592014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1593d954db57SHong Zhang #endif
1594d954db57SHong Zhang 
1595fb785984SHong Zhang /*
1596fb785984SHong Zhang    PETSc interface to SUPERLU_DIST
1597fb785984SHong Zhang */
1598fb785984SHong Zhang #ifdef PETSC_HAVE_SUPERLU_DIST
1599fb785984SHong Zhang PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*);
1600fb785984SHong Zhang #endif
1601fb785984SHong Zhang 
1602575704cbSPieter Ghysels /*
1603575704cbSPieter Ghysels    PETSc interface to STRUMPACK
1604575704cbSPieter Ghysels */
1605575704cbSPieter Ghysels #ifdef PETSC_HAVE_STRUMPACK
1606575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat,PetscReal);
1607575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat,PetscInt);
1608575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool);
1609575704cbSPieter Ghysels #endif
1610575704cbSPieter Ghysels 
1611575704cbSPieter Ghysels 
1612aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA
16133f9c0db1SPaul Mullowney /*E
1614e057df02SPaul Mullowney     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
16152692e278SPaul Mullowney     matrices.
1616e057df02SPaul Mullowney 
1617e057df02SPaul Mullowney     Not Collective
1618e057df02SPaul Mullowney 
16198468deeeSKarl Rupp +   MAT_CUSPARSE_CSR - Compressed Sparse Row
16202692e278SPaul Mullowney .   MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later).
16212692e278SPaul Mullowney -   MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).
1622e057df02SPaul Mullowney 
1623e057df02SPaul Mullowney     Level: intermediate
1624e057df02SPaul Mullowney 
1625af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1626e057df02SPaul Mullowney 
1627e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1628e057df02SPaul Mullowney E*/
1629e057df02SPaul Mullowney 
16303f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;
1631e057df02SPaul Mullowney 
1632e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1633e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
1634e057df02SPaul Mullowney 
16353f9c0db1SPaul Mullowney /*E
1636e057df02SPaul Mullowney     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
16372692e278SPaul Mullowney     matrices whose operation should use a particular storage format.
1638e057df02SPaul Mullowney 
1639e057df02SPaul Mullowney     Not Collective
1640e057df02SPaul Mullowney 
16418468deeeSKarl Rupp +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
16428468deeeSKarl Rupp .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
16438468deeeSKarl Rupp .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
16448468deeeSKarl Rupp -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices
1645e057df02SPaul Mullowney 
1646e057df02SPaul Mullowney     Level: intermediate
1647e057df02SPaul Mullowney 
1648e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1649e057df02SPaul Mullowney E*/
165036d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;
1651e057df02SPaul Mullowney 
165210e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
165310e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1654e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
16559ae82921SPaul Mullowney #endif
16569ae82921SPaul Mullowney 
165790c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
1658014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1659014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1660e057df02SPaul Mullowney 
16613f9c0db1SPaul Mullowney /*E
1662e057df02SPaul Mullowney     MatCUSPStorageFormat - indicates the storage format for CUSP (GPU)
166336d62e41SPaul Mullowney     matrices.
1664e057df02SPaul Mullowney 
1665e057df02SPaul Mullowney     Not Collective
1666e057df02SPaul Mullowney 
16678468deeeSKarl Rupp +   MAT_CUSP_CSR - Compressed Sparse Row
16688468deeeSKarl Rupp .   MAT_CUSP_DIA - Diagonal
16698468deeeSKarl Rupp -   MAT_CUSP_ELL - Ellpack
1670e057df02SPaul Mullowney 
1671e057df02SPaul Mullowney     Level: intermediate
1672e057df02SPaul Mullowney 
1673af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1674e057df02SPaul Mullowney 
1675e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
1676e057df02SPaul Mullowney E*/
16773f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat;
1678e057df02SPaul Mullowney 
1679e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1680e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[];
1681e057df02SPaul Mullowney 
16823f9c0db1SPaul Mullowney /*E
1683e057df02SPaul Mullowney     MatCUSPFormatOperation - indicates the operation of CUSP (GPU)
168436d62e41SPaul Mullowney     matrices whose operation should use a particular storage format.
1685e057df02SPaul Mullowney 
1686e057df02SPaul Mullowney     Not Collective
1687e057df02SPaul Mullowney 
16888468deeeSKarl Rupp +   MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
16898468deeeSKarl Rupp .   MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
16908468deeeSKarl Rupp .   MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
16918468deeeSKarl Rupp -   MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices
1692e057df02SPaul Mullowney 
1693e057df02SPaul Mullowney     Level: intermediate
1694e057df02SPaul Mullowney 
1695af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1696e057df02SPaul Mullowney 
1697e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat
1698e057df02SPaul Mullowney E*/
169936d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation;
1700e057df02SPaul Mullowney 
1701e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat);
1702e057df02SPaul Mullowney #endif
170390c53902SBarry Smith 
1704d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL)
1705d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1706d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1707d67ff14aSKarl Rupp #endif
1708d67ff14aSKarl Rupp 
170954efbe56SHong Zhang /*
171054efbe56SHong Zhang    PETSc interface to FFTW
171154efbe56SHong Zhang */
171254efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
1713014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1714014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
17152a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*);
171654efbe56SHong Zhang #endif
171754efbe56SHong Zhang 
1718382fd914SHong Zhang /*
1719382fd914SHong Zhang    PETSc interface to ELEMENTAL
1720382fd914SHong Zhang */
1721db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL)
17222ef15aa2SHong Zhang #if defined(PETSC_USE_COMPLEX)
17232ef15aa2SHong Zhang typedef El::Complex<PetscReal> PetscElemScalar;
17242ef15aa2SHong Zhang #else
17252ef15aa2SHong Zhang typedef PetscScalar PetscElemScalar;
1726382fd914SHong Zhang #endif
1727607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
1728db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
1729db31f6deSJed Brown #endif
1730db31f6deSJed Brown 
1731014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1732014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1733014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1734014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1735014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1736014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
173719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1738014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1739014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1740d8588912SDave May 
17414325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
1742e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);
17434325cce7SMatthew G Knepley 
1744b541e6a4SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**);
174553dd7562SDmitry Karpeev 
1746c094ef40SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat);
1747c094ef40SMatthew G. Knepley 
1748539c167fSBarry Smith PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*);
1749539c167fSBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*);
1750539c167fSBarry Smith 
175123a5497aSJed Brown #endif
1752